20void PAHcurlMassAssembleDiagonal2D(
const int D1D,
24 const Array<real_t> &bo,
25 const Array<real_t> &bc,
26 const Vector &pa_data,
29 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
30 auto Bc =
Reshape(bc.Read(), Q1D, D1D);
31 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, symmetric ? 3 : 4, NE);
32 auto D =
Reshape(diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
36 constexpr static int VDIM = 2;
37 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
41 for (
int c = 0; c < VDIM; ++c)
43 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
44 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
48 for (
int dy = 0; dy < D1Dy; ++dy)
50 for (
int qx = 0; qx < Q1D; ++qx)
53 for (
int qy = 0; qy < Q1D; ++qy)
55 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
57 mass[qx] += wy * wy * ((c == 0) ? op(qx,qy,0,e) :
58 op(qx,qy,symmetric ? 2 : 3, e));
62 for (
int dx = 0; dx < D1Dx; ++dx)
64 for (
int qx = 0; qx < Q1D; ++qx)
66 const real_t wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
67 D(dx + (dy * D1Dx) + osc, e) += mass[qx] * wx * wx;
77void PAHcurlMassAssembleDiagonal3D(
const int D1D,
81 const Array<real_t> &bo,
82 const Array<real_t> &bc,
83 const Vector &pa_data,
87 "Error: D1D > MAX_D1D");
89 "Error: Q1D > MAX_Q1D");
90 constexpr static int VDIM = 3;
92 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
93 auto Bc =
Reshape(bc.Read(), Q1D, D1D);
94 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, Q1D, symmetric ? 6 : 9, NE);
95 auto D =
Reshape(diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
99 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
103 for (
int c = 0; c < VDIM; ++c)
105 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
106 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
107 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
109 const int opc = (c == 0) ? 0 : ((c == 1) ? (symmetric ? 3 : 4) :
110 (symmetric ? 5 : 8));
114 for (
int dz = 0; dz < D1Dz; ++dz)
116 for (
int dy = 0; dy < D1Dy; ++dy)
118 for (
int qx = 0; qx < Q1D; ++qx)
121 for (
int qy = 0; qy < Q1D; ++qy)
123 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
125 for (
int qz = 0; qz < Q1D; ++qz)
127 const real_t wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
129 mass[qx] += wy * wy * wz * wz * op(qx,qy,qz,opc,e);
134 for (
int dx = 0; dx < D1Dx; ++dx)
136 for (
int qx = 0; qx < Q1D; ++qx)
138 const real_t wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
139 D(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += mass[qx] * wx * wx;
145 osc += D1Dx * D1Dy * D1Dz;
150void PAHcurlMassApply2D(
const int D1D,
153 const bool symmetric,
154 const Array<real_t> &bo,
155 const Array<real_t> &bc,
156 const Array<real_t> &bot,
157 const Array<real_t> &bct,
158 const Vector &pa_data,
162 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
163 auto Bc =
Reshape(bc.Read(), Q1D, D1D);
164 auto Bot =
Reshape(bot.Read(), D1D-1, Q1D);
165 auto Bct =
Reshape(bct.Read(), D1D, Q1D);
166 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, symmetric ? 3 : 4, NE);
167 auto X =
Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
168 auto Y =
Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
172 constexpr static int VDIM = 2;
173 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
174 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
176 real_t mass[MAX_Q1D][MAX_Q1D][VDIM];
178 for (
int qy = 0; qy < Q1D; ++qy)
180 for (
int qx = 0; qx < Q1D; ++qx)
182 for (
int c = 0; c < VDIM; ++c)
184 mass[qy][qx][c] = 0.0;
191 for (
int c = 0; c < VDIM; ++c)
193 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
194 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
196 for (
int dy = 0; dy < D1Dy; ++dy)
199 for (
int qx = 0; qx < Q1D; ++qx)
204 for (
int dx = 0; dx < D1Dx; ++dx)
206 const real_t t = X(dx + (dy * D1Dx) + osc, e);
207 for (
int qx = 0; qx < Q1D; ++qx)
209 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
213 for (
int qy = 0; qy < Q1D; ++qy)
215 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
216 for (
int qx = 0; qx < Q1D; ++qx)
218 mass[qy][qx][c] += massX[qx] * wy;
227 for (
int qy = 0; qy < Q1D; ++qy)
229 for (
int qx = 0; qx < Q1D; ++qx)
231 const real_t O11 = op(qx,qy,0,e);
232 const real_t O21 = op(qx,qy,1,e);
233 const real_t O12 = symmetric ? O21 : op(qx,qy,2,e);
234 const real_t O22 = symmetric ? op(qx,qy,2,e) : op(qx,qy,3,e);
235 const real_t massX = mass[qy][qx][0];
236 const real_t massY = mass[qy][qx][1];
237 mass[qy][qx][0] = (O11*massX)+(O12*massY);
238 mass[qy][qx][1] = (O21*massX)+(O22*massY);
242 for (
int qy = 0; qy < Q1D; ++qy)
246 for (
int c = 0; c < VDIM; ++c)
248 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
249 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
252 for (
int dx = 0; dx < D1Dx; ++dx)
256 for (
int qx = 0; qx < Q1D; ++qx)
258 for (
int dx = 0; dx < D1Dx; ++dx)
260 massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
264 for (
int dy = 0; dy < D1Dy; ++dy)
266 const real_t wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
268 for (
int dx = 0; dx < D1Dx; ++dx)
270 Y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
280void PAHcurlMassApply3D(
const int D1D,
283 const bool symmetric,
284 const Array<real_t> &bo,
285 const Array<real_t> &bc,
286 const Array<real_t> &bot,
287 const Array<real_t> &bct,
288 const Vector &pa_data,
293 "Error: D1D > MAX_D1D");
295 "Error: Q1D > MAX_Q1D");
296 constexpr static int VDIM = 3;
298 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
299 auto Bc =
Reshape(bc.Read(), Q1D, D1D);
300 auto Bot =
Reshape(bot.Read(), D1D-1, Q1D);
301 auto Bct =
Reshape(bct.Read(), D1D, Q1D);
302 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, Q1D, symmetric ? 6 : 9, NE);
303 auto X =
Reshape(x.Read(), 3*(D1D-1)*D1D*D1D, NE);
304 auto Y =
Reshape(y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
308 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
309 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
311 real_t mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
313 for (
int qz = 0; qz < Q1D; ++qz)
315 for (
int qy = 0; qy < Q1D; ++qy)
317 for (
int qx = 0; qx < Q1D; ++qx)
319 for (
int c = 0; c < VDIM; ++c)
321 mass[qz][qy][qx][c] = 0.0;
329 for (
int c = 0; c < VDIM; ++c)
331 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
332 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
333 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
335 for (
int dz = 0; dz < D1Dz; ++dz)
337 real_t massXY[MAX_Q1D][MAX_Q1D];
338 for (
int qy = 0; qy < Q1D; ++qy)
340 for (
int qx = 0; qx < Q1D; ++qx)
342 massXY[qy][qx] = 0.0;
346 for (
int dy = 0; dy < D1Dy; ++dy)
349 for (
int qx = 0; qx < Q1D; ++qx)
354 for (
int dx = 0; dx < D1Dx; ++dx)
356 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
357 for (
int qx = 0; qx < Q1D; ++qx)
359 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
363 for (
int qy = 0; qy < Q1D; ++qy)
365 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
366 for (
int qx = 0; qx < Q1D; ++qx)
368 const real_t wx = massX[qx];
369 massXY[qy][qx] += wx * wy;
374 for (
int qz = 0; qz < Q1D; ++qz)
376 const real_t wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
377 for (
int qy = 0; qy < Q1D; ++qy)
379 for (
int qx = 0; qx < Q1D; ++qx)
381 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
387 osc += D1Dx * D1Dy * D1Dz;
391 for (
int qz = 0; qz < Q1D; ++qz)
393 for (
int qy = 0; qy < Q1D; ++qy)
395 for (
int qx = 0; qx < Q1D; ++qx)
397 const real_t O11 = op(qx,qy,qz,0,e);
398 const real_t O12 = op(qx,qy,qz,1,e);
399 const real_t O13 = op(qx,qy,qz,2,e);
400 const real_t O21 = symmetric ? O12 : op(qx,qy,qz,3,e);
401 const real_t O22 = symmetric ? op(qx,qy,qz,3,e) : op(qx,qy,qz,4,e);
402 const real_t O23 = symmetric ? op(qx,qy,qz,4,e) : op(qx,qy,qz,5,e);
403 const real_t O31 = symmetric ? O13 : op(qx,qy,qz,6,e);
404 const real_t O32 = symmetric ? O23 : op(qx,qy,qz,7,e);
405 const real_t O33 = symmetric ? op(qx,qy,qz,5,e) : op(qx,qy,qz,8,e);
406 const real_t massX = mass[qz][qy][qx][0];
407 const real_t massY = mass[qz][qy][qx][1];
408 const real_t massZ = mass[qz][qy][qx][2];
409 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
410 mass[qz][qy][qx][1] = (O21*massX)+(O22*massY)+(O23*massZ);
411 mass[qz][qy][qx][2] = (O31*massX)+(O32*massY)+(O33*massZ);
416 for (
int qz = 0; qz < Q1D; ++qz)
418 real_t massXY[MAX_D1D][MAX_D1D];
422 for (
int c = 0; c < VDIM; ++c)
424 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
425 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
426 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
428 for (
int dy = 0; dy < D1Dy; ++dy)
430 for (
int dx = 0; dx < D1Dx; ++dx)
432 massXY[dy][dx] = 0.0;
435 for (
int qy = 0; qy < Q1D; ++qy)
438 for (
int dx = 0; dx < D1Dx; ++dx)
442 for (
int qx = 0; qx < Q1D; ++qx)
444 for (
int dx = 0; dx < D1Dx; ++dx)
446 massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
449 for (
int dy = 0; dy < D1Dy; ++dy)
451 const real_t wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
452 for (
int dx = 0; dx < D1Dx; ++dx)
454 massXY[dy][dx] += massX[dx] * wy;
459 for (
int dz = 0; dz < D1Dz; ++dz)
461 const real_t wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
462 for (
int dy = 0; dy < D1Dy; ++dy)
464 for (
int dx = 0; dx < D1Dx; ++dx)
466 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
471 osc += D1Dx * D1Dy * D1Dz;
477void PACurlCurlSetup2D(
const int Q1D,
479 const Array<real_t> &w,
484 const int NQ = Q1D*Q1D;
486 auto J =
Reshape(j.Read(), NQ, 2, 2, NE);
487 auto C =
Reshape(coeff.Read(), NQ, NE);
488 auto y =
Reshape(op.Write(), NQ, NE);
491 for (
int q = 0; q < NQ; ++q)
493 const real_t J11 = J(q,0,0,e);
494 const real_t J21 = J(q,1,0,e);
495 const real_t J12 = J(q,0,1,e);
496 const real_t J22 = J(q,1,1,e);
497 const real_t detJ = (J11*J22)-(J21*J12);
498 y(q,e) = W[q] * C(q,e) / detJ;
503void PACurlCurlSetup3D(
const int Q1D,
506 const Array<real_t> &w,
511 const int NQ = Q1D*Q1D*Q1D;
512 const bool symmetric = (coeffDim != 9);
514 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
515 auto C =
Reshape(coeff.Read(), coeffDim, NQ, NE);
516 auto y =
Reshape(op.Write(), NQ, symmetric ? 6 : 9, NE);
520 for (
int q = 0; q < NQ; ++q)
522 const real_t J11 = J(q,0,0,e);
523 const real_t J21 = J(q,1,0,e);
524 const real_t J31 = J(q,2,0,e);
525 const real_t J12 = J(q,0,1,e);
526 const real_t J22 = J(q,1,1,e);
527 const real_t J32 = J(q,2,1,e);
528 const real_t J13 = J(q,0,2,e);
529 const real_t J23 = J(q,1,2,e);
530 const real_t J33 = J(q,2,2,e);
531 const real_t detJ = J11 * (J22 * J33 - J32 * J23) -
532 J21 * (J12 * J33 - J32 * J13) +
533 J31 * (J12 * J23 - J22 * J13);
535 const real_t c_detJ = W[q] / detJ;
537 if (coeffDim == 6 || coeffDim == 9)
540 const real_t M11 = C(0, q, e);
541 const real_t M12 = C(1, q, e);
542 const real_t M13 = C(2, q, e);
543 const real_t M21 = (!symmetric) ? C(3, q, e) : M12;
544 const real_t M22 = (!symmetric) ? C(4, q, e) : C(3, q, e);
545 const real_t M23 = (!symmetric) ? C(5, q, e) : C(4, q, e);
546 const real_t M31 = (!symmetric) ? C(6, q, e) : M13;
547 const real_t M32 = (!symmetric) ? C(7, q, e) : M23;
548 const real_t M33 = (!symmetric) ? C(8, q, e) : C(5, q, e);
551 const real_t R11 = M11*J11 + M12*J21 + M13*J31;
552 const real_t R12 = M11*J12 + M12*J22 + M13*J32;
553 const real_t R13 = M11*J13 + M12*J23 + M13*J33;
554 const real_t R21 = M21*J11 + M22*J21 + M23*J31;
555 const real_t R22 = M21*J12 + M22*J22 + M23*J32;
556 const real_t R23 = M21*J13 + M22*J23 + M23*J33;
557 const real_t R31 = M31*J11 + M32*J21 + M33*J31;
558 const real_t R32 = M31*J12 + M32*J22 + M33*J32;
559 const real_t R33 = M31*J13 + M32*J23 + M33*J33;
562 y(q,0,e) = c_detJ * (J11*R11 + J21*R21 + J31*R31);
563 const real_t Y12 = c_detJ * (J11*R12 + J21*R22 + J31*R32);
565 y(q,2,e) = c_detJ * (J11*R13 + J21*R23 + J31*R33);
567 const real_t Y21 = c_detJ * (J12*R11 + J22*R21 + J32*R31);
568 const real_t Y22 = c_detJ * (J12*R12 + J22*R22 + J32*R32);
569 const real_t Y23 = c_detJ * (J12*R13 + J22*R23 + J32*R33);
571 const real_t Y33 = c_detJ * (J13*R13 + J23*R23 + J33*R33);
573 y(q,3,e) = symmetric ? Y22 : Y21;
574 y(q,4,e) = symmetric ? Y23 : Y22;
575 y(q,5,e) = symmetric ? Y33 : Y23;
579 y(q,6,e) = c_detJ * (J13*R11 + J23*R21 + J33*R31);
580 y(q,7,e) = c_detJ * (J13*R12 + J23*R22 + J33*R32);
587 const real_t D1 = C(0, q, e);
588 const real_t D2 = coeffDim == 3 ? C(1, q, e) : D1;
589 const real_t D3 = coeffDim == 3 ? C(2, q, e) : D1;
591 y(q,0,e) = c_detJ * (D1*J11*J11 + D2*J21*J21 + D3*J31*J31);
592 y(q,1,e) = c_detJ * (D1*J11*J12 + D2*J21*J22 + D3*J31*J32);
593 y(q,2,e) = c_detJ * (D1*J11*J13 + D2*J21*J23 + D3*J31*J33);
594 y(q,3,e) = c_detJ * (D1*J12*J12 + D2*J22*J22 + D3*J32*J32);
595 y(q,4,e) = c_detJ * (D1*J12*J13 + D2*J22*J23 + D3*J32*J33);
596 y(q,5,e) = c_detJ * (D1*J13*J13 + D2*J23*J23 + D3*J33*J33);
602void PACurlCurlAssembleDiagonal2D(
const int D1D,
605 const Array<real_t> &bo,
606 const Array<real_t> &gc,
607 const Vector &pa_data,
610 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
611 auto Gc =
Reshape(gc.Read(), Q1D, D1D);
612 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, NE);
613 auto D =
Reshape(diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
617 constexpr static int VDIM = 2;
618 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
622 for (
int c = 0; c < VDIM; ++c)
624 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
625 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
629 for (
int dy = 0; dy < D1Dy; ++dy)
631 for (
int qx = 0; qx < Q1D; ++qx)
634 for (
int qy = 0; qy < Q1D; ++qy)
636 const real_t wy = (c == 1) ? Bo(qy,dy) : -Gc(qy,dy);
637 t[qx] += wy * wy * op(qx,qy,e);
641 for (
int dx = 0; dx < D1Dx; ++dx)
643 for (
int qx = 0; qx < Q1D; ++qx)
645 const real_t wx = ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
646 D(dx + (dy * D1Dx) + osc, e) += t[qx] * wx * wx;
656void PACurlCurlApply2D(
const int D1D,
659 const Array<real_t> &bo,
660 const Array<real_t> &bot,
661 const Array<real_t> &gc,
662 const Array<real_t> &gct,
663 const Vector &pa_data,
668 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
669 auto Bot =
Reshape(bot.Read(), D1D-1, Q1D);
670 auto Gc =
Reshape(gc.Read(), Q1D, D1D);
671 auto Gct =
Reshape(gct.Read(), D1D, Q1D);
672 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, NE);
673 auto X =
Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
674 auto Y =
Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
678 constexpr static int VDIM = 2;
679 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
680 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
682 real_t curl[MAX_Q1D][MAX_Q1D];
686 for (
int qy = 0; qy < Q1D; ++qy)
688 for (
int qx = 0; qx < Q1D; ++qx)
696 for (
int c = 0; c < VDIM; ++c)
698 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
699 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
701 for (
int dy = 0; dy < D1Dy; ++dy)
704 for (
int qx = 0; qx < Q1D; ++qx)
709 for (
int dx = 0; dx < D1Dx; ++dx)
711 const real_t t = X(dx + (dy * D1Dx) + osc, e);
712 for (
int qx = 0; qx < Q1D; ++qx)
714 gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
718 for (
int qy = 0; qy < Q1D; ++qy)
720 const real_t wy = (c == 0) ? -Gc(qy,dy) : Bo(qy,dy);
721 for (
int qx = 0; qx < Q1D; ++qx)
723 curl[qy][qx] += gradX[qx] * wy;
732 for (
int qy = 0; qy < Q1D; ++qy)
734 for (
int qx = 0; qx < Q1D; ++qx)
736 curl[qy][qx] *= op(qx,qy,e);
740 for (
int qy = 0; qy < Q1D; ++qy)
744 for (
int c = 0; c < VDIM; ++c)
746 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
747 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
750 for (
int dx = 0; dx < D1Dx; ++dx)
754 for (
int qx = 0; qx < Q1D; ++qx)
756 for (
int dx = 0; dx < D1Dx; ++dx)
758 gradX[dx] += curl[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
761 for (
int dy = 0; dy < D1Dy; ++dy)
763 const real_t wy = (c == 0) ? -Gct(dy,qy) : Bot(dy,qy);
765 for (
int dx = 0; dx < D1Dx; ++dx)
767 Y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
777void PAHcurlL2Setup2D(
const int Q1D,
779 const Array<real_t> &w,
783 const int NQ = Q1D*Q1D;
785 auto C =
Reshape(coeff.Read(), NQ, NE);
786 auto y =
Reshape(op.Write(), NQ, NE);
789 for (
int q = 0; q < NQ; ++q)
791 y(q,e) = W[q] * C(q,e);
796void PAHcurlL2Setup3D(
const int NQ,
799 const Array<real_t> &w,
804 auto C =
Reshape(coeff.Read(), coeffDim, NQ, NE);
805 auto y =
Reshape(op.Write(), coeffDim, NQ, NE);
809 for (
int q = 0; q < NQ; ++q)
811 for (
int c=0; c<coeffDim; ++c)
813 y(c,q,e) = W[q] * C(c,q,e);
819void PAHcurlL2Apply2D(
const int D1D,
823 const Array<real_t> &bo,
824 const Array<real_t> &bot,
825 const Array<real_t> &bt,
826 const Array<real_t> &gc,
827 const Vector &pa_data,
831 const int H1 = (D1Dtest == D1D);
833 MFEM_VERIFY(y.Size() == NE*D1Dtest*D1Dtest,
"Test vector of wrong dimension");
835 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
836 auto Bot =
Reshape(bot.Read(), D1D-1, Q1D);
837 auto Bt =
Reshape(bt.Read(), D1D, Q1D);
838 auto Gc =
Reshape(gc.Read(), Q1D, D1D);
839 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, NE);
840 auto X =
Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
841 auto Y =
Reshape(y.ReadWrite(), D1Dtest, D1Dtest, NE);
845 constexpr static int VDIM = 2;
846 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
847 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
849 real_t curl[MAX_Q1D][MAX_Q1D];
853 for (
int qy = 0; qy < Q1D; ++qy)
855 for (
int qx = 0; qx < Q1D; ++qx)
863 for (
int c = 0; c < VDIM; ++c)
865 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
866 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
868 for (
int dy = 0; dy < D1Dy; ++dy)
871 for (
int qx = 0; qx < Q1D; ++qx)
876 for (
int dx = 0; dx < D1Dx; ++dx)
878 const real_t t = X(dx + (dy * D1Dx) + osc, e);
879 for (
int qx = 0; qx < Q1D; ++qx)
881 gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
885 for (
int qy = 0; qy < Q1D; ++qy)
887 const real_t wy = (c == 0) ? -Gc(qy,dy) : Bo(qy,dy);
888 for (
int qx = 0; qx < Q1D; ++qx)
890 curl[qy][qx] += gradX[qx] * wy;
899 for (
int qy = 0; qy < Q1D; ++qy)
901 for (
int qx = 0; qx < Q1D; ++qx)
903 curl[qy][qx] *= op(qx,qy,e);
907 for (
int qy = 0; qy < Q1D; ++qy)
910 for (
int dx = 0; dx < D1Dtest; ++dx)
914 for (
int qx = 0; qx < Q1D; ++qx)
916 const real_t s = curl[qy][qx];
917 for (
int dx = 0; dx < D1Dtest; ++dx)
919 sol_x[dx] += s * ((H1 == 1) ? Bt(dx,qx) : Bot(dx,qx));
922 for (
int dy = 0; dy < D1Dtest; ++dy)
924 const real_t wy = (H1 == 1) ? Bt(dy,qy) : Bot(dy,qy);
926 for (
int dx = 0; dx < D1Dtest; ++dx)
928 Y(dx,dy,e) += sol_x[dx] * wy;
935void PAHcurlL2ApplyTranspose2D(
const int D1D,
939 const Array<real_t> &bo,
940 const Array<real_t> &bot,
941 const Array<real_t> &
b,
942 const Array<real_t> &gct,
943 const Vector &pa_data,
947 const int H1 = (D1Dtest == D1D);
949 MFEM_VERIFY(x.Size() == NE*D1Dtest*D1Dtest,
"Test vector of wrong dimension");
951 auto Bo =
Reshape(bo.Read(), Q1D, D1D-1);
952 auto B =
Reshape(
b.Read(), Q1D, D1D);
953 auto Bot =
Reshape(bot.Read(), D1D-1, Q1D);
954 auto Gct =
Reshape(gct.Read(), D1D, Q1D);
955 auto op =
Reshape(pa_data.Read(), Q1D, Q1D, NE);
956 auto X =
Reshape(x.Read(), D1Dtest, D1Dtest, NE);
957 auto Y =
Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
961 constexpr static int VDIM = 2;
962 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
963 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
965 real_t mass[MAX_Q1D][MAX_Q1D];
969 for (
int qy = 0; qy < Q1D; ++qy)
971 for (
int qx = 0; qx < Q1D; ++qx)
977 for (
int dy = 0; dy < D1Dtest; ++dy)
980 for (
int qy = 0; qy < Q1D; ++qy)
984 for (
int dx = 0; dx < D1Dtest; ++dx)
986 const real_t s = X(dx,dy,e);
987 for (
int qx = 0; qx < Q1D; ++qx)
989 sol_x[qx] += s * ((H1 == 1) ? B(qx,dx) : Bo(qx,dx));
992 for (
int qy = 0; qy < Q1D; ++qy)
994 const real_t d2q = (H1 == 1) ? B(qy,dy) : Bo(qy,dy);
995 for (
int qx = 0; qx < Q1D; ++qx)
997 mass[qy][qx] += d2q * sol_x[qx];
1003 for (
int qy = 0; qy < Q1D; ++qy)
1005 for (
int qx = 0; qx < Q1D; ++qx)
1007 mass[qy][qx] *= op(qx,qy,e);
1011 for (
int qy = 0; qy < Q1D; ++qy)
1015 for (
int c = 0; c < VDIM; ++c)
1017 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1018 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1021 for (
int dx = 0; dx < D1Dx; ++dx)
1025 for (
int qx = 0; qx < Q1D; ++qx)
1027 for (
int dx = 0; dx < D1Dx; ++dx)
1029 gradX[dx] += mass[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
1032 for (
int dy = 0; dy < D1Dy; ++dy)
1034 const real_t wy = (c == 0) ? -Gct(dy,qy) : Bot(dy,qy);
1036 for (
int dx = 0; dx < D1Dx; ++dx)
1038 Y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void forall(int N, lambda &&body)
static const DeviceDofQuadLimits & Get()
Return a const reference to the DeviceDofQuadLimits singleton.