12 #include "../general/forall.hpp"
26 static void PAVectorDiffusionSetup2D(
const int Q1D,
28 const Array<double> &w,
33 const int NQ = Q1D*Q1D;
36 auto J =
Reshape(j.Read(), NQ, 2, 2, NE);
37 auto y =
Reshape(op.Write(), NQ, 3, NE);
39 const bool const_c = c.Size() == 1;
40 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
46 for (
int q = 0; q < NQ; ++q)
48 const double J11 = J(q,0,0,e);
49 const double J21 = J(q,1,0,e);
50 const double J12 = J(q,0,1,e);
51 const double J22 = J(q,1,1,e);
53 const double C1 = const_c ? C(0,0) : C(q,e);
54 const double c_detJ = W[q] * C1 / ((J11*J22)-(J21*J12));
55 y(q,0,e) = c_detJ * (J12*J12 + J22*J22);
56 y(q,1,e) = -c_detJ * (J12*J11 + J22*J21);
57 y(q,2,e) = c_detJ * (J11*J11 + J21*J21);
63 static void PAVectorDiffusionSetup3D(
const int Q1D,
65 const Array<double> &w,
70 const int NQ = Q1D*Q1D*Q1D;
72 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
73 auto y =
Reshape(op.Write(), NQ, 6, NE);
75 const bool const_c = c.Size() == 1;
76 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
82 for (
int q = 0; q < NQ; ++q)
84 const double J11 = J(q,0,0,e);
85 const double J21 = J(q,1,0,e);
86 const double J31 = J(q,2,0,e);
87 const double J12 = J(q,0,1,e);
88 const double J22 = J(q,1,1,e);
89 const double J32 = J(q,2,1,e);
90 const double J13 = J(q,0,2,e);
91 const double J23 = J(q,1,2,e);
92 const double J33 = J(q,2,2,e);
93 const double detJ = J11 * (J22 * J33 - J32 * J23) -
94 J21 * (J12 * J33 - J32 * J13) +
95 J31 * (J12 * J23 - J22 * J13);
97 const double C1 = const_c ? C(0,0) : C(q,e);
99 const double c_detJ = W[q] * C1 / detJ;
101 const double A11 = (J22 * J33) - (J23 * J32);
102 const double A12 = (J32 * J13) - (J12 * J33);
103 const double A13 = (J12 * J23) - (J22 * J13);
104 const double A21 = (J31 * J23) - (J21 * J33);
105 const double A22 = (J11 * J33) - (J13 * J31);
106 const double A23 = (J21 * J13) - (J11 * J23);
107 const double A31 = (J21 * J32) - (J31 * J22);
108 const double A32 = (J31 * J12) - (J11 * J32);
109 const double A33 = (J11 * J22) - (J12 * J21);
111 y(q,0,e) = c_detJ * (A11*A11 + A12*A12 + A13*A13);
112 y(q,1,e) = c_detJ * (A11*A21 + A12*A22 + A13*A23);
113 y(q,2,e) = c_detJ * (A11*A31 + A12*A32 + A13*A33);
114 y(q,3,e) = c_detJ * (A21*A21 + A22*A22 + A23*A23);
115 y(q,4,e) = c_detJ * (A21*A31 + A22*A32 + A23*A33);
116 y(q,5,e) = c_detJ * (A31*A31 + A32*A32 + A33*A33);
121 static void PAVectorDiffusionSetup(
const int dim,
124 const Array<double> &W,
129 if (!(dim == 2 || dim == 3))
131 MFEM_ABORT(
"Dimension not supported.");
135 PAVectorDiffusionSetup2D(Q1D, NE, W, J, C, op);
139 PAVectorDiffusionSetup3D(Q1D, NE, W, J, C, op);
165 const int dims = el.
GetDim();
166 const int symmDims = (dims * (dims + 1)) / 2;
175 pa_data.SetSize(symmDims * nq * ne, Device::GetDeviceMemoryType());
177 MFEM_VERIFY(!VQ && !MQ,
178 "Only scalar coefficient supported for partial assembly for VectorDiffusionIntegrator");
184 const Vector &j = geom->J;
186 if (dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAVectorDiffusionSetup"); }
187 if (dim == 2 && sdim == 3)
189 constexpr
int DIM = 2;
190 constexpr
int SDIM = 3;
191 const int NQ = quad1D*quad1D;
196 const bool const_c = coeff.
Size() == 1;
197 const auto C = const_c ?
Reshape(coeff.
Read(), 1,1) :
202 for (
int q = 0; q < NQ; ++q)
204 const double wq = W[q];
205 const double J11 = J(q,0,0,e);
206 const double J21 = J(q,1,0,e);
207 const double J31 = J(q,2,0,e);
208 const double J12 = J(q,0,1,e);
209 const double J22 = J(q,1,1,e);
210 const double J32 = J(q,2,1,e);
211 const double E = J11*J11 + J21*J21 + J31*J31;
212 const double G = J12*J12 + J22*J22 + J32*J32;
213 const double F = J11*J12 + J21*J22 + J31*J32;
214 const double iw = 1.0 / sqrt(E*G - F*F);
215 const double C1 = const_c ? C(0,0) : C(q,e);
216 const double alpha = wq * C1 * iw;
217 D(q,0,e) = alpha * G;
218 D(q,1,e) = -alpha * F;
219 D(q,2,e) = alpha * E;
225 PAVectorDiffusionSetup(dim, quad1D, ne, w, j, coeff, d);
230 template<
int T_D1D = 0,
int T_Q1D = 0,
int T_VDIM = 0>
static
231 void PAVectorDiffusionApply2D(
const int NE,
243 const int D1D = T_D1D ? T_D1D : d1d;
244 const int Q1D = T_Q1D ? T_Q1D : q1d;
245 const int VDIM = T_VDIM ? T_VDIM : vdim;
246 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
247 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
257 const int D1D = T_D1D ? T_D1D : d1d;
258 const int Q1D = T_Q1D ? T_Q1D : q1d;
259 const int VDIM = T_VDIM ? T_VDIM : vdim;
260 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
261 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
263 double grad[max_Q1D][max_Q1D][2];
264 for (
int c = 0; c < VDIM; c++)
266 for (
int qy = 0; qy < Q1D; ++qy)
268 for (
int qx = 0; qx < Q1D; ++qx)
270 grad[qy][qx][0] = 0.0;
271 grad[qy][qx][1] = 0.0;
274 for (
int dy = 0; dy < D1D; ++dy)
276 double gradX[max_Q1D][2];
277 for (
int qx = 0; qx < Q1D; ++qx)
282 for (
int dx = 0; dx < D1D; ++dx)
284 const double s = x(dx,dy,c,e);
285 for (
int qx = 0; qx < Q1D; ++qx)
287 gradX[qx][0] += s * B(qx,dx);
288 gradX[qx][1] += s * G(qx,dx);
291 for (
int qy = 0; qy < Q1D; ++qy)
293 const double wy = B(qy,dy);
294 const double wDy = G(qy,dy);
295 for (
int qx = 0; qx < Q1D; ++qx)
297 grad[qy][qx][0] += gradX[qx][1] * wy;
298 grad[qy][qx][1] += gradX[qx][0] * wDy;
303 for (
int qy = 0; qy < Q1D; ++qy)
305 for (
int qx = 0; qx < Q1D; ++qx)
307 const int q = qx + qy * Q1D;
308 const double O11 = D(q,0,e);
309 const double O12 = D(q,1,e);
310 const double O22 = D(q,2,e);
311 const double gradX = grad[qy][qx][0];
312 const double gradY = grad[qy][qx][1];
313 grad[qy][qx][0] = (O11 * gradX) + (O12 * gradY);
314 grad[qy][qx][1] = (O12 * gradX) + (O22 * gradY);
317 for (
int qy = 0; qy < Q1D; ++qy)
319 double gradX[max_D1D][2];
320 for (
int dx = 0; dx < D1D; ++dx)
325 for (
int qx = 0; qx < Q1D; ++qx)
327 const double gX = grad[qy][qx][0];
328 const double gY = grad[qy][qx][1];
329 for (
int dx = 0; dx < D1D; ++dx)
331 const double wx = Bt(dx,qx);
332 const double wDx = Gt(dx,qx);
333 gradX[dx][0] += gX * wDx;
334 gradX[dx][1] += gY * wx;
337 for (
int dy = 0; dy < D1D; ++dy)
339 const double wy = Bt(dy,qy);
340 const double wDy = Gt(dy,qy);
341 for (
int dx = 0; dx < D1D; ++dx)
343 y(dx,dy,c,e) += ((gradX[dx][0] * wy) + (gradX[dx][1] * wDy));
352 template<
const int T_D1D = 0,
353 const int T_Q1D = 0>
static
354 void PAVectorDiffusionApply3D(
const int NE,
355 const Array<double> &b,
356 const Array<double> &g,
357 const Array<double> &bt,
358 const Array<double> >,
362 int d1d = 0,
int q1d = 0)
364 const int D1D = T_D1D ? T_D1D : d1d;
365 const int Q1D = T_Q1D ? T_Q1D : q1d;
366 constexpr
int VDIM = 3;
367 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
368 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
369 auto B =
Reshape(b.Read(), Q1D, D1D);
370 auto G =
Reshape(g.Read(), Q1D, D1D);
371 auto Bt =
Reshape(bt.Read(), D1D, Q1D);
372 auto Gt =
Reshape(gt.Read(), D1D, Q1D);
373 auto op =
Reshape(op_.Read(), Q1D*Q1D*Q1D, 6, NE);
374 auto x =
Reshape(x_.Read(), D1D, D1D, D1D, VDIM, NE);
375 auto y =
Reshape(y_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
378 const int D1D = T_D1D ? T_D1D : d1d;
379 const int Q1D = T_Q1D ? T_Q1D : q1d;
380 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
381 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
382 for (
int c = 0; c < VDIM; ++ c)
384 double grad[max_Q1D][max_Q1D][max_Q1D][3];
385 for (
int qz = 0; qz < Q1D; ++qz)
387 for (
int qy = 0; qy < Q1D; ++qy)
389 for (
int qx = 0; qx < Q1D; ++qx)
391 grad[qz][qy][qx][0] = 0.0;
392 grad[qz][qy][qx][1] = 0.0;
393 grad[qz][qy][qx][2] = 0.0;
397 for (
int dz = 0; dz < D1D; ++dz)
399 double gradXY[max_Q1D][max_Q1D][3];
400 for (
int qy = 0; qy < Q1D; ++qy)
402 for (
int qx = 0; qx < Q1D; ++qx)
404 gradXY[qy][qx][0] = 0.0;
405 gradXY[qy][qx][1] = 0.0;
406 gradXY[qy][qx][2] = 0.0;
409 for (
int dy = 0; dy < D1D; ++dy)
411 double gradX[max_Q1D][2];
412 for (
int qx = 0; qx < Q1D; ++qx)
417 for (
int dx = 0; dx < D1D; ++dx)
419 const double s = x(dx,dy,dz,c,e);
420 for (
int qx = 0; qx < Q1D; ++qx)
422 gradX[qx][0] += s * B(qx,dx);
423 gradX[qx][1] += s * G(qx,dx);
426 for (
int qy = 0; qy < Q1D; ++qy)
428 const double wy = B(qy,dy);
429 const double wDy = G(qy,dy);
430 for (
int qx = 0; qx < Q1D; ++qx)
432 const double wx = gradX[qx][0];
433 const double wDx = gradX[qx][1];
434 gradXY[qy][qx][0] += wDx * wy;
435 gradXY[qy][qx][1] += wx * wDy;
436 gradXY[qy][qx][2] += wx * wy;
440 for (
int qz = 0; qz < Q1D; ++qz)
442 const double wz = B(qz,dz);
443 const double wDz = G(qz,dz);
444 for (
int qy = 0; qy < Q1D; ++qy)
446 for (
int qx = 0; qx < Q1D; ++qx)
448 grad[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
449 grad[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
450 grad[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
456 for (
int qz = 0; qz < Q1D; ++qz)
458 for (
int qy = 0; qy < Q1D; ++qy)
460 for (
int qx = 0; qx < Q1D; ++qx)
462 const int q = qx + (qy + qz * Q1D) * Q1D;
463 const double O11 = op(q,0,e);
464 const double O12 = op(q,1,e);
465 const double O13 = op(q,2,e);
466 const double O22 = op(q,3,e);
467 const double O23 = op(q,4,e);
468 const double O33 = op(q,5,e);
469 const double gradX = grad[qz][qy][qx][0];
470 const double gradY = grad[qz][qy][qx][1];
471 const double gradZ = grad[qz][qy][qx][2];
472 grad[qz][qy][qx][0] = (O11*gradX)+(O12*gradY)+(O13*gradZ);
473 grad[qz][qy][qx][1] = (O12*gradX)+(O22*gradY)+(O23*gradZ);
474 grad[qz][qy][qx][2] = (O13*gradX)+(O23*gradY)+(O33*gradZ);
478 for (
int qz = 0; qz < Q1D; ++qz)
480 double gradXY[max_D1D][max_D1D][3];
481 for (
int dy = 0; dy < D1D; ++dy)
483 for (
int dx = 0; dx < D1D; ++dx)
485 gradXY[dy][dx][0] = 0;
486 gradXY[dy][dx][1] = 0;
487 gradXY[dy][dx][2] = 0;
490 for (
int qy = 0; qy < Q1D; ++qy)
492 double gradX[max_D1D][3];
493 for (
int dx = 0; dx < D1D; ++dx)
499 for (
int qx = 0; qx < Q1D; ++qx)
501 const double gX = grad[qz][qy][qx][0];
502 const double gY = grad[qz][qy][qx][1];
503 const double gZ = grad[qz][qy][qx][2];
504 for (
int dx = 0; dx < D1D; ++dx)
506 const double wx = Bt(dx,qx);
507 const double wDx = Gt(dx,qx);
508 gradX[dx][0] += gX * wDx;
509 gradX[dx][1] += gY * wx;
510 gradX[dx][2] += gZ * wx;
513 for (
int dy = 0; dy < D1D; ++dy)
515 const double wy = Bt(dy,qy);
516 const double wDy = Gt(dy,qy);
517 for (
int dx = 0; dx < D1D; ++dx)
519 gradXY[dy][dx][0] += gradX[dx][0] * wy;
520 gradXY[dy][dx][1] += gradX[dx][1] * wDy;
521 gradXY[dy][dx][2] += gradX[dx][2] * wy;
525 for (
int dz = 0; dz < D1D; ++dz)
527 const double wz = Bt(dz,qz);
528 const double wDz = Gt(dz,qz);
529 for (
int dy = 0; dy < D1D; ++dy)
531 for (
int dx = 0; dx < D1D; ++dx)
534 ((gradXY[dy][dx][0] * wz) +
535 (gradXY[dy][dx][1] * wz) +
536 (gradXY[dy][dx][2] * wDz));
546 void VectorDiffusionIntegrator::AddMultPA(
const Vector &x,
Vector &y)
const
550 ceedOp->AddMult(x, y);
554 const int D1D = dofs1D;
555 const int Q1D = quad1D;
560 const Vector &D = pa_data;
562 if (dim == 2 && sdim == 3)
564 switch ((dofs1D << 4 ) | quad1D)
566 case 0x22:
return PAVectorDiffusionApply2D<2,2,3>(ne,B,G,Bt,Gt,D,x,y);
567 case 0x33:
return PAVectorDiffusionApply2D<3,3,3>(ne,B,G,Bt,Gt,D,x,y);
568 case 0x44:
return PAVectorDiffusionApply2D<4,4,3>(ne,B,G,Bt,Gt,D,x,y);
569 case 0x55:
return PAVectorDiffusionApply2D<5,5,3>(ne,B,G,Bt,Gt,D,x,y);
571 return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim);
574 if (dim == 2 && sdim == 2)
575 {
return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim); }
577 if (dim == 3 && sdim == 3)
578 {
return PAVectorDiffusionApply3D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D); }
580 MFEM_ABORT(
"Unknown kernel.");
584 template<
int T_D1D = 0,
int T_Q1D = 0>
585 static void PAVectorDiffusionDiagonal2D(
const int NE,
593 const int D1D = T_D1D ? T_D1D : d1d;
594 const int Q1D = T_Q1D ? T_Q1D : q1d;
595 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
596 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
605 const int D1D = T_D1D ? T_D1D : d1d;
606 const int Q1D = T_Q1D ? T_Q1D : q1d;
607 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
608 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
610 double QD0[MQ1][MD1];
611 double QD1[MQ1][MD1];
612 double QD2[MQ1][MD1];
613 for (
int qx = 0; qx < Q1D; ++qx)
615 for (
int dy = 0; dy < D1D; ++dy)
620 for (
int qy = 0; qy < Q1D; ++qy)
622 const int q = qx + qy * Q1D;
623 const double D0 = D(q,0,e);
624 const double D1 = D(q,1,e);
625 const double D2 = D(q,2,e);
626 QD0[qx][dy] += B(qy, dy) * B(qy, dy) * D0;
627 QD1[qx][dy] += B(qy, dy) * G(qy, dy) * D1;
628 QD2[qx][dy] += G(qy, dy) * G(qy, dy) * D2;
632 for (
int dy = 0; dy < D1D; ++dy)
634 for (
int dx = 0; dx < D1D; ++dx)
637 for (
int qx = 0; qx < Q1D; ++qx)
639 temp += G(qx, dx) * G(qx, dx) * QD0[qx][dy];
640 temp += G(qx, dx) * B(qx, dx) * QD1[qx][dy];
641 temp += B(qx, dx) * G(qx, dx) * QD1[qx][dy];
642 temp += B(qx, dx) * B(qx, dx) * QD2[qx][dy];
644 Y(dx,dy,0,e) += temp;
645 Y(dx,dy,1,e) += temp;
651 template<
int T_D1D = 0,
int T_Q1D = 0>
652 static void PAVectorDiffusionDiagonal3D(
const int NE,
653 const Array<double> &b,
654 const Array<double> &g,
660 constexpr
int DIM = 3;
661 const int D1D = T_D1D ? T_D1D : d1d;
662 const int Q1D = T_Q1D ? T_Q1D : q1d;
663 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
664 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
665 MFEM_VERIFY(D1D <= MD1,
"");
666 MFEM_VERIFY(Q1D <= MQ1,
"");
667 auto B =
Reshape(b.Read(), Q1D, D1D);
668 auto G =
Reshape(g.Read(), Q1D, D1D);
669 auto Q =
Reshape(d.Read(), Q1D*Q1D*Q1D, 6, NE);
670 auto Y =
Reshape(y.ReadWrite(), D1D, D1D, D1D, 3, NE);
673 const int D1D = T_D1D ? T_D1D : d1d;
674 const int Q1D = T_Q1D ? T_Q1D : q1d;
675 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
676 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
677 double QQD[MQ1][MQ1][MD1];
678 double QDD[MQ1][MD1][MD1];
679 for (
int i = 0; i <
DIM; ++i)
681 for (
int j = 0; j <
DIM; ++j)
684 for (
int qx = 0; qx < Q1D; ++qx)
686 for (
int qy = 0; qy < Q1D; ++qy)
688 for (
int dz = 0; dz < D1D; ++dz)
690 QQD[qx][qy][dz] = 0.0;
691 for (
int qz = 0; qz < Q1D; ++qz)
693 const int q = qx + (qy + qz * Q1D) * Q1D;
694 const int k = j >= i ?
695 3 - (3-i)*(2-i)/2 + j:
696 3 - (3-j)*(2-j)/2 + i;
697 const double O = Q(q,k,e);
698 const double Bz = B(qz,dz);
699 const double Gz = G(qz,dz);
700 const double L = i==2 ? Gz : Bz;
701 const double R = j==2 ? Gz : Bz;
702 QQD[qx][qy][dz] += L * O * R;
708 for (
int qx = 0; qx < Q1D; ++qx)
710 for (
int dz = 0; dz < D1D; ++dz)
712 for (
int dy = 0; dy < D1D; ++dy)
714 QDD[qx][dy][dz] = 0.0;
715 for (
int qy = 0; qy < Q1D; ++qy)
717 const double By = B(qy,dy);
718 const double Gy = G(qy,dy);
719 const double L = i==1 ? Gy : By;
720 const double R = j==1 ? Gy : By;
721 QDD[qx][dy][dz] += L * QQD[qx][qy][dz] * R;
727 for (
int dz = 0; dz < D1D; ++dz)
729 for (
int dy = 0; dy < D1D; ++dy)
731 for (
int dx = 0; dx < D1D; ++dx)
734 for (
int qx = 0; qx < Q1D; ++qx)
736 const double Bx = B(qx,dx);
737 const double Gx = G(qx,dx);
738 const double L = i==0 ? Gx : Bx;
739 const double R = j==0 ? Gx : Bx;
740 temp += L * QDD[qx][dy][dz] * R;
742 Y(dx, dy, dz, 0, e) += temp;
743 Y(dx, dy, dz, 1, e) += temp;
744 Y(dx, dy, dz, 2, e) += temp;
753 static void PAVectorDiffusionAssembleDiagonal(
const int dim,
757 const Array<double> &B,
758 const Array<double> &G,
764 return PAVectorDiffusionDiagonal2D(NE, B, G, op, y, D1D, Q1D);
768 return PAVectorDiffusionDiagonal3D(NE, B, G, op, y, D1D, Q1D);
770 MFEM_ABORT(
"Dimension not implemented.");
773 void VectorDiffusionIntegrator::AssembleDiagonalPA(
Vector &diag)
777 ceedOp->GetDiagonal(diag);
781 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.
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
int Size() const
Returns the size of the vector.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
Class to represent a coefficient evaluated at quadrature points.
int GetNE() const
Returns number of elements in the mesh.
Mesh * GetMesh() const
Returns the mesh.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
int SpaceDimension() const
Represent a DiffusionIntegrator with AssemblyLevel::Partial using libCEED.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
const IntegrationRule & GetRule(const Integrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Class representing the storage layout of a QuadratureFunction.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.