12 #include "../general/forall.hpp"
21 const Array<double> &w,
30 const Array<double> &bo,
31 const Array<double> &bc,
32 const Vector &pa_data,
39 const Array<double> &bo,
40 const Array<double> &bc,
41 const Vector &pa_data,
44 template<
int T_D1D = 0,
int T_Q1D = 0>
49 const Array<double> &bo,
50 const Array<double> &bc,
51 const Vector &pa_data,
58 const Array<double> &bo,
59 const Array<double> &bc,
60 const Array<double> &bot,
61 const Array<double> &bct,
62 const Vector &pa_data,
70 const Array<double> &bo,
71 const Array<double> &bc,
72 const Array<double> &bot,
73 const Array<double> &bct,
74 const Vector &pa_data,
78 template<
int T_D1D = 0,
int T_Q1D = 0>
83 const Array<double> &bo,
84 const Array<double> &bc,
85 const Array<double> &bot,
86 const Array<double> &bct,
87 const Vector &pa_data,
93 const Array<double> &w,
100 const Array<double> &w,
108 const Array<double> &bc,
109 const Array<double> &gc,
110 const Array<double> &bot,
111 const Array<double> &bct,
112 const Vector &pa_data,
119 const Array<double> &bc,
120 const Array<double> &gc,
121 const Array<double> &bot,
122 const Array<double> &bct,
123 const Vector &pa_data,
130 const Array<double> &_Bo,
131 const Array<double> &_Bc,
138 const Array<double> &_Bo,
139 const Array<double> &_Bc,
146 const Array<double> &_Bo,
147 const Array<double> &_Bc,
148 const Array<double> &_Bot,
149 const Array<double> &_Bct,
157 const Array<double> &_Bo,
158 const Array<double> &_Bc,
159 const Array<double> &_Bot,
160 const Array<double> &_Bct,
168 const Array<double> &w,
178 const bool transpose,
184 const bool symmetric = (coeffDim != 9);
186 auto J =
Reshape(j.
Read(), Q1D, Q1D, Q1D, 3, 3, NE);
187 auto coeff =
Reshape(_coeff.
Read(), coeffDim, Q1D, Q1D, Q1D, NE);
191 const int i12 = transpose ? 3 : 1;
192 const int i13 = transpose ? 6 : 2;
193 const int i21 = transpose ? 1 : 3;
195 const int i23 = transpose ? 7 : 5;
196 const int i31 = transpose ? 2 : 6;
197 const int i32 = transpose ? 5 : 7;
200 MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
202 MFEM_FOREACH_THREAD(qx,x,Q1D)
204 MFEM_FOREACH_THREAD(qy,y,Q1D)
206 MFEM_FOREACH_THREAD(qz,z,Q1D)
208 const double J11 = J(qx,qy,qz,0,0,e);
209 const double J21 = J(qx,qy,qz,1,0,e);
210 const double J31 = J(qx,qy,qz,2,0,e);
211 const double J12 = J(qx,qy,qz,0,1,e);
212 const double J22 = J(qx,qy,qz,1,1,e);
213 const double J32 = J(qx,qy,qz,2,1,e);
214 const double J13 = J(qx,qy,qz,0,2,e);
215 const double J23 = J(qx,qy,qz,1,2,e);
216 const double J33 = J(qx,qy,qz,2,2,e);
217 const double detJ = J11 * (J22 * J33 - J32 * J23) -
218 J21 * (J12 * J33 - J32 * J13) +
219 J31 * (J12 * J23 - J22 * J13);
220 const double w_detJ = W(qx,qy,qz) / detJ;
222 const double A11 = (J22 * J33) - (J23 * J32);
223 const double A12 = (J32 * J13) - (J12 * J33);
224 const double A13 = (J12 * J23) - (J22 * J13);
225 const double A21 = (J31 * J23) - (J21 * J33);
226 const double A22 = (J11 * J33) - (J13 * J31);
227 const double A23 = (J21 * J13) - (J11 * J23);
228 const double A31 = (J21 * J32) - (J31 * J22);
229 const double A32 = (J31 * J12) - (J11 * J32);
230 const double A33 = (J11 * J22) - (J12 * J21);
232 if (coeffDim == 6 || coeffDim == 9)
235 const double M11 = (!symmetric) ? coeff(i11,qx,qy,qz,e) : coeff(0,qx,qy,qz,e);
236 const double M12 = (!symmetric) ? coeff(i12,qx,qy,qz,e) : coeff(1,qx,qy,qz,e);
237 const double M13 = (!symmetric) ? coeff(i13,qx,qy,qz,e) : coeff(2,qx,qy,qz,e);
238 const double M21 = (!symmetric) ? coeff(i21,qx,qy,qz,e) : M12;
239 const double M22 = (!symmetric) ? coeff(i22,qx,qy,qz,e) : coeff(3,qx,qy,qz,e);
240 const double M23 = (!symmetric) ? coeff(i23,qx,qy,qz,e) : coeff(4,qx,qy,qz,e);
241 const double M31 = (!symmetric) ? coeff(i31,qx,qy,qz,e) : M13;
242 const double M32 = (!symmetric) ? coeff(i32,qx,qy,qz,e) : M23;
243 const double M33 = (!symmetric) ? coeff(i33,qx,qy,qz,e) : coeff(5,qx,qy,qz,e);
245 const double R11 = M11*J11 + M12*J12 + M13*J13;
246 const double R12 = M11*J21 + M12*J22 + M13*J23;
247 const double R13 = M11*J31 + M12*J32 + M13*J33;
248 const double R21 = M21*J11 + M22*J12 + M23*J13;
249 const double R22 = M21*J21 + M22*J22 + M23*J23;
250 const double R23 = M21*J31 + M22*J32 + M23*J33;
251 const double R31 = M31*J11 + M32*J12 + M33*J13;
252 const double R32 = M31*J21 + M32*J22 + M33*J23;
253 const double R33 = M31*J31 + M32*J32 + M33*J33;
256 y(i11,qx,qy,qz,e) = w_detJ * (A11*R11 + A12*R21 + A13*R31);
257 y(i12,qx,qy,qz,e) = w_detJ * (A11*R12 + A12*R22 + A13*R32);
258 y(i13,qx,qy,qz,e) = w_detJ * (A11*R13 + A12*R23 + A13*R33);
259 y(i21,qx,qy,qz,e) = w_detJ * (A21*R11 + A22*R21 + A23*R31);
260 y(i22,qx,qy,qz,e) = w_detJ * (A21*R12 + A22*R22 + A23*R32);
261 y(i23,qx,qy,qz,e) = w_detJ * (A21*R13 + A22*R23 + A23*R33);
262 y(i31,qx,qy,qz,e) = w_detJ * (A31*R11 + A32*R21 + A33*R31);
263 y(i32,qx,qy,qz,e) = w_detJ * (A31*R12 + A32*R22 + A33*R32);
264 y(i33,qx,qy,qz,e) = w_detJ * (A31*R13 + A32*R23 + A33*R33);
266 else if (coeffDim == 3)
268 const double D1 = coeff(0,qx,qy,qz,e);
269 const double D2 = coeff(1,qx,qy,qz,e);
270 const double D3 = coeff(2,qx,qy,qz,e);
272 y(i11,qx,qy,qz,e) = w_detJ * (D1*A11*J11 + D2*A12*J21 + D3*A13*J31);
273 y(i12,qx,qy,qz,e) = w_detJ * (D1*A11*J12 + D2*A12*J22 + D3*A13*J32);
274 y(i13,qx,qy,qz,e) = w_detJ * (D1*A11*J13 + D2*A12*J23 + D3*A13*J33);
275 y(i21,qx,qy,qz,e) = w_detJ * (D1*A21*J11 + D2*A22*J21 + D3*A23*J31);
276 y(i22,qx,qy,qz,e) = w_detJ * (D1*A21*J12 + D2*A22*J22 + D3*A23*J32);
277 y(i23,qx,qy,qz,e) = w_detJ * (D1*A21*J13 + D2*A22*J23 + D3*A23*J33);
278 y(i31,qx,qy,qz,e) = w_detJ * (D1*A31*J11 + D2*A32*J21 + D3*A33*J31);
279 y(i32,qx,qy,qz,e) = w_detJ * (D1*A31*J12 + D2*A32*J22 + D3*A33*J32);
280 y(i33,qx,qy,qz,e) = w_detJ * (D1*A31*J13 + D2*A32*J23 + D3*A33*J33);
294 const bool transpose,
300 const bool symmetric = (coeffDim != 4);
303 auto coeff =
Reshape(_coeff.
Read(), coeffDim, Q1D, Q1D, NE);
307 const int i12 = transpose ? 2 : 1;
308 const int i21 = transpose ? 1 : 2;
311 MFEM_FORALL_2D(e, NE, Q1D, Q1D, 1,
313 MFEM_FOREACH_THREAD(qx,x,Q1D)
315 MFEM_FOREACH_THREAD(qy,y,Q1D)
317 const double J11 = J(qx,qy,0,0,e);
318 const double J21 = J(qx,qy,1,0,e);
319 const double J12 = J(qx,qy,0,1,e);
320 const double J22 = J(qx,qy,1,1,e);
321 const double w_detJ = W(qx,qy) / (J11*J22) - (J21*J12);
323 if (coeffDim == 3 || coeffDim == 4)
326 const double M11 = coeff(i11,qx,qy,e);
327 const double M12 = (!symmetric) ? coeff(i12,qx,qy,e) : coeff(1,qx,qy,e);
328 const double M21 = (!symmetric) ? coeff(i21,qx,qy,e) : M12;
329 const double M22 = (!symmetric) ? coeff(i22,qx,qy,e) : coeff(2,qx,qy,e);
331 const double R11 = M11*J11 + M12*J21;
332 const double R12 = M11*J12 + M12*J22;
333 const double R21 = M21*J11 + M22*J21;
334 const double R22 = M21*J12 + M22*J22;
337 y(i11,qx,qy,e) = w_detJ * ( J22*R11 - J12*R21);
338 y(i12,qx,qy,e) = w_detJ * ( J22*R12 - J12*R22);
339 y(i21,qx,qy,e) = w_detJ * (-J21*R11 + J11*R21);
340 y(i22,qx,qy,e) = w_detJ * (-J21*R12 + J11*R22);
342 else if (coeffDim == 2)
344 const double D1 = coeff(0,qx,qy,e);
345 const double D2 = coeff(1,qx,qy,e);
346 const double R11 = D1*J11;
347 const double R12 = D1*J12;
348 const double R21 = D2*J21;
349 const double R22 = D2*J22;
350 y(i11,qx,qy,e) = w_detJ * ( J22*R11 - J12*R21);
351 y(i12,qx,qy,e) = w_detJ * ( J22*R12 - J12*R22);
352 y(i21,qx,qy,e) = w_detJ * (-J21*R11 + J11*R21);
353 y(i22,qx,qy,e) = w_detJ * (-J21*R12 + J11*R22);
366 const bool scalarCoeff,
367 const bool trialHcurl,
379 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D
");
380 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D
");
381 constexpr static int VDIM = 3;
383 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
384 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
385 auto Bot = Reshape(_Bot.Read(), D1Dtest-1, Q1D);
386 auto Bct = Reshape(_Bct.Read(), D1Dtest, Q1D);
387 auto op = Reshape(_op.Read(), scalarCoeff ? 1 : 9, Q1D, Q1D, Q1D, NE);
388 auto x = Reshape(_x.Read(), 3*(D1D-1)*D1D*(trialHcurl ? D1D : D1D-1), NE);
389 auto y = Reshape(_y.ReadWrite(), 3*(D1Dtest-1)*D1Dtest*
390 (trialHcurl ? D1Dtest-1 : D1Dtest), NE);
394 double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
396 for (int qz = 0; qz < Q1D; ++qz)
398 for (int qy = 0; qy < Q1D; ++qy)
400 for (int qx = 0; qx < Q1D; ++qx)
402 for (int c = 0; c < VDIM; ++c)
404 mass[qz][qy][qx][c] = 0.0;
411 for (int c = 0; c < VDIM; ++c) // loop over x, y, z trial components
413 const int D1Dz = trialHcurl ? ((c == 2) ? D1D - 1 : D1D) :
414 ((c == 2) ? D1D : D1D - 1);
415 const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
416 ((c == 1) ? D1D : D1D - 1);
417 const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
418 ((c == 0) ? D1D : D1D - 1);
420 for (int dz = 0; dz < D1Dz; ++dz)
422 double massXY[MAX_Q1D][MAX_Q1D];
423 for (int qy = 0; qy < Q1D; ++qy)
425 for (int qx = 0; qx < Q1D; ++qx)
427 massXY[qy][qx] = 0.0;
431 for (int dy = 0; dy < D1Dy; ++dy)
433 double massX[MAX_Q1D];
434 for (int qx = 0; qx < Q1D; ++qx)
439 for (int dx = 0; dx < D1Dx; ++dx)
441 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
442 for (int qx = 0; qx < Q1D; ++qx)
444 massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
445 ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
449 for (int qy = 0; qy < Q1D; ++qy)
451 const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
452 ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
453 for (int qx = 0; qx < Q1D; ++qx)
455 const double wx = massX[qx];
456 massXY[qy][qx] += wx * wy;
461 for (int qz = 0; qz < Q1D; ++qz)
463 const double wz = trialHcurl ? ((c == 2) ? Bo(qz,dz) : Bc(qz,dz)) :
464 ((c == 2) ? Bc(qz,dz) : Bo(qz,dz));
465 for (int qy = 0; qy < Q1D; ++qy)
467 for (int qx = 0; qx < Q1D; ++qx)
469 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
475 osc += D1Dx * D1Dy * D1Dz;
476 } // loop (c) over components
479 for (int qz = 0; qz < Q1D; ++qz)
481 for (int qy = 0; qy < Q1D; ++qy)
483 for (int qx = 0; qx < Q1D; ++qx)
485 const double O11 = op(0,qx,qy,qz,e);
486 const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,qz,e);
487 const double O13 = scalarCoeff ? 0.0 : op(2,qx,qy,qz,e);
488 const double O21 = scalarCoeff ? 0.0 : op(3,qx,qy,qz,e);
489 const double O22 = scalarCoeff ? O11 : op(4,qx,qy,qz,e);
490 const double O23 = scalarCoeff ? 0.0 : op(5,qx,qy,qz,e);
491 const double O31 = scalarCoeff ? 0.0 : op(6,qx,qy,qz,e);
492 const double O32 = scalarCoeff ? 0.0 : op(7,qx,qy,qz,e);
493 const double O33 = scalarCoeff ? O11 : op(8,qx,qy,qz,e);
494 const double massX = mass[qz][qy][qx][0];
495 const double massY = mass[qz][qy][qx][1];
496 const double massZ = mass[qz][qy][qx][2];
497 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
498 mass[qz][qy][qx][1] = (O21*massX)+(O22*massY)+(O23*massZ);
499 mass[qz][qy][qx][2] = (O31*massX)+(O32*massY)+(O33*massZ);
504 for (int qz = 0; qz < Q1D; ++qz)
506 double massXY[HDIV_MAX_D1D][HDIV_MAX_D1D];
509 for (int c = 0; c < VDIM; ++c) // loop over x, y, z test components
511 const int D1Dz = trialHcurl ? ((c == 2) ? D1Dtest : D1Dtest - 1) :
512 ((c == 2) ? D1Dtest - 1 : D1Dtest);
513 const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
514 ((c == 1) ? D1Dtest - 1 : D1Dtest);
515 const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
516 ((c == 0) ? D1Dtest - 1 : D1Dtest);
518 for (int dy = 0; dy < D1Dy; ++dy)
520 for (int dx = 0; dx < D1Dx; ++dx)
522 massXY[dy][dx] = 0.0;
525 for (int qy = 0; qy < Q1D; ++qy)
527 double massX[HDIV_MAX_D1D];
528 for (int dx = 0; dx < D1Dx; ++dx)
532 for (int qx = 0; qx < Q1D; ++qx)
534 for (int dx = 0; dx < D1Dx; ++dx)
536 massX[dx] += mass[qz][qy][qx][c] * (trialHcurl ?
537 ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
538 ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
541 for (int dy = 0; dy < D1Dy; ++dy)
543 const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
544 ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
545 for (int dx = 0; dx < D1Dx; ++dx)
547 massXY[dy][dx] += massX[dx] * wy;
552 for (int dz = 0; dz < D1Dz; ++dz)
554 const double wz = trialHcurl ? ((c == 2) ? Bct(dz,qz) : Bot(dz,qz)) :
555 ((c == 2) ? Bot(dz,qz) : Bct(dz,qz));
556 for (int dy = 0; dy < D1Dy; ++dy)
558 for (int dx = 0; dx < D1Dx; ++dx)
560 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) +=
566 osc += D1Dx * D1Dy * D1Dz;
569 }); // end of element loop
572 // Mass operator for H(curl) and H(div) functions, using Piola transformations
573 // u = dF^{-T} \hat{u} in H(curl), v = (1 / det dF) dF \hat{v} in H(div).
574 void PAHcurlHdivMassApply2D(const int D1D,
578 const bool scalarCoeff,
579 const bool trialHcurl,
580 const Array<double> &_Bo,
581 const Array<double> &_Bc,
582 const Array<double> &_Bot,
583 const Array<double> &_Bct,
588 constexpr static int MAX_D1D = HCURL_MAX_D1D;
589 constexpr static int MAX_Q1D = HCURL_MAX_Q1D;
591 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D
");
592 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D
");
593 constexpr static int VDIM = 2;
595 auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
596 auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
597 auto Bot = Reshape(_Bot.Read(), D1Dtest-1, Q1D);
598 auto Bct = Reshape(_Bct.Read(), D1Dtest, Q1D);
599 auto op = Reshape(_op.Read(), scalarCoeff ? 1 : 4, Q1D, Q1D, NE);
600 auto x = Reshape(_x.Read(), 2*(D1D-1)*D1D, NE);
601 auto y = Reshape(_y.ReadWrite(), 2*(D1Dtest-1)*D1Dtest, NE);
605 double mass[MAX_Q1D][MAX_Q1D][VDIM];
607 for (int qy = 0; qy < Q1D; ++qy)
609 for (int qx = 0; qx < Q1D; ++qx)
611 for (int c = 0; c < VDIM; ++c)
613 mass[qy][qx][c] = 0.0;
619 for (int c = 0; c < VDIM; ++c) // loop over x, y trial components
621 const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
622 ((c == 1) ? D1D : D1D - 1);
623 const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
624 ((c == 0) ? D1D : D1D - 1);
626 for (int dy = 0; dy < D1Dy; ++dy)
628 double massX[MAX_Q1D];
629 for (int qx = 0; qx < Q1D; ++qx)
634 for (int dx = 0; dx < D1Dx; ++dx)
636 const double t = x(dx + (dy * D1Dx) + osc, e);
637 for (int qx = 0; qx < Q1D; ++qx)
639 massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
640 ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
644 for (int qy = 0; qy < Q1D; ++qy)
646 const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
647 ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
648 for (int qx = 0; qx < Q1D; ++qx)
650 mass[qy][qx][c] += massX[qx] * wy;
656 } // loop (c) over components
659 for (int qy = 0; qy < Q1D; ++qy)
661 for (int qx = 0; qx < Q1D; ++qx)
663 const double O11 = op(0,qx,qy,e);
664 const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,e);
665 const double O21 = scalarCoeff ? 0.0 : op(2,qx,qy,e);
666 const double O22 = scalarCoeff ? O11 : op(3,qx,qy,e);
667 const double massX = mass[qy][qx][0];
668 const double massY = mass[qy][qx][1];
669 mass[qy][qx][0] = (O11*massX)+(O12*massY);
670 mass[qy][qx][1] = (O21*massX)+(O22*massY);
675 for (int c = 0; c < VDIM; ++c) // loop over x, y test components
677 const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
678 ((c == 1) ? D1Dtest - 1 : D1Dtest);
679 const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
680 ((c == 0) ? D1Dtest - 1 : D1Dtest);
682 for (int qy = 0; qy < Q1D; ++qy)
684 double massX[HDIV_MAX_D1D];
685 for (int dx = 0; dx < D1Dx; ++dx)
689 for (int qx = 0; qx < Q1D; ++qx)
691 for (int dx = 0; dx < D1Dx; ++dx)
693 massX[dx] += mass[qy][qx][c] * (trialHcurl ?
694 ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
695 ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
698 for (int dy = 0; dy < D1Dy; ++dy)
700 const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
701 ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
702 for (int dx = 0; dx < D1Dx; ++dx)
704 y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
711 }); // end of element loop
714 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &fes)
716 AssemblePA(fes, fes);
719 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &trial_fes,
720 const FiniteElementSpace &test_fes)
722 // Assumes tensor-product elements
723 Mesh *mesh = trial_fes.GetMesh();
725 const FiniteElement *trial_fel = trial_fes.GetFE(0);
726 const VectorTensorFiniteElement *trial_el =
727 dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
730 const FiniteElement *test_fel = test_fes.GetFE(0);
731 const VectorTensorFiniteElement *test_el =
732 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
735 const IntegrationRule *ir
736 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
737 *mesh->GetElementTransformation(0));
738 const int dims = trial_el->GetDim();
739 MFEM_VERIFY(dims == 2 || dims == 3, "");
741 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
742 const int nq = ir->GetNPoints();
743 dim = mesh->Dimension();
744 MFEM_VERIFY(dim == 2 || dim == 3, "");
746 ne = trial_fes.GetNE();
747 MFEM_VERIFY(ne == test_fes.GetNE(),
748 "Different meshes
for test and trial spaces
");
749 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
750 mapsC = &trial_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
751 mapsO = &trial_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
752 dofs1D = mapsC->ndof;
753 quad1D = mapsC->nqpt;
755 mapsCtest = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
756 mapsOtest = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
757 dofs1Dtest = mapsCtest->ndof;
759 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
761 trial_fetype = trial_el->GetDerivType();
762 test_fetype = test_el->GetDerivType();
764 const int MQsymmDim = MQ ? (MQ->GetWidth() * (MQ->GetWidth() + 1)) / 2 : 0;
765 const int MQfullDim = MQ ? (MQ->GetHeight() * MQ->GetWidth()) : 0;
766 const int MQdim = MQ ? (MQ->IsSymmetric() ? MQsymmDim : MQfullDim) : 0;
767 const int coeffDim = MQ ? MQdim : (VQ ? VQ->GetVDim() : 1);
769 symmetric = MQ ? MQ->IsSymmetric() : true;
771 const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
772 const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
773 const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
774 const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
776 if ((trial_curl && test_div) || (trial_div && test_curl))
777 pa_data.SetSize((coeffDim == 1 ? 1 : dim*dim) * nq * ne,
778 Device::GetMemoryType());
780 pa_data.SetSize((symmetric ? symmDims : MQfullDim) * nq * ne,
781 Device::GetMemoryType());
783 Vector coeff(coeffDim * ne * nq);
785 auto coeffh = Reshape(coeff.HostWrite(), coeffDim, nq, ne);
788 Vector D(VQ ? coeffDim : 0);
795 Msymm.SetSize(MQsymmDim);
805 MFEM_VERIFY(coeffDim == dim, "");
809 MFEM_VERIFY(coeffDim == MQdim, "");
810 MFEM_VERIFY(MQ->GetHeight() == dim && MQ->GetWidth() == dim, "");
813 for (int e=0; e<ne; ++e)
815 ElementTransformation *tr = mesh->GetElementTransformation(e);
816 for (int p=0; p<nq; ++p)
820 if (MQ->IsSymmetric())
822 MQ->EvalSymmetric(Msymm, *tr, ir->IntPoint(p));
824 for (int i=0; i<MQsymmDim; ++i)
826 coeffh(i, p, e) = Msymm[i];
831 MQ->Eval(M, *tr, ir->IntPoint(p));
833 for (int i=0; i<dim; ++i)
834 for (int j=0; j<dim; ++j)
836 coeffh(j+(i*dim), p, e) = M(i,j);
842 VQ->Eval(D, *tr, ir->IntPoint(p));
843 for (int i=0; i<coeffDim; ++i)
845 coeffh(i, p, e) = D[i];
850 coeffh(0, p, e) = Q->Eval(*tr, ir->IntPoint(p));
856 if (trial_curl && test_curl && dim == 3)
858 PADiffusionSetup3D(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
861 else if (trial_curl && test_curl && dim == 2)
863 PADiffusionSetup2D<2>(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
866 else if (trial_div && test_div && dim == 3)
868 PAHdivSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
871 else if (trial_div && test_div && dim == 2)
873 PAHdivSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
876 else if (((trial_curl && test_div) || (trial_div && test_curl)) &&
877 test_fel->GetOrder() == trial_fel->GetOrder())
881 PAHcurlL2Setup(nq, coeffDim, ne, ir->GetWeights(), coeff, pa_data);
885 const bool tr = (trial_div && test_curl);
887 PAHcurlHdivSetup3D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
888 geom->J, coeff, pa_data);
890 PAHcurlHdivSetup2D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
891 geom->J, coeff, pa_data);
896 MFEM_ABORT("Unknown kernel.
");
900 void VectorFEMassIntegrator::AssembleDiagonalPA(Vector& diag)
904 if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
906 if (Device::Allows(Backend::DEVICE_MASK))
908 const int ID = (dofs1D << 4) | quad1D;
911 case 0x23: return SmemPAHcurlMassAssembleDiagonal3D<2,3>(dofs1D, quad1D, ne,
913 mapsO->B, mapsC->B, pa_data, diag);
914 case 0x34: return SmemPAHcurlMassAssembleDiagonal3D<3,4>(dofs1D, quad1D, ne,
916 mapsO->B, mapsC->B, pa_data, diag);
917 case 0x45: return SmemPAHcurlMassAssembleDiagonal3D<4,5>(dofs1D, quad1D, ne,
919 mapsO->B, mapsC->B, pa_data, diag);
920 case 0x56: return SmemPAHcurlMassAssembleDiagonal3D<5,6>(dofs1D, quad1D, ne,
922 mapsO->B, mapsC->B, pa_data, diag);
923 default: return SmemPAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
924 mapsO->B, mapsC->B, pa_data, diag);
928 PAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
929 mapsO->B, mapsC->B, pa_data, diag);
931 else if (trial_fetype == mfem::FiniteElement::DIV &&
932 test_fetype == trial_fetype)
934 PAHdivMassAssembleDiagonal3D(dofs1D, quad1D, ne,
935 mapsO->B, mapsC->B, pa_data, diag);
939 MFEM_ABORT("Unknown kernel.
");
944 if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
946 PAHcurlMassAssembleDiagonal2D(dofs1D, quad1D, ne, symmetric,
947 mapsO->B, mapsC->B, pa_data, diag);
949 else if (trial_fetype == mfem::FiniteElement::DIV &&
950 test_fetype == trial_fetype)
952 PAHdivMassAssembleDiagonal2D(dofs1D, quad1D, ne,
953 mapsO->B, mapsC->B, pa_data, diag);
957 MFEM_ABORT("Unknown kernel.
");
962 void VectorFEMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
964 const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
965 const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
966 const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
967 const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
971 if (trial_curl && test_curl)
973 if (Device::Allows(Backend::DEVICE_MASK))
975 const int ID = (dofs1D << 4) | quad1D;
978 case 0x23: return SmemPAHcurlMassApply3D<2,3>(dofs1D, quad1D, ne, symmetric,
981 mapsC->Bt, pa_data, x, y);
982 case 0x34: return SmemPAHcurlMassApply3D<3,4>(dofs1D, quad1D, ne, symmetric,
985 mapsC->Bt, pa_data, x, y);
986 case 0x45: return SmemPAHcurlMassApply3D<4,5>(dofs1D, quad1D, ne, symmetric,
989 mapsC->Bt, pa_data, x, y);
990 case 0x56: return SmemPAHcurlMassApply3D<5,6>(dofs1D, quad1D, ne, symmetric,
993 mapsC->Bt, pa_data, x, y);
994 default: return SmemPAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B,
996 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1000 PAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B, mapsO->Bt,
1001 mapsC->Bt, pa_data, x, y);
1003 else if (trial_div && test_div)
1005 PAHdivMassApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1006 mapsC->Bt, pa_data, x, y);
1008 else if (trial_curl && test_div)
1010 const bool scalarCoeff = !(VQ || MQ);
1011 PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1012 true, mapsO->B, mapsC->B, mapsOtest->Bt,
1013 mapsCtest->Bt, pa_data, x, y);
1015 else if (trial_div && test_curl)
1017 const bool scalarCoeff = !(VQ || MQ);
1018 PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1019 false, mapsO->B, mapsC->B, mapsOtest->Bt,
1020 mapsCtest->Bt, pa_data, x, y);
1024 MFEM_ABORT("Unknown kernel.
");
1029 if (trial_curl && test_curl)
1031 PAHcurlMassApply2D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B,
1032 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1034 else if (trial_div && test_div)
1036 PAHdivMassApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1037 mapsC->Bt, pa_data, x, y);
1039 else if ((trial_curl && test_div) || (trial_div && test_curl))
1041 const bool scalarCoeff = !(VQ || MQ);
1042 PAHcurlHdivMassApply2D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1043 trial_curl, mapsO->B, mapsC->B, mapsOtest->Bt,
1044 mapsCtest->Bt, pa_data, x, y);
1048 MFEM_ABORT("Unknown kernel.
");
1053 void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace
1055 const FiniteElementSpace &test_fes)
1057 // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1058 Mesh *mesh = trial_fes.GetMesh();
1059 const FiniteElement *trial_fel = trial_fes.GetFE(0);
1060 const FiniteElement *test_fel = test_fes.GetFE(0);
1062 const NodalTensorFiniteElement *trial_el =
1063 dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1066 const VectorTensorFiniteElement *test_el =
1067 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1070 const IntegrationRule *ir
1071 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
1072 *mesh->GetElementTransformation(0));
1073 const int dims = trial_el->GetDim();
1074 MFEM_VERIFY(dims == 2 || dims == 3, "");
1076 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
1077 const int nq = ir->GetNPoints();
1078 dim = mesh->Dimension();
1079 MFEM_VERIFY(dim == 2 || dim == 3, "");
1081 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1083 ne = trial_fes.GetNE();
1084 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
1085 mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
1086 mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
1087 dofs1D = mapsC->ndof;
1088 quad1D = mapsC->nqpt;
1090 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
1092 pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
1094 Vector coeff(ne * nq);
1098 for (int e=0; e<ne; ++e)
1100 ElementTransformation *tr = mesh->GetElementTransformation(e);
1101 for (int p=0; p<nq; ++p)
1103 coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
1108 // Use the same setup functions as VectorFEMassIntegrator.
1109 if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
1111 PADiffusionSetup3D(quad1D, 1, ne, ir->GetWeights(), geom->J,
1114 else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
1116 PADiffusionSetup2D<2>(quad1D, 1, ne, ir->GetWeights(), geom->J,
1121 MFEM_ABORT("Unknown kernel.
");
1125 void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const
1128 PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1129 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1131 PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1132 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1135 MFEM_ABORT("Unsupported dimension!
");
void PAHcurlH1Apply3D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &gc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHcurlMassApply3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHcurlHdivSetup2D(const int Q1D, const int coeffDim, const int NE, const bool transpose, const Array< double > &_w, const Vector &j, Vector &_coeff, Vector &op)
void PAHdivMassAssembleDiagonal2D(const int D1D, const int Q1D, const int NE, const Array< double > &_Bo, const Array< double > &_Bc, const Vector &_op, Vector &_diag)
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void PAHcurlL2Setup(const int NQ, const int coeffDim, const int NE, const Array< double > &w, Vector &coeff, Vector &op)
void PAHcurlHdivMassApply3D(const int D1D, const int D1Dtest, const int Q1D, const int NE, const bool scalarCoeff, const bool trialHcurl, const Array< double > &_Bo, const Array< double > &_Bc, const Array< double > &_Bot, const Array< double > &_Bct, const Vector &_op, const Vector &_x, Vector &_y)
void PAHdivMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const Array< double > &_Bo, const Array< double > &_Bc, const Vector &_op, Vector &_diag)
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void SmemPAHcurlMassApply3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHdivMassApply3D(const int D1D, const int Q1D, const int NE, const Array< double > &_Bo, const Array< double > &_Bc, const Array< double > &_Bot, const Array< double > &_Bct, const Vector &_op, const Vector &_x, Vector &_y)
void PAHcurlMassAssembleDiagonal2D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
constexpr int HCURL_MAX_D1D
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
void PAHdivSetup3D(const int Q1D, const int NE, const Array< double > &w, const Vector &j, Vector &_coeff, Vector &op)
void PAHdivMassApply2D(const int D1D, const int Q1D, const int NE, const Array< double > &_Bo, const Array< double > &_Bc, const Array< double > &_Bot, const Array< double > &_Bct, const Vector &_op, const Vector &_x, Vector &_y)
void PAHdivSetup2D(const int Q1D, const int NE, const Array< double > &w, const Vector &j, Vector &_coeff, Vector &op)
void PAHcurlMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
void PADiffusionSetup3D(const int Q1D, const int coeffDim, const int NE, const Array< double > &w, const Vector &j, const Vector &c, Vector &d)
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void PAHcurlMassApply2D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void SmemPAHcurlMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
void PAHcurlHdivSetup3D(const int Q1D, const int coeffDim, const int NE, const bool transpose, const Array< double > &_w, const Vector &j, Vector &_coeff, Vector &op)
void PAHcurlH1Apply2D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &gc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
constexpr int HCURL_MAX_Q1D