12 #include "../general/forall.hpp"
73 static void PAGradientSetup2D(
const int Q1D,
75 const Array<double> &w,
80 const int NQ = Q1D*Q1D;
82 auto J =
Reshape(j.Read(), NQ, 2, 2, NE);
83 auto y =
Reshape(op.Write(), NQ, 2, 2, NE);
85 const bool const_c = c.Size() == 1;
86 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
91 for (
int q = 0; q < NQ; ++q)
93 const double J11 = J(q,0,0,e);
94 const double J12 = J(q,0,1,e);
95 const double J21 = J(q,1,0,e);
96 const double J22 = J(q,1,1,e);
98 const double Co = const_c ? C(0,0) : C(q,e);
99 y(q,0,0,e) = W[q] * Co * J22;
100 y(q,0,1,e) = W[q] * Co * -J12;
101 y(q,1,0,e) = W[q] * Co * -J21;
102 y(q,1,1,e) = W[q] * Co * J11;
108 static void PAGradientSetup3D(
const int Q1D,
110 const Array<double> &w,
115 const int NQ = Q1D*Q1D*Q1D;
117 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
118 auto y =
Reshape(op.Write(), NQ, 3, 3, NE);
120 const bool const_c = c.Size() == 1;
121 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
126 for (
int q = 0; q < NQ; ++q)
128 const double Co = const_c ? C(0,0) : C(q,e);
130 const double J11 = J(q,0,0,e);
131 const double J21 = J(q,1,0,e);
132 const double J31 = J(q,2,0,e);
133 const double J12 = J(q,0,1,e);
134 const double J22 = J(q,1,1,e);
135 const double J32 = J(q,2,1,e);
136 const double J13 = J(q,0,2,e);
137 const double J23 = J(q,1,2,e);
138 const double J33 = J(q,2,2,e);
139 const double cw = W[q] * Co;
141 const double A11 = (J22 * J33) - (J23 * J32);
142 const double A12 = (J32 * J13) - (J12 * J33);
143 const double A13 = (J12 * J23) - (J22 * J13);
144 const double A21 = (J31 * J23) - (J21 * J33);
145 const double A22 = (J11 * J33) - (J13 * J31);
146 const double A23 = (J21 * J13) - (J11 * J23);
147 const double A31 = (J21 * J32) - (J31 * J22);
148 const double A32 = (J31 * J12) - (J11 * J32);
149 const double A33 = (J11 * J22) - (J12 * J21);
151 y(q,0,0,e) = cw * A11;
152 y(q,0,1,e) = cw * A12;
153 y(q,0,2,e) = cw * A13;
154 y(q,1,0,e) = cw * A21;
155 y(q,1,1,e) = cw * A22;
156 y(q,1,2,e) = cw * A23;
157 y(q,2,0,e) = cw * A31;
158 y(q,2,1,e) = cw * A32;
159 y(q,2,2,e) = cw * A33;
164 static void PAGradientSetup(
const int dim,
169 const Array<double> &W,
174 if (dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAGradientSetup"); }
177 PAGradientSetup2D(Q1D, NE, W, J, COEFF, op);
181 PAGradientSetup3D(Q1D, NE, W, J, COEFF, op);
189 MFEM_ASSERT(trial_fes.
GetOrdering() == Ordering::byNODES,
190 "PA Only supports Ordering::byNODES!");
195 const IntegrationRule *ir = IntRule ? IntRule : &GetRule(trial_fe, test_fe,
197 const int dims = trial_fe.
GetDim();
198 const int dimsToStore = dims * dims;
200 dim = mesh->Dimension();
201 ne = trial_fes.
GetNE();
202 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
203 trial_maps = &trial_fe.
GetDofToQuad(*ir, DofToQuad::TENSOR);
204 trial_dofs1D = trial_maps->
ndof;
205 quad1D = trial_maps->nqpt;
206 test_maps = &test_fe.
GetDofToQuad(*ir, DofToQuad::TENSOR);
207 test_dofs1D = test_maps->
ndof;
208 MFEM_ASSERT(quad1D == test_maps->nqpt,
209 "PA requires test and trial space to have same number of quadrature points!");
222 coeff(0) = cQ->constant;
225 dynamic_cast<QuadratureFunctionCoefficient*>(Q))
228 MFEM_VERIFY(qFun.
Size() == ne*nq,
229 "Incompatible QuadratureFunction dimension \n");
232 "IntegrationRule used within integrator and in"
233 " QuadratureFunction appear to be different");
235 coeff.
MakeRef(const_cast<QuadratureFunction &>(qFun),0);
241 for (
int e = 0; e < ne; ++e)
244 for (
int q = 0; q < nq; ++q)
246 C(q,e) = Q->Eval(T, ir->
IntPoint(q));
251 PAGradientSetup(dim, trial_dofs1D, test_dofs1D, quad1D,
256 template<
int T_TR_D1D = 0,
int T_TE_D1D = 0,
int T_Q1D = 0>
257 static void PAGradientApply2D(
const int NE,
264 const int tr_d1d = 0,
265 const int te_d1d = 0,
268 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
269 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
270 const int Q1D = T_Q1D ? T_Q1D : q1d;
271 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
272 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
273 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
282 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
283 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
284 const int Q1D = T_Q1D ? T_Q1D : q1d;
287 constexpr
int max_TE_D1D = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
288 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
290 double grad[max_Q1D][max_Q1D][VDIM];
291 for (
int qy = 0; qy < Q1D; ++qy)
293 for (
int qx = 0; qx < Q1D; ++qx)
295 grad[qy][qx][0] = 0.0;
296 grad[qy][qx][1] = 0.0;
299 for (
int dy = 0; dy < TR_D1D; ++dy)
301 double gradX[max_Q1D][VDIM];
302 for (
int qx = 0; qx < Q1D; ++qx)
307 for (
int dx = 0; dx < TR_D1D; ++dx)
309 const double s = x(dx,dy,e);
310 for (
int qx = 0; qx < Q1D; ++qx)
312 gradX[qx][0] += s * G(qx,dx);
313 gradX[qx][1] += s * B(qx,dx);
316 for (
int qy = 0; qy < Q1D; ++qy)
318 const double wy = B(qy,dy);
319 const double wDy = G(qy,dy);
320 for (
int qx = 0; qx < Q1D; ++qx)
322 grad[qy][qx][0] += gradX[qx][0] * wy;
323 grad[qy][qx][1] += gradX[qx][1] * wDy;
328 for (
int qy = 0; qy < Q1D; ++qy)
330 for (
int qx = 0; qx < Q1D; ++qx)
332 const int q = qx + qy * Q1D;
333 const double gradX = grad[qy][qx][0];
334 const double gradY = grad[qy][qx][1];
336 grad[qy][qx][0] = gradX*op(q,0,0,e) + gradY*op(q,1,0,e);
337 grad[qy][qx][1] = gradX*op(q,0,1,e) + gradY*op(q,1,1,e);
341 for (
int qy = 0; qy < Q1D; ++qy)
343 double opX[max_TE_D1D][VDIM];
344 for (
int dx = 0; dx < TE_D1D; ++dx)
348 for (
int qx = 0; qx < Q1D; ++qx)
350 opX[dx][0] += Bt(dx,qx)*grad[qy][qx][0];
351 opX[dx][1] += Bt(dx,qx)*grad[qy][qx][1];
354 for (
int dy = 0; dy < TE_D1D; ++dy)
356 for (
int dx = 0; dx < TE_D1D; ++dx)
358 y(dx,dy,0,e) += Bt(dy,qy)*opX[dx][0];
359 y(dx,dy,1,e) += Bt(dy,qy)*opX[dx][1];
369 template<
int T_TR_D1D = 0,
int T_TE_D1D = 0,
int T_Q1D = 0>
370 static void PAGradientApplyTranspose2D(
const int NE,
371 const Array<double> &bt,
372 const Array<double> >,
373 const Array<double> &b,
377 const int tr_d1d = 0,
378 const int te_d1d = 0,
382 MFEM_ASSERT(
false,
"PAGradientApplyTranspose2D not implemented.");
386 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
387 static void PAGradientApply3D(
const int NE,
388 const Array<double> &b,
389 const Array<double> &g,
390 const Array<double> &bt,
398 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
399 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
400 const int Q1D = T_Q1D ? T_Q1D : q1d;
401 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
402 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
403 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
404 auto B =
Reshape(b.Read(), Q1D, TR_D1D);
405 auto G =
Reshape(g.Read(), Q1D, TR_D1D);
406 auto Bt =
Reshape(bt.Read(), TE_D1D, Q1D);
407 auto op =
Reshape(op_.Read(), Q1D*Q1D*Q1D, 3,3, NE);
408 auto x =
Reshape(x_.Read(), TR_D1D, TR_D1D, TR_D1D, NE);
409 auto y =
Reshape(y_.ReadWrite(), TE_D1D, TE_D1D, TE_D1D, 3, NE);
412 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
413 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
414 const int Q1D = T_Q1D ? T_Q1D : q1d;
417 constexpr
int max_TE_D1D = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
418 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
420 double grad[max_Q1D][max_Q1D][max_Q1D][VDIM];
421 for (
int qz = 0; qz < Q1D; ++qz)
423 for (
int qy = 0; qy < Q1D; ++qy)
425 for (
int qx = 0; qx < Q1D; ++qx)
427 grad[qz][qy][qx][0] = 0.0;
428 grad[qz][qy][qx][1] = 0.0;
429 grad[qz][qy][qx][2] = 0.0;
433 for (
int dz = 0; dz < TR_D1D; ++dz)
435 double gradXY[max_Q1D][max_Q1D][3];
436 for (
int qy = 0; qy < Q1D; ++qy)
438 for (
int qx = 0; qx < Q1D; ++qx)
440 gradXY[qy][qx][0] = 0.0;
441 gradXY[qy][qx][1] = 0.0;
442 gradXY[qy][qx][2] = 0.0;
445 for (
int dy = 0; dy < TR_D1D; ++dy)
447 double gradX[max_Q1D][2];
448 for (
int qx = 0; qx < Q1D; ++qx)
453 for (
int dx = 0; dx < TR_D1D; ++dx)
455 const double s = x(dx,dy,dz,e);
456 for (
int qx = 0; qx < Q1D; ++qx)
458 gradX[qx][0] += s * B(qx,dx);
459 gradX[qx][1] += s * G(qx,dx);
462 for (
int qy = 0; qy < Q1D; ++qy)
464 const double wy = B(qy,dy);
465 const double wDy = G(qy,dy);
466 for (
int qx = 0; qx < Q1D; ++qx)
468 const double wx = gradX[qx][0];
469 const double wDx = gradX[qx][1];
470 gradXY[qy][qx][0] += wDx * wy;
471 gradXY[qy][qx][1] += wx * wDy;
472 gradXY[qy][qx][2] += wx * wy;
476 for (
int qz = 0; qz < Q1D; ++qz)
478 const double wz = B(qz,dz);
479 const double wDz = G(qz,dz);
480 for (
int qy = 0; qy < Q1D; ++qy)
482 for (
int qx = 0; qx < Q1D; ++qx)
484 grad[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
485 grad[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
486 grad[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
492 for (
int qz = 0; qz < Q1D; ++qz)
494 for (
int qy = 0; qy < Q1D; ++qy)
496 for (
int qx = 0; qx < Q1D; ++qx)
498 const int q = qx + (qy + qz * Q1D) * Q1D;
499 const double gradX = grad[qz][qy][qx][0];
500 const double gradY = grad[qz][qy][qx][1];
501 const double gradZ = grad[qz][qy][qx][2];
503 grad[qz][qy][qx][0] = gradX*op(q,0,0,e) + gradY*op(q,1,0,e) + gradZ*op(q,2,0,e);
504 grad[qz][qy][qx][1] = gradX*op(q,0,1,e) + gradY*op(q,1,1,e) + gradZ*op(q,2,1,e);
505 grad[qz][qy][qx][2] = gradX*op(q,0,2,e) + gradY*op(q,1,2,e) + gradZ*op(q,2,2,e);
510 for (
int qz = 0; qz < Q1D; ++qz)
512 double opXY[max_TE_D1D][max_TE_D1D][VDIM];
513 for (
int dy = 0; dy < TE_D1D; ++dy)
515 for (
int dx = 0; dx < TE_D1D; ++dx)
517 opXY[dy][dx][0] = 0.0;
518 opXY[dy][dx][1] = 0.0;
519 opXY[dy][dx][2] = 0.0;
522 for (
int qy = 0; qy < Q1D; ++qy)
524 double opX[max_TE_D1D][VDIM];
525 for (
int dx = 0; dx < TE_D1D; ++dx)
530 for (
int qx = 0; qx < Q1D; ++qx)
532 opX[dx][0] += Bt(dx,qx)*grad[qz][qy][qx][0];
533 opX[dx][1] += Bt(dx,qx)*grad[qz][qy][qx][1];
534 opX[dx][2] += Bt(dx,qx)*grad[qz][qy][qx][2];
537 for (
int dy = 0; dy < TE_D1D; ++dy)
539 for (
int dx = 0; dx < TE_D1D; ++dx)
541 opXY[dy][dx][0] += Bt(dy,qy)*opX[dx][0];
542 opXY[dy][dx][1] += Bt(dy,qy)*opX[dx][1];
543 opXY[dy][dx][2] += Bt(dy,qy)*opX[dx][2];
547 for (
int dz = 0; dz < TE_D1D; ++dz)
549 for (
int dy = 0; dy < TE_D1D; ++dy)
551 for (
int dx = 0; dx < TE_D1D; ++dx)
553 y(dx,dy,dz,0,e) += Bt(dz,qz)*opXY[dy][dx][0];
554 y(dx,dy,dz,1,e) += Bt(dz,qz)*opXY[dy][dx][1];
555 y(dx,dy,dz,2,e) += Bt(dz,qz)*opXY[dy][dx][2];
565 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
566 static void PAGradientApplyTranspose3D(
const int NE,
567 const Array<double> &bt,
568 const Array<double> >,
569 const Array<double> &b,
577 MFEM_ASSERT(
false,
"Gradient PA Apply Transpose 3D not implemented.");
581 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
582 static void SmemPAGradientApply3D(
const int NE,
583 const Array<double> &b_,
584 const Array<double> &g_,
585 const Array<double> &bt_,
589 const int tr_d1d = 0,
590 const int te_d1d = 0,
593 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
594 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
595 const int Q1D = T_Q1D ? T_Q1D : q1d;
597 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
598 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
599 MFEM_VERIFY(TR_D1D <= Q1D,
"");
600 MFEM_VERIFY(TE_D1D <= Q1D,
"");
601 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
603 auto b =
Reshape(b_.Read(), Q1D, TR_D1D);
604 auto g =
Reshape(g_.Read(), Q1D, TR_D1D);
605 auto bt =
Reshape(bt_.Read(), TE_D1D, Q1D);
606 auto D =
Reshape(d_.Read(), Q1D*Q1D*Q1D, 3, 3, NE);
607 auto x =
Reshape(x_.Read(), TR_D1D, TR_D1D, TR_D1D, NE);
608 auto y =
Reshape(y_.ReadWrite(), TE_D1D, TE_D1D, TE_D1D, 3, NE);
610 MFEM_FORALL_3D(e, NE, (Q1D>8)?8:Q1D, (Q1D>8)?8:Q1D, (Q1D>8)?8:Q1D,
612 const int tidz = MFEM_THREAD_ID(z);
613 const int D1DR = T_TR_D1D ? T_TR_D1D : tr_d1d;
614 const int D1DE = T_TE_D1D ? T_TE_D1D : te_d1d;
615 const int Q1D = T_Q1D ? T_Q1D : q1d;
616 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
617 constexpr
int MD1R = T_TR_D1D ? T_TR_D1D :
MAX_D1D;
618 constexpr
int MD1E = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
619 constexpr
int MD1 = MD1E > MD1R ? MD1E : MD1R;
620 constexpr
int MDQ = MQ1 > MD1 ? MQ1 : MD1;
621 MFEM_SHARED
double sBG[2][MQ1*MD1];
622 double (*B)[MD1] = (double (*)[MD1]) (sBG+0);
623 double (*G)[MD1] = (double (*)[MD1]) (sBG+1);
624 double (*Bt)[MQ1] = (double (*)[MQ1]) (sBG+0);
625 MFEM_SHARED
double sm0[3][MDQ*MDQ*MDQ];
626 MFEM_SHARED
double sm1[3][MDQ*MDQ*MDQ];
627 double (*X)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+2);
628 double (*DDQ0)[MD1][MQ1] = (double (*)[MD1][MQ1]) (sm0+0);
629 double (*DDQ1)[MD1][MQ1] = (double (*)[MD1][MQ1]) (sm0+1);
631 double (*DQQ0)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+0);
632 double (*DQQ1)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+1);
633 double (*DQQ2)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+2);
635 double (*QQQ0)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+0);
636 double (*QQQ1)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+1);
637 double (*QQQ2)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+2);
639 double (*QQD0)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+0);
640 double (*QQD1)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+1);
641 double (*QQD2)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+2);
643 double (*QDD0)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+0);
644 double (*QDD1)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+1);
645 double (*QDD2)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+2);
646 MFEM_FOREACH_THREAD(dz,z,D1DR)
648 MFEM_FOREACH_THREAD(dy,y,D1DR)
650 MFEM_FOREACH_THREAD(dx,x,D1DR)
652 X[dz][dy][dx] = x(dx,dy,dz,e);
658 MFEM_FOREACH_THREAD(d,y,D1DR)
660 MFEM_FOREACH_THREAD(q,x,Q1D)
668 MFEM_FOREACH_THREAD(dz,z,D1DR)
670 MFEM_FOREACH_THREAD(dy,y,D1DR)
672 MFEM_FOREACH_THREAD(qx,x,Q1D)
676 for (
int dx = 0; dx < D1DR; ++dx)
678 const double coord = X[dz][dy][dx];
679 u += coord * B[qx][dx];
680 v += coord * G[qx][dx];
682 DDQ0[dz][dy][qx] =
u;
683 DDQ1[dz][dy][qx] = v;
688 MFEM_FOREACH_THREAD(dz,z,D1DR)
690 MFEM_FOREACH_THREAD(qy,y,Q1D)
692 MFEM_FOREACH_THREAD(qx,x,Q1D)
697 for (
int dy = 0; dy < D1DR; ++dy)
699 u += DDQ1[dz][dy][qx] * B[qy][dy];
700 v += DDQ0[dz][dy][qx] * G[qy][dy];
701 w += DDQ0[dz][dy][qx] * B[qy][dy];
703 DQQ0[dz][qy][qx] =
u;
704 DQQ1[dz][qy][qx] = v;
705 DQQ2[dz][qy][qx] = w;
710 MFEM_FOREACH_THREAD(qz,z,Q1D)
712 MFEM_FOREACH_THREAD(qy,y,Q1D)
714 MFEM_FOREACH_THREAD(qx,x,Q1D)
719 for (
int dz = 0; dz < D1DR; ++dz)
721 u += DQQ0[dz][qy][qx] * B[qz][dz];
722 v += DQQ1[dz][qy][qx] * B[qz][dz];
723 w += DQQ2[dz][qy][qx] * G[qz][dz];
725 QQQ0[qz][qy][qx] =
u;
726 QQQ1[qz][qy][qx] = v;
727 QQQ2[qz][qy][qx] = w;
732 MFEM_FOREACH_THREAD(qz,z,Q1D)
734 MFEM_FOREACH_THREAD(qy,y,Q1D)
736 MFEM_FOREACH_THREAD(qx,x,Q1D)
738 const int q = qx + (qy + qz * Q1D) * Q1D;
739 const double gX = QQQ0[qz][qy][qx];
740 const double gY = QQQ1[qz][qy][qx];
741 const double gZ = QQQ2[qz][qy][qx];
742 QQQ0[qz][qy][qx] = (D(q,0,0,e)*gX) + (D(q,1,0,e)*gY) + (D(q,2,0,e)*gZ);
743 QQQ1[qz][qy][qx] = (D(q,0,1,e)*gX) + (D(q,1,1,e)*gY) + (D(q,2,1,e)*gZ);
744 QQQ2[qz][qy][qx] = (D(q,0,2,e)*gX) + (D(q,1,2,e)*gY) + (D(q,2,2,e)*gZ);
751 MFEM_FOREACH_THREAD(d,y,D1DE)
753 MFEM_FOREACH_THREAD(q,x,Q1D)
760 MFEM_FOREACH_THREAD(qz,z,Q1D)
762 MFEM_FOREACH_THREAD(qy,y,Q1D)
764 MFEM_FOREACH_THREAD(dx,x,D1DE)
769 for (
int qx = 0; qx < Q1D; ++qx)
771 u += QQQ0[qz][qy][qx] * Bt[dx][qx];
772 v += QQQ1[qz][qy][qx] * Bt[dx][qx];
773 w += QQQ2[qz][qy][qx] * Bt[dx][qx];
775 QQD0[qz][qy][dx] =
u;
776 QQD1[qz][qy][dx] = v;
777 QQD2[qz][qy][dx] = w;
782 MFEM_FOREACH_THREAD(qz,z,Q1D)
784 MFEM_FOREACH_THREAD(dy,y,D1DE)
786 MFEM_FOREACH_THREAD(dx,x,D1DE)
791 for (
int qy = 0; qy < Q1D; ++qy)
793 u += QQD0[qz][qy][dx] * Bt[dy][qy];
794 v += QQD1[qz][qy][dx] * Bt[dy][qy];
795 w += QQD2[qz][qy][dx] * Bt[dy][qy];
797 QDD0[qz][dy][dx] =
u;
798 QDD1[qz][dy][dx] = v;
799 QDD2[qz][dy][dx] = w;
804 MFEM_FOREACH_THREAD(dz,z,D1DE)
806 MFEM_FOREACH_THREAD(dy,y,D1DE)
808 MFEM_FOREACH_THREAD(dx,x,D1DE)
813 for (
int qz = 0; qz < Q1D; ++qz)
815 u += QDD0[qz][dy][dx] * Bt[dz][qz];
816 v += QDD1[qz][dy][dx] * Bt[dz][qz];
817 w += QDD2[qz][dy][dx] * Bt[dz][qz];
819 y(dx,dy,dz,0,e) +=
u;
820 y(dx,dy,dz,1,e) += v;
821 y(dx,dy,dz,2,e) += w;
828 static void PAGradientApply(
const int dim,
833 const Array<double> &B,
834 const Array<double> &G,
835 const Array<double> &Bt,
839 bool transpose=
false)
844 return PAGradientApply2D(NE,B,G,Bt,op,x,y,TR_D1D,TE_D1D,Q1D);
848 return PAGradientApply3D(NE,B,G,Bt,op,x,y,TR_D1D,TE_D1D,Q1D);
850 MFEM_ABORT(
"Unknown kernel.");
856 PAGradientApply(dim, trial_dofs1D, test_dofs1D, quad1D, ne,
857 trial_maps->B, trial_maps->G, test_maps->Bt, pa_data, x, y,
862 void GradientIntegrator::AddMultTransposePA(
const Vector &x,
Vector &y)
const
864 MFEM_ABORT(
"PA Gradient AddMultTransposePA not implemented.");
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
Ordering::Type GetOrdering() const
Return the ordering method.
int GetDim() const
Returns the reference space dimension for the finite element.
void trans(const Vector &u, Vector &x)
Class for an integration rule - an Array of IntegrationPoint.
A coefficient that is constant across space and time.
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
void SetSize(int s)
Resize the vector to size s.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
int Size() const
Returns the size of the vector.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
int GetNE() const
Returns number of elements in the mesh.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
Mesh * GetMesh() const
Returns the mesh.
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
double u(const Vector &xvec)
Class representing a function through its values (scalar or vector) at quadrature points...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).