12 #include "../general/forall.hpp"
74 static void PAGradientSetup2D(
const int Q1D,
76 const Array<double> &w,
81 const int NQ = Q1D*Q1D;
83 auto J =
Reshape(j.Read(), NQ, 2, 2, NE);
84 auto y =
Reshape(op.Write(), NQ, 2, 2, NE);
86 const bool const_c = c.Size() == 1;
87 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
92 for (
int q = 0; q < NQ; ++q)
94 const double J11 = J(q,0,0,e);
95 const double J12 = J(q,0,1,e);
96 const double J21 = J(q,1,0,e);
97 const double J22 = J(q,1,1,e);
99 const double Co = const_c ? C(0,0) : C(q,e);
100 y(q,0,0,e) = W[q] * Co * J22;
101 y(q,0,1,e) = W[q] * Co * -J12;
102 y(q,1,0,e) = W[q] * Co * -J21;
103 y(q,1,1,e) = W[q] * Co * J11;
109 static void PAGradientSetup3D(
const int Q1D,
111 const Array<double> &w,
116 const int NQ = Q1D*Q1D*Q1D;
118 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
119 auto y =
Reshape(op.Write(), NQ, 3, 3, NE);
121 const bool const_c = c.Size() == 1;
122 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
127 for (
int q = 0; q < NQ; ++q)
129 const double Co = const_c ? C(0,0) : C(q,e);
131 const double J11 = J(q,0,0,e);
132 const double J21 = J(q,1,0,e);
133 const double J31 = J(q,2,0,e);
134 const double J12 = J(q,0,1,e);
135 const double J22 = J(q,1,1,e);
136 const double J32 = J(q,2,1,e);
137 const double J13 = J(q,0,2,e);
138 const double J23 = J(q,1,2,e);
139 const double J33 = J(q,2,2,e);
140 const double cw = W[q] * Co;
142 const double A11 = (J22 * J33) - (J23 * J32);
143 const double A12 = (J32 * J13) - (J12 * J33);
144 const double A13 = (J12 * J23) - (J22 * J13);
145 const double A21 = (J31 * J23) - (J21 * J33);
146 const double A22 = (J11 * J33) - (J13 * J31);
147 const double A23 = (J21 * J13) - (J11 * J23);
148 const double A31 = (J21 * J32) - (J31 * J22);
149 const double A32 = (J31 * J12) - (J11 * J32);
150 const double A33 = (J11 * J22) - (J12 * J21);
152 y(q,0,0,e) = cw * A11;
153 y(q,0,1,e) = cw * A12;
154 y(q,0,2,e) = cw * A13;
155 y(q,1,0,e) = cw * A21;
156 y(q,1,1,e) = cw * A22;
157 y(q,1,2,e) = cw * A23;
158 y(q,2,0,e) = cw * A31;
159 y(q,2,1,e) = cw * A32;
160 y(q,2,2,e) = cw * A33;
165 static void PAGradientSetup(
const int dim,
170 const Array<double> &W,
175 if (dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAGradientSetup"); }
178 PAGradientSetup2D(Q1D, NE, W, J, COEFF, op);
182 PAGradientSetup3D(Q1D, NE, W, J, COEFF, op);
190 MFEM_ASSERT(trial_fes.
GetOrdering() == Ordering::byNODES,
191 "PA Only supports Ordering::byNODES!");
198 const int dims = trial_fe.
GetDim();
199 const int dimsToStore = dims * dims;
201 dim = mesh->Dimension();
202 ne = trial_fes.
GetNE();
203 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
204 trial_maps = &trial_fe.
GetDofToQuad(*ir, DofToQuad::TENSOR);
205 trial_dofs1D = trial_maps->
ndof;
206 quad1D = trial_maps->nqpt;
207 test_maps = &test_fe.
GetDofToQuad(*ir, DofToQuad::TENSOR);
208 test_dofs1D = test_maps->
ndof;
209 MFEM_ASSERT(quad1D == test_maps->nqpt,
210 "PA requires test and trial space to have same number of quadrature points!");
216 PAGradientSetup(dim, trial_dofs1D, test_dofs1D, quad1D,
217 ne, ir->
GetWeights(), geom->J, coeff, pa_data);
221 template<
int T_TR_D1D = 0,
int T_TE_D1D = 0,
int T_Q1D = 0>
222 static void PAGradientApply2D(
const int NE,
229 const int tr_d1d = 0,
230 const int te_d1d = 0,
233 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
234 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
235 const int Q1D = T_Q1D ? T_Q1D : q1d;
236 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
237 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
238 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
247 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
248 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
249 const int Q1D = T_Q1D ? T_Q1D : q1d;
252 constexpr
int max_TE_D1D = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
253 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
255 double grad[max_Q1D][max_Q1D][VDIM];
256 for (
int qy = 0; qy < Q1D; ++qy)
258 for (
int qx = 0; qx < Q1D; ++qx)
260 grad[qy][qx][0] = 0.0;
261 grad[qy][qx][1] = 0.0;
264 for (
int dy = 0; dy < TR_D1D; ++dy)
266 double gradX[max_Q1D][VDIM];
267 for (
int qx = 0; qx < Q1D; ++qx)
272 for (
int dx = 0; dx < TR_D1D; ++dx)
274 const double s = x(dx,dy,e);
275 for (
int qx = 0; qx < Q1D; ++qx)
277 gradX[qx][0] += s * G(qx,dx);
278 gradX[qx][1] += s * B(qx,dx);
281 for (
int qy = 0; qy < Q1D; ++qy)
283 const double wy = B(qy,dy);
284 const double wDy = G(qy,dy);
285 for (
int qx = 0; qx < Q1D; ++qx)
287 grad[qy][qx][0] += gradX[qx][0] * wy;
288 grad[qy][qx][1] += gradX[qx][1] * wDy;
293 for (
int qy = 0; qy < Q1D; ++qy)
295 for (
int qx = 0; qx < Q1D; ++qx)
297 const int q = qx + qy * Q1D;
298 const double gradX = grad[qy][qx][0];
299 const double gradY = grad[qy][qx][1];
301 grad[qy][qx][0] = gradX*op(q,0,0,e) + gradY*op(q,1,0,e);
302 grad[qy][qx][1] = gradX*op(q,0,1,e) + gradY*op(q,1,1,e);
306 for (
int qy = 0; qy < Q1D; ++qy)
308 double opX[max_TE_D1D][VDIM];
309 for (
int dx = 0; dx < TE_D1D; ++dx)
313 for (
int qx = 0; qx < Q1D; ++qx)
315 opX[dx][0] += Bt(dx,qx)*grad[qy][qx][0];
316 opX[dx][1] += Bt(dx,qx)*grad[qy][qx][1];
319 for (
int dy = 0; dy < TE_D1D; ++dy)
321 for (
int dx = 0; dx < TE_D1D; ++dx)
323 y(dx,dy,0,e) += Bt(dy,qy)*opX[dx][0];
324 y(dx,dy,1,e) += Bt(dy,qy)*opX[dx][1];
334 template<
int T_TR_D1D = 0,
int T_TE_D1D = 0,
int T_Q1D = 0>
335 static void PAGradientApplyTranspose2D(
const int NE,
336 const Array<double> &bt,
337 const Array<double> >,
338 const Array<double> &b,
342 const int tr_d1d = 0,
343 const int te_d1d = 0,
347 MFEM_ASSERT(
false,
"PAGradientApplyTranspose2D not implemented.");
351 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
352 static void PAGradientApply3D(
const int NE,
353 const Array<double> &b,
354 const Array<double> &g,
355 const Array<double> &bt,
363 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
364 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
365 const int Q1D = T_Q1D ? T_Q1D : q1d;
366 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
367 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
368 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
369 auto B =
Reshape(b.Read(), Q1D, TR_D1D);
370 auto G =
Reshape(g.Read(), Q1D, TR_D1D);
371 auto Bt =
Reshape(bt.Read(), TE_D1D, Q1D);
372 auto op =
Reshape(op_.Read(), Q1D*Q1D*Q1D, 3,3, NE);
373 auto x =
Reshape(x_.Read(), TR_D1D, TR_D1D, TR_D1D, NE);
374 auto y =
Reshape(y_.ReadWrite(), TE_D1D, TE_D1D, TE_D1D, 3, NE);
377 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
378 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
379 const int Q1D = T_Q1D ? T_Q1D : q1d;
382 constexpr
int max_TE_D1D = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
383 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
385 double grad[max_Q1D][max_Q1D][max_Q1D][VDIM];
386 for (
int qz = 0; qz < Q1D; ++qz)
388 for (
int qy = 0; qy < Q1D; ++qy)
390 for (
int qx = 0; qx < Q1D; ++qx)
392 grad[qz][qy][qx][0] = 0.0;
393 grad[qz][qy][qx][1] = 0.0;
394 grad[qz][qy][qx][2] = 0.0;
398 for (
int dz = 0; dz < TR_D1D; ++dz)
400 double gradXY[max_Q1D][max_Q1D][3];
401 for (
int qy = 0; qy < Q1D; ++qy)
403 for (
int qx = 0; qx < Q1D; ++qx)
405 gradXY[qy][qx][0] = 0.0;
406 gradXY[qy][qx][1] = 0.0;
407 gradXY[qy][qx][2] = 0.0;
410 for (
int dy = 0; dy < TR_D1D; ++dy)
412 double gradX[max_Q1D][2];
413 for (
int qx = 0; qx < Q1D; ++qx)
418 for (
int dx = 0; dx < TR_D1D; ++dx)
420 const double s = x(dx,dy,dz,e);
421 for (
int qx = 0; qx < Q1D; ++qx)
423 gradX[qx][0] += s * B(qx,dx);
424 gradX[qx][1] += s * G(qx,dx);
427 for (
int qy = 0; qy < Q1D; ++qy)
429 const double wy = B(qy,dy);
430 const double wDy = G(qy,dy);
431 for (
int qx = 0; qx < Q1D; ++qx)
433 const double wx = gradX[qx][0];
434 const double wDx = gradX[qx][1];
435 gradXY[qy][qx][0] += wDx * wy;
436 gradXY[qy][qx][1] += wx * wDy;
437 gradXY[qy][qx][2] += wx * wy;
441 for (
int qz = 0; qz < Q1D; ++qz)
443 const double wz = B(qz,dz);
444 const double wDz = G(qz,dz);
445 for (
int qy = 0; qy < Q1D; ++qy)
447 for (
int qx = 0; qx < Q1D; ++qx)
449 grad[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
450 grad[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
451 grad[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
457 for (
int qz = 0; qz < Q1D; ++qz)
459 for (
int qy = 0; qy < Q1D; ++qy)
461 for (
int qx = 0; qx < Q1D; ++qx)
463 const int q = qx + (qy + qz * Q1D) * Q1D;
464 const double gradX = grad[qz][qy][qx][0];
465 const double gradY = grad[qz][qy][qx][1];
466 const double gradZ = grad[qz][qy][qx][2];
468 grad[qz][qy][qx][0] = gradX*op(q,0,0,e) + gradY*op(q,1,0,e) + gradZ*op(q,2,0,e);
469 grad[qz][qy][qx][1] = gradX*op(q,0,1,e) + gradY*op(q,1,1,e) + gradZ*op(q,2,1,e);
470 grad[qz][qy][qx][2] = gradX*op(q,0,2,e) + gradY*op(q,1,2,e) + gradZ*op(q,2,2,e);
475 for (
int qz = 0; qz < Q1D; ++qz)
477 double opXY[max_TE_D1D][max_TE_D1D][VDIM];
478 for (
int dy = 0; dy < TE_D1D; ++dy)
480 for (
int dx = 0; dx < TE_D1D; ++dx)
482 opXY[dy][dx][0] = 0.0;
483 opXY[dy][dx][1] = 0.0;
484 opXY[dy][dx][2] = 0.0;
487 for (
int qy = 0; qy < Q1D; ++qy)
489 double opX[max_TE_D1D][VDIM];
490 for (
int dx = 0; dx < TE_D1D; ++dx)
495 for (
int qx = 0; qx < Q1D; ++qx)
497 opX[dx][0] += Bt(dx,qx)*grad[qz][qy][qx][0];
498 opX[dx][1] += Bt(dx,qx)*grad[qz][qy][qx][1];
499 opX[dx][2] += Bt(dx,qx)*grad[qz][qy][qx][2];
502 for (
int dy = 0; dy < TE_D1D; ++dy)
504 for (
int dx = 0; dx < TE_D1D; ++dx)
506 opXY[dy][dx][0] += Bt(dy,qy)*opX[dx][0];
507 opXY[dy][dx][1] += Bt(dy,qy)*opX[dx][1];
508 opXY[dy][dx][2] += Bt(dy,qy)*opX[dx][2];
512 for (
int dz = 0; dz < TE_D1D; ++dz)
514 for (
int dy = 0; dy < TE_D1D; ++dy)
516 for (
int dx = 0; dx < TE_D1D; ++dx)
518 y(dx,dy,dz,0,e) += Bt(dz,qz)*opXY[dy][dx][0];
519 y(dx,dy,dz,1,e) += Bt(dz,qz)*opXY[dy][dx][1];
520 y(dx,dy,dz,2,e) += Bt(dz,qz)*opXY[dy][dx][2];
530 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
531 static void PAGradientApplyTranspose3D(
const int NE,
532 const Array<double> &bt,
533 const Array<double> >,
534 const Array<double> &b,
542 MFEM_ASSERT(
false,
"Gradient PA Apply Transpose 3D not implemented.");
546 template<const
int T_TR_D1D = 0, const
int T_TE_D1D = 0, const
int T_Q1D = 0>
547 static void SmemPAGradientApply3D(
const int NE,
548 const Array<double> &b_,
549 const Array<double> &g_,
550 const Array<double> &bt_,
554 const int tr_d1d = 0,
555 const int te_d1d = 0,
558 const int TR_D1D = T_TR_D1D ? T_TR_D1D : tr_d1d;
559 const int TE_D1D = T_TE_D1D ? T_TE_D1D : te_d1d;
560 const int Q1D = T_Q1D ? T_Q1D : q1d;
562 MFEM_VERIFY(TR_D1D <=
MAX_D1D,
"");
563 MFEM_VERIFY(TE_D1D <=
MAX_D1D,
"");
564 MFEM_VERIFY(TR_D1D <= Q1D,
"");
565 MFEM_VERIFY(TE_D1D <= Q1D,
"");
566 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
568 auto b =
Reshape(b_.Read(), Q1D, TR_D1D);
569 auto g =
Reshape(g_.Read(), Q1D, TR_D1D);
570 auto bt =
Reshape(bt_.Read(), TE_D1D, Q1D);
571 auto D =
Reshape(d_.Read(), Q1D*Q1D*Q1D, 3, 3, NE);
572 auto x =
Reshape(x_.Read(), TR_D1D, TR_D1D, TR_D1D, NE);
573 auto y =
Reshape(y_.ReadWrite(), TE_D1D, TE_D1D, TE_D1D, 3, NE);
575 MFEM_FORALL_3D(e, NE, (Q1D>8)?8:Q1D, (Q1D>8)?8:Q1D, (Q1D>8)?8:Q1D,
577 const int tidz = MFEM_THREAD_ID(z);
578 const int D1DR = T_TR_D1D ? T_TR_D1D : tr_d1d;
579 const int D1DE = T_TE_D1D ? T_TE_D1D : te_d1d;
580 const int Q1D = T_Q1D ? T_Q1D : q1d;
581 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
582 constexpr
int MD1R = T_TR_D1D ? T_TR_D1D :
MAX_D1D;
583 constexpr
int MD1E = T_TE_D1D ? T_TE_D1D :
MAX_D1D;
584 constexpr
int MD1 = MD1E > MD1R ? MD1E : MD1R;
585 constexpr
int MDQ = MQ1 > MD1 ? MQ1 : MD1;
586 MFEM_SHARED
double sBG[2][MQ1*MD1];
587 double (*B)[MD1] = (double (*)[MD1]) (sBG+0);
588 double (*G)[MD1] = (double (*)[MD1]) (sBG+1);
589 double (*Bt)[MQ1] = (double (*)[MQ1]) (sBG+0);
590 MFEM_SHARED
double sm0[3][MDQ*MDQ*MDQ];
591 MFEM_SHARED
double sm1[3][MDQ*MDQ*MDQ];
592 double (*X)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+2);
593 double (*DDQ0)[MD1][MQ1] = (double (*)[MD1][MQ1]) (sm0+0);
594 double (*DDQ1)[MD1][MQ1] = (double (*)[MD1][MQ1]) (sm0+1);
596 double (*DQQ0)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+0);
597 double (*DQQ1)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+1);
598 double (*DQQ2)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm1+2);
600 double (*QQQ0)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+0);
601 double (*QQQ1)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+1);
602 double (*QQQ2)[MQ1][MQ1] = (double (*)[MQ1][MQ1]) (sm0+2);
604 double (*QQD0)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+0);
605 double (*QQD1)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+1);
606 double (*QQD2)[MQ1][MD1] = (double (*)[MQ1][MD1]) (sm1+2);
608 double (*QDD0)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+0);
609 double (*QDD1)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+1);
610 double (*QDD2)[MD1][MD1] = (double (*)[MD1][MD1]) (sm0+2);
611 MFEM_FOREACH_THREAD(dz,z,D1DR)
613 MFEM_FOREACH_THREAD(dy,y,D1DR)
615 MFEM_FOREACH_THREAD(dx,x,D1DR)
617 X[dz][dy][dx] = x(dx,dy,dz,e);
623 MFEM_FOREACH_THREAD(d,y,D1DR)
625 MFEM_FOREACH_THREAD(q,x,Q1D)
633 MFEM_FOREACH_THREAD(dz,z,D1DR)
635 MFEM_FOREACH_THREAD(dy,y,D1DR)
637 MFEM_FOREACH_THREAD(qx,x,Q1D)
641 for (
int dx = 0; dx < D1DR; ++dx)
643 const double coord = X[dz][dy][dx];
644 u += coord * B[qx][dx];
645 v += coord * G[qx][dx];
647 DDQ0[dz][dy][qx] =
u;
648 DDQ1[dz][dy][qx] = v;
653 MFEM_FOREACH_THREAD(dz,z,D1DR)
655 MFEM_FOREACH_THREAD(qy,y,Q1D)
657 MFEM_FOREACH_THREAD(qx,x,Q1D)
662 for (
int dy = 0; dy < D1DR; ++dy)
664 u += DDQ1[dz][dy][qx] * B[qy][dy];
665 v += DDQ0[dz][dy][qx] * G[qy][dy];
666 w += DDQ0[dz][dy][qx] * B[qy][dy];
668 DQQ0[dz][qy][qx] =
u;
669 DQQ1[dz][qy][qx] = v;
670 DQQ2[dz][qy][qx] = w;
675 MFEM_FOREACH_THREAD(qz,z,Q1D)
677 MFEM_FOREACH_THREAD(qy,y,Q1D)
679 MFEM_FOREACH_THREAD(qx,x,Q1D)
684 for (
int dz = 0; dz < D1DR; ++dz)
686 u += DQQ0[dz][qy][qx] * B[qz][dz];
687 v += DQQ1[dz][qy][qx] * B[qz][dz];
688 w += DQQ2[dz][qy][qx] * G[qz][dz];
690 QQQ0[qz][qy][qx] =
u;
691 QQQ1[qz][qy][qx] = v;
692 QQQ2[qz][qy][qx] = w;
697 MFEM_FOREACH_THREAD(qz,z,Q1D)
699 MFEM_FOREACH_THREAD(qy,y,Q1D)
701 MFEM_FOREACH_THREAD(qx,x,Q1D)
703 const int q = qx + (qy + qz * Q1D) * Q1D;
704 const double gX = QQQ0[qz][qy][qx];
705 const double gY = QQQ1[qz][qy][qx];
706 const double gZ = QQQ2[qz][qy][qx];
707 QQQ0[qz][qy][qx] = (D(q,0,0,e)*gX) + (D(q,1,0,e)*gY) + (D(q,2,0,e)*gZ);
708 QQQ1[qz][qy][qx] = (D(q,0,1,e)*gX) + (D(q,1,1,e)*gY) + (D(q,2,1,e)*gZ);
709 QQQ2[qz][qy][qx] = (D(q,0,2,e)*gX) + (D(q,1,2,e)*gY) + (D(q,2,2,e)*gZ);
716 MFEM_FOREACH_THREAD(d,y,D1DE)
718 MFEM_FOREACH_THREAD(q,x,Q1D)
725 MFEM_FOREACH_THREAD(qz,z,Q1D)
727 MFEM_FOREACH_THREAD(qy,y,Q1D)
729 MFEM_FOREACH_THREAD(dx,x,D1DE)
734 for (
int qx = 0; qx < Q1D; ++qx)
736 u += QQQ0[qz][qy][qx] * Bt[dx][qx];
737 v += QQQ1[qz][qy][qx] * Bt[dx][qx];
738 w += QQQ2[qz][qy][qx] * Bt[dx][qx];
740 QQD0[qz][qy][dx] =
u;
741 QQD1[qz][qy][dx] = v;
742 QQD2[qz][qy][dx] = w;
747 MFEM_FOREACH_THREAD(qz,z,Q1D)
749 MFEM_FOREACH_THREAD(dy,y,D1DE)
751 MFEM_FOREACH_THREAD(dx,x,D1DE)
756 for (
int qy = 0; qy < Q1D; ++qy)
758 u += QQD0[qz][qy][dx] * Bt[dy][qy];
759 v += QQD1[qz][qy][dx] * Bt[dy][qy];
760 w += QQD2[qz][qy][dx] * Bt[dy][qy];
762 QDD0[qz][dy][dx] =
u;
763 QDD1[qz][dy][dx] = v;
764 QDD2[qz][dy][dx] = w;
769 MFEM_FOREACH_THREAD(dz,z,D1DE)
771 MFEM_FOREACH_THREAD(dy,y,D1DE)
773 MFEM_FOREACH_THREAD(dx,x,D1DE)
778 for (
int qz = 0; qz < Q1D; ++qz)
780 u += QDD0[qz][dy][dx] * Bt[dz][qz];
781 v += QDD1[qz][dy][dx] * Bt[dz][qz];
782 w += QDD2[qz][dy][dx] * Bt[dz][qz];
784 y(dx,dy,dz,0,e) +=
u;
785 y(dx,dy,dz,1,e) += v;
786 y(dx,dy,dz,2,e) += w;
793 static void PAGradientApply(
const int dim,
798 const Array<double> &B,
799 const Array<double> &G,
800 const Array<double> &Bt,
804 bool transpose=
false)
809 return PAGradientApply2D(NE,B,G,Bt,op,x,y,TR_D1D,TE_D1D,Q1D);
813 return PAGradientApply3D(NE,B,G,Bt,op,x,y,TR_D1D,TE_D1D,Q1D);
815 MFEM_ABORT(
"Unknown kernel.");
821 PAGradientApply(dim, trial_dofs1D, test_dofs1D, quad1D, ne,
822 trial_maps->B, trial_maps->G, test_maps->Bt, pa_data, x, y,
827 void GradientIntegrator::AddMultTransposePA(
const Vector &x,
Vector &y)
const
829 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.
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.
Class to represent a coefficient evaluated at quadrature points.
int GetNE() const
Returns number of elements in the mesh.
Mesh * GetMesh() const
Returns the mesh.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
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...
const IntegrationRule & GetRule(const Integrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
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...
Class representing the storage layout of a QuadratureFunction.
double u(const Vector &xvec)
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.