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> &bo,
121 const Array<double> &bct,
122 const Array<double> &gct,
123 const Vector &pa_data,
130 const Array<double> &bc,
131 const Array<double> &gc,
132 const Array<double> &bot,
133 const Array<double> &bct,
134 const Vector &pa_data,
141 const Array<double> &bc,
142 const Array<double> &bo,
143 const Array<double> &bct,
144 const Array<double> &gct,
145 const Vector &pa_data,
152 const Array<double> &Bo_,
153 const Array<double> &Bc_,
160 const Array<double> &Bo_,
161 const Array<double> &Bc_,
168 const Array<double> &Bo_,
169 const Array<double> &Bc_,
170 const Array<double> &Bot_,
171 const Array<double> &Bct_,
179 const Array<double> &Bo_,
180 const Array<double> &Bc_,
181 const Array<double> &Bot_,
182 const Array<double> &Bct_,
190 const Array<double> &w,
200 const bool transpose,
206 const bool symmetric = (coeffDim != 9);
208 auto J =
Reshape(j.
Read(), Q1D, Q1D, Q1D, 3, 3, NE);
209 auto coeff =
Reshape(coeff_.
Read(), coeffDim, Q1D, Q1D, Q1D, NE);
213 const int i12 = transpose ? 3 : 1;
214 const int i13 = transpose ? 6 : 2;
215 const int i21 = transpose ? 1 : 3;
217 const int i23 = transpose ? 7 : 5;
218 const int i31 = transpose ? 2 : 6;
219 const int i32 = transpose ? 5 : 7;
222 MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
224 MFEM_FOREACH_THREAD(qx,x,Q1D)
226 MFEM_FOREACH_THREAD(qy,y,Q1D)
228 MFEM_FOREACH_THREAD(qz,z,Q1D)
230 const double J11 = J(qx,qy,qz,0,0,e);
231 const double J21 = J(qx,qy,qz,1,0,e);
232 const double J31 = J(qx,qy,qz,2,0,e);
233 const double J12 = J(qx,qy,qz,0,1,e);
234 const double J22 = J(qx,qy,qz,1,1,e);
235 const double J32 = J(qx,qy,qz,2,1,e);
236 const double J13 = J(qx,qy,qz,0,2,e);
237 const double J23 = J(qx,qy,qz,1,2,e);
238 const double J33 = J(qx,qy,qz,2,2,e);
239 const double detJ = J11 * (J22 * J33 - J32 * J23) -
240 J21 * (J12 * J33 - J32 * J13) +
241 J31 * (J12 * J23 - J22 * J13);
242 const double w_detJ = W(qx,qy,qz) / detJ;
244 const double A11 = (J22 * J33) - (J23 * J32);
245 const double A12 = (J32 * J13) - (J12 * J33);
246 const double A13 = (J12 * J23) - (J22 * J13);
247 const double A21 = (J31 * J23) - (J21 * J33);
248 const double A22 = (J11 * J33) - (J13 * J31);
249 const double A23 = (J21 * J13) - (J11 * J23);
250 const double A31 = (J21 * J32) - (J31 * J22);
251 const double A32 = (J31 * J12) - (J11 * J32);
252 const double A33 = (J11 * J22) - (J12 * J21);
254 if (coeffDim == 6 || coeffDim == 9)
257 const double M11 = (!symmetric) ? coeff(i11,qx,qy,qz,e) : coeff(0,qx,qy,qz,e);
258 const double M12 = (!symmetric) ? coeff(i12,qx,qy,qz,e) : coeff(1,qx,qy,qz,e);
259 const double M13 = (!symmetric) ? coeff(i13,qx,qy,qz,e) : coeff(2,qx,qy,qz,e);
260 const double M21 = (!symmetric) ? coeff(i21,qx,qy,qz,e) : M12;
261 const double M22 = (!symmetric) ? coeff(i22,qx,qy,qz,e) : coeff(3,qx,qy,qz,e);
262 const double M23 = (!symmetric) ? coeff(i23,qx,qy,qz,e) : coeff(4,qx,qy,qz,e);
263 const double M31 = (!symmetric) ? coeff(i31,qx,qy,qz,e) : M13;
264 const double M32 = (!symmetric) ? coeff(i32,qx,qy,qz,e) : M23;
265 const double M33 = (!symmetric) ? coeff(i33,qx,qy,qz,e) : coeff(5,qx,qy,qz,e);
267 const double R11 = M11*J11 + M21*J21 + M31*J31;
268 const double R12 = M11*J12 + M21*J22 + M31*J32;
269 const double R13 = M11*J13 + M21*J23 + M31*J33;
270 const double R21 = M12*J11 + M22*J21 + M32*J31;
271 const double R22 = M12*J12 + M22*J22 + M32*J32;
272 const double R23 = M12*J13 + M22*J23 + M32*J33;
273 const double R31 = M13*J11 + M23*J21 + M33*J31;
274 const double R32 = M13*J12 + M23*J22 + M33*J32;
275 const double R33 = M13*J13 + M23*J23 + M33*J33;
278 y(i11,qx,qy,qz,e) = w_detJ * (A11*R11 + A12*R21 + A13*R31);
279 y(i21,qx,qy,qz,e) = w_detJ * (A11*R12 + A12*R22 + A13*R32);
280 y(i31,qx,qy,qz,e) = w_detJ * (A11*R13 + A12*R23 + A13*R33);
281 y(i12,qx,qy,qz,e) = w_detJ * (A21*R11 + A22*R21 + A23*R31);
282 y(i22,qx,qy,qz,e) = w_detJ * (A21*R12 + A22*R22 + A23*R32);
283 y(i32,qx,qy,qz,e) = w_detJ * (A21*R13 + A22*R23 + A23*R33);
284 y(i13,qx,qy,qz,e) = w_detJ * (A31*R11 + A32*R21 + A33*R31);
285 y(i23,qx,qy,qz,e) = w_detJ * (A31*R12 + A32*R22 + A33*R32);
286 y(i33,qx,qy,qz,e) = w_detJ * (A31*R13 + A32*R23 + A33*R33);
288 else if (coeffDim == 3)
290 const double D1 = coeff(0,qx,qy,qz,e);
291 const double D2 = coeff(1,qx,qy,qz,e);
292 const double D3 = coeff(2,qx,qy,qz,e);
295 y(i11,qx,qy,qz,e) = w_detJ * (D1*A11*J11 + D2*A12*J21 + D3*A13*J31);
296 y(i21,qx,qy,qz,e) = w_detJ * (D1*A11*J12 + D2*A12*J22 + D3*A13*J32);
297 y(i31,qx,qy,qz,e) = w_detJ * (D1*A11*J13 + D2*A12*J23 + D3*A13*J33);
298 y(i12,qx,qy,qz,e) = w_detJ * (D1*A21*J11 + D2*A22*J21 + D3*A23*J31);
299 y(i22,qx,qy,qz,e) = w_detJ * (D1*A21*J12 + D2*A22*J22 + D3*A23*J32);
300 y(i32,qx,qy,qz,e) = w_detJ * (D1*A21*J13 + D2*A22*J23 + D3*A23*J33);
301 y(i13,qx,qy,qz,e) = w_detJ * (D1*A31*J11 + D2*A32*J21 + D3*A33*J31);
302 y(i23,qx,qy,qz,e) = w_detJ * (D1*A31*J12 + D2*A32*J22 + D3*A33*J32);
303 y(i33,qx,qy,qz,e) = w_detJ * (D1*A31*J13 + D2*A32*J23 + D3*A33*J33);
317 const bool transpose,
323 const bool symmetric = (coeffDim != 4);
326 auto coeff =
Reshape(coeff_.
Read(), coeffDim, Q1D, Q1D, NE);
330 const int i12 = transpose ? 2 : 1;
331 const int i21 = transpose ? 1 : 2;
334 MFEM_FORALL_2D(e, NE, Q1D, Q1D, 1,
336 MFEM_FOREACH_THREAD(qx,x,Q1D)
338 MFEM_FOREACH_THREAD(qy,y,Q1D)
340 const double J11 = J(qx,qy,0,0,e);
341 const double J21 = J(qx,qy,1,0,e);
342 const double J12 = J(qx,qy,0,1,e);
343 const double J22 = J(qx,qy,1,1,e);
344 const double w_detJ = W(qx,qy) / ((J11*J22) - (J21*J12));
346 if (coeffDim == 3 || coeffDim == 4)
349 const double M11 = coeff(i11,qx,qy,e);
350 const double M12 = (!symmetric) ? coeff(i12,qx,qy,e) : coeff(1,qx,qy,e);
351 const double M21 = (!symmetric) ? coeff(i21,qx,qy,e) : M12;
352 const double M22 = (!symmetric) ? coeff(i22,qx,qy,e) : coeff(2,qx,qy,e);
355 const double R11 = ( J22*M11 - J12*M12);
356 const double R12 = ( J22*M21 - J12*M22);
357 const double R21 = (-J21*M11 + J11*M12);
358 const double R22 = (-J21*M21 + J11*M22);
361 y(i11,qx,qy,e) = w_detJ * (R11*J11 + R12*J21);
362 y(i21,qx,qy,e) = w_detJ * (R11*J12 + R12*J22);
363 y(i12,qx,qy,e) = w_detJ * (R21*J11 + R22*J21);
364 y(i22,qx,qy,e) = w_detJ * (R21*J12 + R22*J22);
366 else if (coeffDim == 2)
368 const double D1 = coeff(0,qx,qy,e);
369 const double D2 = coeff(1,qx,qy,e);
370 const double R11 = D1*J11;
371 const double R12 = D1*J12;
372 const double R21 = D2*J21;
373 const double R22 = D2*J22;
374 y(i11,qx,qy,e) = w_detJ * ( J22*R11 - J12*R21);
375 y(i12,qx,qy,e) = w_detJ * ( J22*R12 - J12*R22);
376 y(i21,qx,qy,e) = w_detJ * (-J21*R11 + J11*R21);
377 y(i22,qx,qy,e) = w_detJ * (-J21*R12 + J11*R22);
390 const bool scalarCoeff,
391 const bool trialHcurl,
403 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D
");
404 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D
");
405 constexpr static int VDIM = 3;
407 auto Bo = Reshape(Bo_.Read(), Q1D, D1D-1);
408 auto Bc = Reshape(Bc_.Read(), Q1D, D1D);
409 auto Bot = Reshape(Bot_.Read(), D1Dtest-1, Q1D);
410 auto Bct = Reshape(Bct_.Read(), D1Dtest, Q1D);
411 auto op = Reshape(op_.Read(), scalarCoeff ? 1 : 9, Q1D, Q1D, Q1D, NE);
412 auto x = Reshape(x_.Read(), 3*(D1D-1)*D1D*(trialHcurl ? D1D : D1D-1), NE);
413 auto y = Reshape(y_.ReadWrite(), 3*(D1Dtest-1)*D1Dtest*
414 (trialHcurl ? D1Dtest-1 : D1Dtest), NE);
418 double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
420 for (int qz = 0; qz < Q1D; ++qz)
422 for (int qy = 0; qy < Q1D; ++qy)
424 for (int qx = 0; qx < Q1D; ++qx)
426 for (int c = 0; c < VDIM; ++c)
428 mass[qz][qy][qx][c] = 0.0;
435 for (int c = 0; c < VDIM; ++c) // loop over x, y, z trial components
437 const int D1Dz = trialHcurl ? ((c == 2) ? D1D - 1 : D1D) :
438 ((c == 2) ? D1D : D1D - 1);
439 const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
440 ((c == 1) ? D1D : D1D - 1);
441 const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
442 ((c == 0) ? D1D : D1D - 1);
444 for (int dz = 0; dz < D1Dz; ++dz)
446 double massXY[MAX_Q1D][MAX_Q1D];
447 for (int qy = 0; qy < Q1D; ++qy)
449 for (int qx = 0; qx < Q1D; ++qx)
451 massXY[qy][qx] = 0.0;
455 for (int dy = 0; dy < D1Dy; ++dy)
457 double massX[MAX_Q1D];
458 for (int qx = 0; qx < Q1D; ++qx)
463 for (int dx = 0; dx < D1Dx; ++dx)
465 const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
466 for (int qx = 0; qx < Q1D; ++qx)
468 massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
469 ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
473 for (int qy = 0; qy < Q1D; ++qy)
475 const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
476 ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
477 for (int qx = 0; qx < Q1D; ++qx)
479 const double wx = massX[qx];
480 massXY[qy][qx] += wx * wy;
485 for (int qz = 0; qz < Q1D; ++qz)
487 const double wz = trialHcurl ? ((c == 2) ? Bo(qz,dz) : Bc(qz,dz)) :
488 ((c == 2) ? Bc(qz,dz) : Bo(qz,dz));
489 for (int qy = 0; qy < Q1D; ++qy)
491 for (int qx = 0; qx < Q1D; ++qx)
493 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
499 osc += D1Dx * D1Dy * D1Dz;
500 } // loop (c) over components
503 for (int qz = 0; qz < Q1D; ++qz)
505 for (int qy = 0; qy < Q1D; ++qy)
507 for (int qx = 0; qx < Q1D; ++qx)
509 const double O11 = op(0,qx,qy,qz,e);
510 const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,qz,e);
511 const double O13 = scalarCoeff ? 0.0 : op(2,qx,qy,qz,e);
512 const double O21 = scalarCoeff ? 0.0 : op(3,qx,qy,qz,e);
513 const double O22 = scalarCoeff ? O11 : op(4,qx,qy,qz,e);
514 const double O23 = scalarCoeff ? 0.0 : op(5,qx,qy,qz,e);
515 const double O31 = scalarCoeff ? 0.0 : op(6,qx,qy,qz,e);
516 const double O32 = scalarCoeff ? 0.0 : op(7,qx,qy,qz,e);
517 const double O33 = scalarCoeff ? O11 : op(8,qx,qy,qz,e);
518 const double massX = mass[qz][qy][qx][0];
519 const double massY = mass[qz][qy][qx][1];
520 const double massZ = mass[qz][qy][qx][2];
521 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
522 mass[qz][qy][qx][1] = (O21*massX)+(O22*massY)+(O23*massZ);
523 mass[qz][qy][qx][2] = (O31*massX)+(O32*massY)+(O33*massZ);
528 for (int qz = 0; qz < Q1D; ++qz)
530 double massXY[HDIV_MAX_D1D][HDIV_MAX_D1D];
533 for (int c = 0; c < VDIM; ++c) // loop over x, y, z test components
535 const int D1Dz = trialHcurl ? ((c == 2) ? D1Dtest : D1Dtest - 1) :
536 ((c == 2) ? D1Dtest - 1 : D1Dtest);
537 const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
538 ((c == 1) ? D1Dtest - 1 : D1Dtest);
539 const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
540 ((c == 0) ? D1Dtest - 1 : D1Dtest);
542 for (int dy = 0; dy < D1Dy; ++dy)
544 for (int dx = 0; dx < D1Dx; ++dx)
546 massXY[dy][dx] = 0.0;
549 for (int qy = 0; qy < Q1D; ++qy)
551 double massX[HDIV_MAX_D1D];
552 for (int dx = 0; dx < D1Dx; ++dx)
556 for (int qx = 0; qx < Q1D; ++qx)
558 for (int dx = 0; dx < D1Dx; ++dx)
560 massX[dx] += mass[qz][qy][qx][c] * (trialHcurl ?
561 ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
562 ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
565 for (int dy = 0; dy < D1Dy; ++dy)
567 const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
568 ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
569 for (int dx = 0; dx < D1Dx; ++dx)
571 massXY[dy][dx] += massX[dx] * wy;
576 for (int dz = 0; dz < D1Dz; ++dz)
578 const double wz = trialHcurl ? ((c == 2) ? Bct(dz,qz) : Bot(dz,qz)) :
579 ((c == 2) ? Bot(dz,qz) : Bct(dz,qz));
580 for (int dy = 0; dy < D1Dy; ++dy)
582 for (int dx = 0; dx < D1Dx; ++dx)
584 y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) +=
590 osc += D1Dx * D1Dy * D1Dz;
593 }); // end of element loop
596 // Mass operator for H(curl) and H(div) functions, using Piola transformations
597 // u = dF^{-T} \hat{u} in H(curl), v = (1 / det dF) dF \hat{v} in H(div).
598 void PAHcurlHdivMassApply2D(const int D1D,
602 const bool scalarCoeff,
603 const bool trialHcurl,
604 const Array<double> &Bo_,
605 const Array<double> &Bc_,
606 const Array<double> &Bot_,
607 const Array<double> &Bct_,
612 constexpr static int MAX_D1D = HCURL_MAX_D1D;
613 constexpr static int MAX_Q1D = HCURL_MAX_Q1D;
615 MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D
");
616 MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D
");
617 constexpr static int VDIM = 2;
619 auto Bo = Reshape(Bo_.Read(), Q1D, D1D-1);
620 auto Bc = Reshape(Bc_.Read(), Q1D, D1D);
621 auto Bot = Reshape(Bot_.Read(), D1Dtest-1, Q1D);
622 auto Bct = Reshape(Bct_.Read(), D1Dtest, Q1D);
623 auto op = Reshape(op_.Read(), scalarCoeff ? 1 : 4, Q1D, Q1D, NE);
624 auto x = Reshape(x_.Read(), 2*(D1D-1)*D1D, NE);
625 auto y = Reshape(y_.ReadWrite(), 2*(D1Dtest-1)*D1Dtest, NE);
629 double mass[MAX_Q1D][MAX_Q1D][VDIM];
631 for (int qy = 0; qy < Q1D; ++qy)
633 for (int qx = 0; qx < Q1D; ++qx)
635 for (int c = 0; c < VDIM; ++c)
637 mass[qy][qx][c] = 0.0;
643 for (int c = 0; c < VDIM; ++c) // loop over x, y trial components
645 const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
646 ((c == 1) ? D1D : D1D - 1);
647 const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
648 ((c == 0) ? D1D : D1D - 1);
650 for (int dy = 0; dy < D1Dy; ++dy)
652 double massX[MAX_Q1D];
653 for (int qx = 0; qx < Q1D; ++qx)
658 for (int dx = 0; dx < D1Dx; ++dx)
660 const double t = x(dx + (dy * D1Dx) + osc, e);
661 for (int qx = 0; qx < Q1D; ++qx)
663 massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
664 ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
668 for (int qy = 0; qy < Q1D; ++qy)
670 const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
671 ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
672 for (int qx = 0; qx < Q1D; ++qx)
674 mass[qy][qx][c] += massX[qx] * wy;
680 } // loop (c) over components
683 for (int qy = 0; qy < Q1D; ++qy)
685 for (int qx = 0; qx < Q1D; ++qx)
687 const double O11 = op(0,qx,qy,e);
688 const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,e);
689 const double O21 = scalarCoeff ? 0.0 : op(2,qx,qy,e);
690 const double O22 = scalarCoeff ? O11 : op(3,qx,qy,e);
691 const double massX = mass[qy][qx][0];
692 const double massY = mass[qy][qx][1];
693 mass[qy][qx][0] = (O11*massX)+(O12*massY);
694 mass[qy][qx][1] = (O21*massX)+(O22*massY);
699 for (int c = 0; c < VDIM; ++c) // loop over x, y test components
701 const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
702 ((c == 1) ? D1Dtest - 1 : D1Dtest);
703 const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
704 ((c == 0) ? D1Dtest - 1 : D1Dtest);
706 for (int qy = 0; qy < Q1D; ++qy)
708 double massX[HDIV_MAX_D1D];
709 for (int dx = 0; dx < D1Dx; ++dx)
713 for (int qx = 0; qx < Q1D; ++qx)
715 for (int dx = 0; dx < D1Dx; ++dx)
717 massX[dx] += mass[qy][qx][c] * (trialHcurl ?
718 ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
719 ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
722 for (int dy = 0; dy < D1Dy; ++dy)
724 const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
725 ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
726 for (int dx = 0; dx < D1Dx; ++dx)
728 y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
735 }); // end of element loop
738 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &fes)
740 AssemblePA(fes, fes);
743 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &trial_fes,
744 const FiniteElementSpace &test_fes)
746 // Assumes tensor-product elements
747 Mesh *mesh = trial_fes.GetMesh();
749 const FiniteElement *trial_fel = trial_fes.GetFE(0);
750 const VectorTensorFiniteElement *trial_el =
751 dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
754 const FiniteElement *test_fel = test_fes.GetFE(0);
755 const VectorTensorFiniteElement *test_el =
756 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
759 const IntegrationRule *ir
760 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
761 *mesh->GetElementTransformation(0));
762 const int dims = trial_el->GetDim();
763 MFEM_VERIFY(dims == 2 || dims == 3, "");
765 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
766 nq = ir->GetNPoints();
767 dim = mesh->Dimension();
768 MFEM_VERIFY(dim == 2 || dim == 3, "");
770 ne = trial_fes.GetNE();
771 MFEM_VERIFY(ne == test_fes.GetNE(),
772 "Different meshes
for test and trial spaces
");
773 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
774 mapsC = &trial_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
775 mapsO = &trial_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
776 dofs1D = mapsC->ndof;
777 quad1D = mapsC->nqpt;
779 mapsCtest = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
780 mapsOtest = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
781 dofs1Dtest = mapsCtest->ndof;
783 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
785 trial_fetype = trial_el->GetDerivType();
786 test_fetype = test_el->GetDerivType();
788 const int MQsymmDim = SMQ ? (SMQ->GetSize() * (SMQ->GetSize() + 1)) / 2 : 0;
789 const int MQfullDim = MQ ? (MQ->GetHeight() * MQ->GetWidth()) : 0;
790 const int MQdim = MQ ? MQfullDim : MQsymmDim;
791 const int coeffDim = (MQ || SMQ) ? MQdim : (DQ ? DQ->GetVDim() : 1);
793 symmetric = (MQ == NULL);
795 const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
796 const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
797 const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
798 const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
800 if ((trial_curl && test_div) || (trial_div && test_curl))
801 pa_data.SetSize((coeffDim == 1 ? 1 : dim*dim) * nq * ne,
802 Device::GetMemoryType());
804 pa_data.SetSize((symmetric ? symmDims : MQfullDim) * nq * ne,
805 Device::GetMemoryType());
807 Vector coeff(coeffDim * ne * nq);
809 auto coeffh = Reshape(coeff.HostWrite(), coeffDim, nq, ne);
810 if (Q || DQ || MQ || SMQ)
812 Vector DM(DQ ? coeffDim : 0);
814 DenseSymmetricMatrix SM;
818 MFEM_VERIFY(coeffDim == dim, "");
822 MFEM_VERIFY(coeffDim == MQdim, "");
823 MFEM_VERIFY(MQ->GetHeight() == dim && MQ->GetWidth() == dim, "");
828 MFEM_VERIFY(SMQ->GetSize() == dim, "");
832 for (int e=0; e<ne; ++e)
834 ElementTransformation *tr = mesh->GetElementTransformation(e);
835 for (int p=0; p<nq; ++p)
839 MQ->Eval(M, *tr, ir->IntPoint(p));
841 for (int i=0; i<dim; ++i)
842 for (int j=0; j<dim; ++j)
844 coeffh(j+(i*dim), p, e) = M(i,j);
849 SMQ->Eval(SM, *tr, ir->IntPoint(p));
851 for (int i=0; i<dim; ++i)
852 for (int j=i; j<dim; ++j, ++cnt)
854 coeffh(cnt, p, e) = SM(i,j);
859 DQ->Eval(DM, *tr, ir->IntPoint(p));
860 for (int i=0; i<coeffDim; ++i)
862 coeffh(i, p, e) = DM[i];
867 coeffh(0, p, e) = Q->Eval(*tr, ir->IntPoint(p));
873 if (trial_curl && test_curl && dim == 3)
875 PADiffusionSetup3D(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
878 else if (trial_curl && test_curl && dim == 2)
880 PADiffusionSetup2D<2>(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
883 else if (trial_div && test_div && dim == 3)
885 PAHdivSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
888 else if (trial_div && test_div && dim == 2)
890 PAHdivSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
893 else if (((trial_curl && test_div) || (trial_div && test_curl)) &&
894 test_fel->GetOrder() == trial_fel->GetOrder())
898 PAHcurlL2Setup(nq, coeffDim, ne, ir->GetWeights(), coeff, pa_data);
902 const bool tr = (trial_div && test_curl);
904 PAHcurlHdivSetup3D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
905 geom->J, coeff, pa_data);
907 PAHcurlHdivSetup2D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
908 geom->J, coeff, pa_data);
913 MFEM_ABORT("Unknown kernel.
");
917 void VectorFEMassIntegrator::AssembleDiagonalPA(Vector& diag)
921 if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
923 if (Device::Allows(Backend::DEVICE_MASK))
925 const int ID = (dofs1D << 4) | quad1D;
928 case 0x23: return SmemPAHcurlMassAssembleDiagonal3D<2,3>(dofs1D, quad1D, ne,
930 mapsO->B, mapsC->B, pa_data, diag);
931 case 0x34: return SmemPAHcurlMassAssembleDiagonal3D<3,4>(dofs1D, quad1D, ne,
933 mapsO->B, mapsC->B, pa_data, diag);
934 case 0x45: return SmemPAHcurlMassAssembleDiagonal3D<4,5>(dofs1D, quad1D, ne,
936 mapsO->B, mapsC->B, pa_data, diag);
937 case 0x56: return SmemPAHcurlMassAssembleDiagonal3D<5,6>(dofs1D, quad1D, ne,
939 mapsO->B, mapsC->B, pa_data, diag);
940 default: return SmemPAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
941 mapsO->B, mapsC->B, pa_data, diag);
945 PAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
946 mapsO->B, mapsC->B, pa_data, diag);
948 else if (trial_fetype == mfem::FiniteElement::DIV &&
949 test_fetype == trial_fetype)
951 PAHdivMassAssembleDiagonal3D(dofs1D, quad1D, ne,
952 mapsO->B, mapsC->B, pa_data, diag);
956 MFEM_ABORT("Unknown kernel.
");
961 if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
963 PAHcurlMassAssembleDiagonal2D(dofs1D, quad1D, ne, symmetric,
964 mapsO->B, mapsC->B, pa_data, diag);
966 else if (trial_fetype == mfem::FiniteElement::DIV &&
967 test_fetype == trial_fetype)
969 PAHdivMassAssembleDiagonal2D(dofs1D, quad1D, ne,
970 mapsO->B, mapsC->B, pa_data, diag);
974 MFEM_ABORT("Unknown kernel.
");
979 void VectorFEMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
981 const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
982 const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
983 const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
984 const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
988 if (trial_curl && test_curl)
990 if (Device::Allows(Backend::DEVICE_MASK))
992 const int ID = (dofs1D << 4) | quad1D;
995 case 0x23: return SmemPAHcurlMassApply3D<2,3>(dofs1D, quad1D, ne, symmetric,
998 mapsC->Bt, pa_data, x, y);
999 case 0x34: return SmemPAHcurlMassApply3D<3,4>(dofs1D, quad1D, ne, symmetric,
1001 mapsC->B, mapsO->Bt,
1002 mapsC->Bt, pa_data, x, y);
1003 case 0x45: return SmemPAHcurlMassApply3D<4,5>(dofs1D, quad1D, ne, symmetric,
1005 mapsC->B, mapsO->Bt,
1006 mapsC->Bt, pa_data, x, y);
1007 case 0x56: return SmemPAHcurlMassApply3D<5,6>(dofs1D, quad1D, ne, symmetric,
1009 mapsC->B, mapsO->Bt,
1010 mapsC->Bt, pa_data, x, y);
1011 default: return SmemPAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B,
1013 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1017 PAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B, mapsO->Bt,
1018 mapsC->Bt, pa_data, x, y);
1020 else if (trial_div && test_div)
1022 PAHdivMassApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1023 mapsC->Bt, pa_data, x, y);
1025 else if (trial_curl && test_div)
1027 const bool scalarCoeff = !(DQ || MQ || SMQ);
1028 PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1029 true, mapsO->B, mapsC->B, mapsOtest->Bt,
1030 mapsCtest->Bt, pa_data, x, y);
1032 else if (trial_div && test_curl)
1034 const bool scalarCoeff = !(DQ || MQ || SMQ);
1035 PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1036 false, mapsO->B, mapsC->B, mapsOtest->Bt,
1037 mapsCtest->Bt, pa_data, x, y);
1041 MFEM_ABORT("Unknown kernel.
");
1046 if (trial_curl && test_curl)
1048 PAHcurlMassApply2D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B,
1049 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1051 else if (trial_div && test_div)
1053 PAHdivMassApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1054 mapsC->Bt, pa_data, x, y);
1056 else if ((trial_curl && test_div) || (trial_div && test_curl))
1058 const bool scalarCoeff = !(DQ || MQ || SMQ);
1059 PAHcurlHdivMassApply2D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1060 trial_curl, mapsO->B, mapsC->B, mapsOtest->Bt,
1061 mapsCtest->Bt, pa_data, x, y);
1065 MFEM_ABORT("Unknown kernel.
");
1070 void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace
1072 const FiniteElementSpace &test_fes)
1074 // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1075 Mesh *mesh = trial_fes.GetMesh();
1076 const FiniteElement *trial_fel = trial_fes.GetFE(0);
1077 const FiniteElement *test_fel = test_fes.GetFE(0);
1079 const NodalTensorFiniteElement *trial_el =
1080 dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1083 const VectorTensorFiniteElement *test_el =
1084 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1087 const IntegrationRule *ir
1088 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
1089 *mesh->GetElementTransformation(0));
1090 const int dims = trial_el->GetDim();
1091 MFEM_VERIFY(dims == 2 || dims == 3, "");
1093 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
1094 const int nq = ir->GetNPoints();
1095 dim = mesh->Dimension();
1096 MFEM_VERIFY(dim == 2 || dim == 3, "");
1098 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1100 ne = trial_fes.GetNE();
1101 geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
1102 mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
1103 mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
1104 dofs1D = mapsC->ndof;
1105 quad1D = mapsC->nqpt;
1107 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
1109 pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
1111 Vector coeff(ne * nq);
1115 for (int e=0; e<ne; ++e)
1117 ElementTransformation *tr = mesh->GetElementTransformation(e);
1118 for (int p=0; p<nq; ++p)
1120 coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
1125 // Use the same setup functions as VectorFEMassIntegrator.
1126 if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
1128 PADiffusionSetup3D(quad1D, 1, ne, ir->GetWeights(), geom->J,
1131 else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
1133 PADiffusionSetup2D<2>(quad1D, 1, ne, ir->GetWeights(), geom->J,
1138 MFEM_ABORT("Unknown kernel.
");
1142 void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const
1145 PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1146 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1148 PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1149 mapsO->Bt, mapsC->Bt, pa_data, x, y);
1152 MFEM_ABORT("Unsupported dimension!
");
1156 void MixedVectorGradientIntegrator::AddMultTransposePA(const Vector &x,
1160 PAHcurlH1ApplyTranspose3D(dofs1D, quad1D, ne, mapsC->B, mapsO->B,
1161 mapsC->Bt, mapsC->Gt, pa_data, x, y);
1163 PAHcurlH1ApplyTranspose2D(dofs1D, quad1D, ne, mapsC->B, mapsO->B,
1164 mapsC->Bt, mapsC->Gt, pa_data, x, y);
1167 MFEM_ABORT("Unsupported dimension!
");
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 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 PAHcurlL2Setup(const int NQ, const int coeffDim, const int NE, const Array< double > &w, Vector &coeff, Vector &op)
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 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 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)
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
constexpr int HCURL_MAX_D1D
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
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 PAHcurlH1ApplyTranspose2D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &bo, const Array< double > &bct, const Array< double > &gct, const Vector &pa_data, 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)
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 PAHdivMassAssembleDiagonal2D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Vector &op_, Vector &diag_)
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 PAHdivMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Vector &op_, Vector &diag_)
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 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
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void PAHcurlH1ApplyTranspose3D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &bo, const Array< double > &bct, const Array< double > &gct, const Vector &pa_data, const Vector &x, Vector &y)