12 #include "../general/forall.hpp"
15 #include "libceed/diffusion.hpp"
25 static void PAVectorDiffusionSetup2D(
const int Q1D,
27 const Array<double> &w,
32 const int NQ = Q1D*Q1D;
35 auto J =
Reshape(j.Read(), NQ, 2, 2, NE);
36 auto y =
Reshape(op.Write(), NQ, 3, NE);
40 for (
int q = 0; q < NQ; ++q)
42 const double J11 = J(q,0,0,e);
43 const double J21 = J(q,1,0,e);
44 const double J12 = J(q,0,1,e);
45 const double J22 = J(q,1,1,e);
46 const double c_detJ = W[q] * COEFF / ((J11*J22)-(J21*J12));
47 y(q,0,e) = c_detJ * (J12*J12 + J22*J22);
48 y(q,1,e) = -c_detJ * (J12*J11 + J22*J21);
49 y(q,2,e) = c_detJ * (J11*J11 + J21*J21);
55 static void PAVectorDiffusionSetup3D(
const int Q1D,
57 const Array<double> &w,
62 const int NQ = Q1D*Q1D*Q1D;
64 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
65 auto y =
Reshape(op.Write(), NQ, 6, NE);
68 for (
int q = 0; q < NQ; ++q)
70 const double J11 = J(q,0,0,e);
71 const double J21 = J(q,1,0,e);
72 const double J31 = J(q,2,0,e);
73 const double J12 = J(q,0,1,e);
74 const double J22 = J(q,1,1,e);
75 const double J32 = J(q,2,1,e);
76 const double J13 = J(q,0,2,e);
77 const double J23 = J(q,1,2,e);
78 const double J33 = J(q,2,2,e);
79 const double detJ = J11 * (J22 * J33 - J32 * J23) -
80 J21 * (J12 * J33 - J32 * J13) +
81 J31 * (J12 * J23 - J22 * J13);
82 const double c_detJ = W[q] * COEFF / detJ;
84 const double A11 = (J22 * J33) - (J23 * J32);
85 const double A12 = (J32 * J13) - (J12 * J33);
86 const double A13 = (J12 * J23) - (J22 * J13);
87 const double A21 = (J31 * J23) - (J21 * J33);
88 const double A22 = (J11 * J33) - (J13 * J31);
89 const double A23 = (J21 * J13) - (J11 * J23);
90 const double A31 = (J21 * J32) - (J31 * J22);
91 const double A32 = (J31 * J12) - (J11 * J32);
92 const double A33 = (J11 * J22) - (J12 * J21);
94 y(q,0,e) = c_detJ * (A11*A11 + A12*A12 + A13*A13);
95 y(q,1,e) = c_detJ * (A11*A21 + A12*A22 + A13*A23);
96 y(q,2,e) = c_detJ * (A11*A31 + A12*A32 + A13*A33);
97 y(q,3,e) = c_detJ * (A21*A21 + A22*A22 + A23*A23);
98 y(q,4,e) = c_detJ * (A21*A31 + A22*A32 + A23*A33);
99 y(q,5,e) = c_detJ * (A31*A31 + A32*A32 + A33*A33);
104 static void PAVectorDiffusionSetup(
const int dim,
107 const Array<double> &W,
112 if (!(dim == 2 || dim == 3))
114 MFEM_ABORT(
"Dimension not supported.");
118 PAVectorDiffusionSetup2D(Q1D, NE, W, J, COEFF, op);
122 PAVectorDiffusionSetup3D(Q1D, NE, W, J, COEFF, op);
132 = IntRule ? IntRule : &DiffusionIntegrator::GetRule(el, el);
133 if (DeviceCanUseCeed())
136 ceedDataPtr =
new CeedData;
137 InitCeedCoeff(Q, *mesh, *ir, ceedDataPtr);
138 return CeedPADiffusionAssemble(fes, *ir, * ceedDataPtr);
140 const int dims = el.
GetDim();
141 const int symmDims = (dims * (dims + 1)) / 2;
150 pa_data.SetSize(symmDims * nq * ne, Device::GetDeviceMemoryType());
155 MFEM_VERIFY(cQ != NULL,
"only ConstantCoefficient is supported!");
161 if (dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAVectorDiffusionSetup"); }
162 if (dim == 2 && sdim == 3)
164 constexpr
int DIM = 2;
165 constexpr
int SDIM = 3;
166 const int NQ = quad1D*quad1D;
172 for (
int q = 0; q < NQ; ++q)
174 const double wq = W[q];
175 const double J11 = J(q,0,0,e);
176 const double J21 = J(q,1,0,e);
177 const double J31 = J(q,2,0,e);
178 const double J12 = J(q,0,1,e);
179 const double J22 = J(q,1,1,e);
180 const double J32 = J(q,2,1,e);
181 const double E = J11*J11 + J21*J21 + J31*J31;
182 const double G = J12*J12 + J22*J22 + J32*J32;
183 const double F = J11*J12 + J21*J22 + J31*J32;
184 const double iw = 1.0 / sqrt(E*G - F*F);
185 const double alpha = wq * coeff * iw;
186 D(q,0,e) = alpha * G;
187 D(q,1,e) = -alpha * F;
188 D(q,2,e) = alpha * E;
194 PAVectorDiffusionSetup(dim, quad1D, ne, w, j, coeff, d);
199 template<
int T_D1D = 0,
int T_Q1D = 0,
int T_VDIM = 0>
static
200 void PAVectorDiffusionApply2D(
const int NE,
212 const int D1D = T_D1D ? T_D1D : d1d;
213 const int Q1D = T_Q1D ? T_Q1D : q1d;
214 const int VDIM = T_VDIM ? T_VDIM : vdim;
215 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
216 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
226 const int D1D = T_D1D ? T_D1D : d1d;
227 const int Q1D = T_Q1D ? T_Q1D : q1d;
228 const int VDIM = T_VDIM ? T_VDIM : vdim;
229 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
230 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
232 double grad[max_Q1D][max_Q1D][2];
233 for (
int c = 0; c < VDIM; c++)
235 for (
int qy = 0; qy < Q1D; ++qy)
237 for (
int qx = 0; qx < Q1D; ++qx)
239 grad[qy][qx][0] = 0.0;
240 grad[qy][qx][1] = 0.0;
243 for (
int dy = 0; dy < D1D; ++dy)
245 double gradX[max_Q1D][2];
246 for (
int qx = 0; qx < Q1D; ++qx)
251 for (
int dx = 0; dx < D1D; ++dx)
253 const double s = x(dx,dy,c,e);
254 for (
int qx = 0; qx < Q1D; ++qx)
256 gradX[qx][0] += s * B(qx,dx);
257 gradX[qx][1] += s * G(qx,dx);
260 for (
int qy = 0; qy < Q1D; ++qy)
262 const double wy = B(qy,dy);
263 const double wDy = G(qy,dy);
264 for (
int qx = 0; qx < Q1D; ++qx)
266 grad[qy][qx][0] += gradX[qx][1] * wy;
267 grad[qy][qx][1] += gradX[qx][0] * wDy;
272 for (
int qy = 0; qy < Q1D; ++qy)
274 for (
int qx = 0; qx < Q1D; ++qx)
276 const int q = qx + qy * Q1D;
277 const double O11 = D(q,0,e);
278 const double O12 = D(q,1,e);
279 const double O22 = D(q,2,e);
280 const double gradX = grad[qy][qx][0];
281 const double gradY = grad[qy][qx][1];
282 grad[qy][qx][0] = (O11 * gradX) + (O12 * gradY);
283 grad[qy][qx][1] = (O12 * gradX) + (O22 * gradY);
286 for (
int qy = 0; qy < Q1D; ++qy)
288 double gradX[max_D1D][2];
289 for (
int dx = 0; dx < D1D; ++dx)
294 for (
int qx = 0; qx < Q1D; ++qx)
296 const double gX = grad[qy][qx][0];
297 const double gY = grad[qy][qx][1];
298 for (
int dx = 0; dx < D1D; ++dx)
300 const double wx = Bt(dx,qx);
301 const double wDx = Gt(dx,qx);
302 gradX[dx][0] += gX * wDx;
303 gradX[dx][1] += gY * wx;
306 for (
int dy = 0; dy < D1D; ++dy)
308 const double wy = Bt(dy,qy);
309 const double wDy = Gt(dy,qy);
310 for (
int dx = 0; dx < D1D; ++dx)
312 y(dx,dy,c,e) += ((gradX[dx][0] * wy) + (gradX[dx][1] * wDy));
321 template<
const int T_D1D = 0,
322 const int T_Q1D = 0>
static
323 void PAVectorDiffusionApply3D(
const int NE,
324 const Array<double> &b,
325 const Array<double> &g,
326 const Array<double> &bt,
327 const Array<double> >,
331 int d1d = 0,
int q1d = 0)
333 const int D1D = T_D1D ? T_D1D : d1d;
334 const int Q1D = T_Q1D ? T_Q1D : q1d;
335 constexpr
int VDIM = 3;
336 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
337 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
338 auto B =
Reshape(b.Read(), Q1D, D1D);
339 auto G =
Reshape(g.Read(), Q1D, D1D);
340 auto Bt =
Reshape(bt.Read(), D1D, Q1D);
341 auto Gt =
Reshape(gt.Read(), D1D, Q1D);
342 auto op =
Reshape(_op.Read(), Q1D*Q1D*Q1D, 6, NE);
343 auto x =
Reshape(_x.Read(), D1D, D1D, D1D, VDIM, NE);
344 auto y =
Reshape(_y.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
347 const int D1D = T_D1D ? T_D1D : d1d;
348 const int Q1D = T_Q1D ? T_Q1D : q1d;
349 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
350 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
351 for (
int c = 0; c < VDIM; ++ c)
353 double grad[max_Q1D][max_Q1D][max_Q1D][3];
354 for (
int qz = 0; qz < Q1D; ++qz)
356 for (
int qy = 0; qy < Q1D; ++qy)
358 for (
int qx = 0; qx < Q1D; ++qx)
360 grad[qz][qy][qx][0] = 0.0;
361 grad[qz][qy][qx][1] = 0.0;
362 grad[qz][qy][qx][2] = 0.0;
366 for (
int dz = 0; dz < D1D; ++dz)
368 double gradXY[max_Q1D][max_Q1D][3];
369 for (
int qy = 0; qy < Q1D; ++qy)
371 for (
int qx = 0; qx < Q1D; ++qx)
373 gradXY[qy][qx][0] = 0.0;
374 gradXY[qy][qx][1] = 0.0;
375 gradXY[qy][qx][2] = 0.0;
378 for (
int dy = 0; dy < D1D; ++dy)
380 double gradX[max_Q1D][2];
381 for (
int qx = 0; qx < Q1D; ++qx)
386 for (
int dx = 0; dx < D1D; ++dx)
388 const double s = x(dx,dy,dz,c,e);
389 for (
int qx = 0; qx < Q1D; ++qx)
391 gradX[qx][0] += s * B(qx,dx);
392 gradX[qx][1] += s * G(qx,dx);
395 for (
int qy = 0; qy < Q1D; ++qy)
397 const double wy = B(qy,dy);
398 const double wDy = G(qy,dy);
399 for (
int qx = 0; qx < Q1D; ++qx)
401 const double wx = gradX[qx][0];
402 const double wDx = gradX[qx][1];
403 gradXY[qy][qx][0] += wDx * wy;
404 gradXY[qy][qx][1] += wx * wDy;
405 gradXY[qy][qx][2] += wx * wy;
409 for (
int qz = 0; qz < Q1D; ++qz)
411 const double wz = B(qz,dz);
412 const double wDz = G(qz,dz);
413 for (
int qy = 0; qy < Q1D; ++qy)
415 for (
int qx = 0; qx < Q1D; ++qx)
417 grad[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
418 grad[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
419 grad[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
425 for (
int qz = 0; qz < Q1D; ++qz)
427 for (
int qy = 0; qy < Q1D; ++qy)
429 for (
int qx = 0; qx < Q1D; ++qx)
431 const int q = qx + (qy + qz * Q1D) * Q1D;
432 const double O11 = op(q,0,e);
433 const double O12 = op(q,1,e);
434 const double O13 = op(q,2,e);
435 const double O22 = op(q,3,e);
436 const double O23 = op(q,4,e);
437 const double O33 = op(q,5,e);
438 const double gradX = grad[qz][qy][qx][0];
439 const double gradY = grad[qz][qy][qx][1];
440 const double gradZ = grad[qz][qy][qx][2];
441 grad[qz][qy][qx][0] = (O11*gradX)+(O12*gradY)+(O13*gradZ);
442 grad[qz][qy][qx][1] = (O12*gradX)+(O22*gradY)+(O23*gradZ);
443 grad[qz][qy][qx][2] = (O13*gradX)+(O23*gradY)+(O33*gradZ);
447 for (
int qz = 0; qz < Q1D; ++qz)
449 double gradXY[max_D1D][max_D1D][3];
450 for (
int dy = 0; dy < D1D; ++dy)
452 for (
int dx = 0; dx < D1D; ++dx)
454 gradXY[dy][dx][0] = 0;
455 gradXY[dy][dx][1] = 0;
456 gradXY[dy][dx][2] = 0;
459 for (
int qy = 0; qy < Q1D; ++qy)
461 double gradX[max_D1D][3];
462 for (
int dx = 0; dx < D1D; ++dx)
468 for (
int qx = 0; qx < Q1D; ++qx)
470 const double gX = grad[qz][qy][qx][0];
471 const double gY = grad[qz][qy][qx][1];
472 const double gZ = grad[qz][qy][qx][2];
473 for (
int dx = 0; dx < D1D; ++dx)
475 const double wx = Bt(dx,qx);
476 const double wDx = Gt(dx,qx);
477 gradX[dx][0] += gX * wDx;
478 gradX[dx][1] += gY * wx;
479 gradX[dx][2] += gZ * wx;
482 for (
int dy = 0; dy < D1D; ++dy)
484 const double wy = Bt(dy,qy);
485 const double wDy = Gt(dy,qy);
486 for (
int dx = 0; dx < D1D; ++dx)
488 gradXY[dy][dx][0] += gradX[dx][0] * wy;
489 gradXY[dy][dx][1] += gradX[dx][1] * wDy;
490 gradXY[dy][dx][2] += gradX[dx][2] * wy;
494 for (
int dz = 0; dz < D1D; ++dz)
496 const double wz = Bt(dz,qz);
497 const double wDz = Gt(dz,qz);
498 for (
int dy = 0; dy < D1D; ++dy)
500 for (
int dx = 0; dx < D1D; ++dx)
503 ((gradXY[dy][dx][0] * wz) +
504 (gradXY[dy][dx][1] * wz) +
505 (gradXY[dy][dx][2] * wDz));
515 void VectorDiffusionIntegrator::AddMultPA(
const Vector &x,
Vector &y)
const
517 if (DeviceCanUseCeed())
519 CeedAddMult(ceedDataPtr, x, y);
523 const int D1D = dofs1D;
524 const int Q1D = quad1D;
529 const Vector &D = pa_data;
531 if (dim == 2 && sdim == 3)
533 switch ((dofs1D << 4 ) | quad1D)
535 case 0x22:
return PAVectorDiffusionApply2D<2,2,3>(ne,B,G,Bt,Gt,D,x,y);
536 case 0x33:
return PAVectorDiffusionApply2D<3,3,3>(ne,B,G,Bt,Gt,D,x,y);
537 case 0x44:
return PAVectorDiffusionApply2D<4,4,3>(ne,B,G,Bt,Gt,D,x,y);
538 case 0x55:
return PAVectorDiffusionApply2D<5,5,3>(ne,B,G,Bt,Gt,D,x,y);
540 return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim);
543 if (dim == 2 && sdim == 2)
544 {
return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim); }
546 if (dim == 3 && sdim == 3)
547 {
return PAVectorDiffusionApply3D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D); }
549 MFEM_ABORT(
"Unknown kernel.");
553 template<
int T_D1D = 0,
int T_Q1D = 0>
554 static void PAVectorDiffusionDiagonal2D(
const int NE,
562 const int D1D = T_D1D ? T_D1D : d1d;
563 const int Q1D = T_Q1D ? T_Q1D : q1d;
564 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
565 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
574 const int D1D = T_D1D ? T_D1D : d1d;
575 const int Q1D = T_Q1D ? T_Q1D : q1d;
576 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
577 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
579 double QD0[MQ1][MD1];
580 double QD1[MQ1][MD1];
581 double QD2[MQ1][MD1];
582 for (
int qx = 0; qx < Q1D; ++qx)
584 for (
int dy = 0; dy < D1D; ++dy)
589 for (
int qy = 0; qy < Q1D; ++qy)
591 const int q = qx + qy * Q1D;
592 const double D0 = D(q,0,e);
593 const double D1 = D(q,1,e);
594 const double D2 = D(q,2,e);
595 QD0[qx][dy] += B(qy, dy) * B(qy, dy) * D0;
596 QD1[qx][dy] += B(qy, dy) * G(qy, dy) * D1;
597 QD2[qx][dy] += G(qy, dy) * G(qy, dy) * D2;
601 for (
int dy = 0; dy < D1D; ++dy)
603 for (
int dx = 0; dx < D1D; ++dx)
606 for (
int qx = 0; qx < Q1D; ++qx)
608 temp += G(qx, dx) * G(qx, dx) * QD0[qx][dy];
609 temp += G(qx, dx) * B(qx, dx) * QD1[qx][dy];
610 temp += B(qx, dx) * G(qx, dx) * QD1[qx][dy];
611 temp += B(qx, dx) * B(qx, dx) * QD2[qx][dy];
613 Y(dx,dy,0,e) += temp;
614 Y(dx,dy,1,e) += temp;
620 template<
int T_D1D = 0,
int T_Q1D = 0>
621 static void PAVectorDiffusionDiagonal3D(
const int NE,
622 const Array<double> &b,
623 const Array<double> &g,
629 constexpr
int DIM = 3;
630 const int D1D = T_D1D ? T_D1D : d1d;
631 const int Q1D = T_Q1D ? T_Q1D : q1d;
632 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
633 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
634 MFEM_VERIFY(D1D <= MD1,
"");
635 MFEM_VERIFY(Q1D <= MQ1,
"");
636 auto B =
Reshape(b.Read(), Q1D, D1D);
637 auto G =
Reshape(g.Read(), Q1D, D1D);
638 auto Q =
Reshape(d.Read(), Q1D*Q1D*Q1D, 6, NE);
639 auto Y =
Reshape(y.ReadWrite(), D1D, D1D, D1D, 3, NE);
642 const int D1D = T_D1D ? T_D1D : d1d;
643 const int Q1D = T_Q1D ? T_Q1D : q1d;
644 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
645 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
646 double QQD[MQ1][MQ1][MD1];
647 double QDD[MQ1][MD1][MD1];
648 for (
int i = 0; i <
DIM; ++i)
650 for (
int j = 0; j <
DIM; ++j)
653 for (
int qx = 0; qx < Q1D; ++qx)
655 for (
int qy = 0; qy < Q1D; ++qy)
657 for (
int dz = 0; dz < D1D; ++dz)
659 QQD[qx][qy][dz] = 0.0;
660 for (
int qz = 0; qz < Q1D; ++qz)
662 const int q = qx + (qy + qz * Q1D) * Q1D;
663 const int k = j >= i ?
664 3 - (3-i)*(2-i)/2 + j:
665 3 - (3-j)*(2-j)/2 + i;
666 const double O = Q(q,k,e);
667 const double Bz = B(qz,dz);
668 const double Gz = G(qz,dz);
669 const double L = i==2 ? Gz : Bz;
670 const double R = j==2 ? Gz : Bz;
671 QQD[qx][qy][dz] += L * O * R;
677 for (
int qx = 0; qx < Q1D; ++qx)
679 for (
int dz = 0; dz < D1D; ++dz)
681 for (
int dy = 0; dy < D1D; ++dy)
683 QDD[qx][dy][dz] = 0.0;
684 for (
int qy = 0; qy < Q1D; ++qy)
686 const double By = B(qy,dy);
687 const double Gy = G(qy,dy);
688 const double L = i==1 ? Gy : By;
689 const double R = j==1 ? Gy : By;
690 QDD[qx][dy][dz] += L * QQD[qx][qy][dz] * R;
696 for (
int dz = 0; dz < D1D; ++dz)
698 for (
int dy = 0; dy < D1D; ++dy)
700 for (
int dx = 0; dx < D1D; ++dx)
703 for (
int qx = 0; qx < Q1D; ++qx)
705 const double Bx = B(qx,dx);
706 const double Gx = G(qx,dx);
707 const double L = i==0 ? Gx : Bx;
708 const double R = j==0 ? Gx : Bx;
709 temp += L * QDD[qx][dy][dz] * R;
711 Y(dx, dy, dz, 0, e) += temp;
712 Y(dx, dy, dz, 1, e) += temp;
713 Y(dx, dy, dz, 2, e) += temp;
722 static void PAVectorDiffusionAssembleDiagonal(
const int dim,
726 const Array<double> &B,
727 const Array<double> &G,
733 return PAVectorDiffusionDiagonal2D(NE, B, G, op, y, D1D, Q1D);
737 return PAVectorDiffusionDiagonal3D(NE, B, G, op, y, D1D, Q1D);
739 MFEM_ABORT(
"Dimension not implemented.");
742 void VectorDiffusionIntegrator::AssembleDiagonalPA(
Vector &diag)
744 if (DeviceCanUseCeed())
746 CeedAssembleDiagonal(ceedDataPtr, diag);
750 PAVectorDiffusionAssembleDiagonal(dim,
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
int GetDim() const
Returns the reference space dimension for the finite element.
Class for an integration rule - an Array of IntegrationPoint.
A coefficient that is constant across space and time.
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.
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
int GetNE() const
Returns number of elements in the mesh.
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags)
Return the mesh geometric factors corresponding to the given integration rule.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Mesh * GetMesh() const
Returns the mesh.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
int SpaceDimension() const
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 double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...