12 #include "../general/forall.hpp"
15 #include "ceed/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);
38 const bool const_c = c.Size() == 1;
39 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
45 for (
int q = 0; q < NQ; ++q)
47 const double J11 = J(q,0,0,e);
48 const double J21 = J(q,1,0,e);
49 const double J12 = J(q,0,1,e);
50 const double J22 = J(q,1,1,e);
52 const double C1 = const_c ? C(0,0) : C(q,e);
53 const double c_detJ = W[q] * C1 / ((J11*J22)-(J21*J12));
54 y(q,0,e) = c_detJ * (J12*J12 + J22*J22);
55 y(q,1,e) = -c_detJ * (J12*J11 + J22*J21);
56 y(q,2,e) = c_detJ * (J11*J11 + J21*J21);
62 static void PAVectorDiffusionSetup3D(
const int Q1D,
64 const Array<double> &w,
69 const int NQ = Q1D*Q1D*Q1D;
71 auto J =
Reshape(j.Read(), NQ, 3, 3, NE);
72 auto y =
Reshape(op.Write(), NQ, 6, NE);
74 const bool const_c = c.Size() == 1;
75 const auto C = const_c ?
Reshape(c.Read(), 1,1) :
81 for (
int q = 0; q < NQ; ++q)
83 const double J11 = J(q,0,0,e);
84 const double J21 = J(q,1,0,e);
85 const double J31 = J(q,2,0,e);
86 const double J12 = J(q,0,1,e);
87 const double J22 = J(q,1,1,e);
88 const double J32 = J(q,2,1,e);
89 const double J13 = J(q,0,2,e);
90 const double J23 = J(q,1,2,e);
91 const double J33 = J(q,2,2,e);
92 const double detJ = J11 * (J22 * J33 - J32 * J23) -
93 J21 * (J12 * J33 - J32 * J13) +
94 J31 * (J12 * J23 - J22 * J13);
96 const double C1 = const_c ? C(0,0) : C(q,e);
98 const double c_detJ = W[q] * C1 / detJ;
100 const double A11 = (J22 * J33) - (J23 * J32);
101 const double A12 = (J32 * J13) - (J12 * J33);
102 const double A13 = (J12 * J23) - (J22 * J13);
103 const double A21 = (J31 * J23) - (J21 * J33);
104 const double A22 = (J11 * J33) - (J13 * J31);
105 const double A23 = (J21 * J13) - (J11 * J23);
106 const double A31 = (J21 * J32) - (J31 * J22);
107 const double A32 = (J31 * J12) - (J11 * J32);
108 const double A33 = (J11 * J22) - (J12 * J21);
110 y(q,0,e) = c_detJ * (A11*A11 + A12*A12 + A13*A13);
111 y(q,1,e) = c_detJ * (A11*A21 + A12*A22 + A13*A23);
112 y(q,2,e) = c_detJ * (A11*A31 + A12*A32 + A13*A33);
113 y(q,3,e) = c_detJ * (A21*A21 + A22*A22 + A23*A23);
114 y(q,4,e) = c_detJ * (A21*A31 + A22*A32 + A23*A33);
115 y(q,5,e) = c_detJ * (A31*A31 + A32*A32 + A33*A33);
120 static void PAVectorDiffusionSetup(
const int dim,
123 const Array<double> &W,
128 if (!(dim == 2 || dim == 3))
130 MFEM_ABORT(
"Dimension not supported.");
134 PAVectorDiffusionSetup2D(Q1D, NE, W, J, C, op);
138 PAVectorDiffusionSetup3D(Q1D, NE, W, J, C, op);
148 = IntRule ? IntRule : &DiffusionIntegrator::GetRule(el, el);
149 if (DeviceCanUseCeed())
152 ceedOp =
new ceed::PADiffusionIntegrator(fes, *ir, Q);
155 const int dims = el.
GetDim();
156 const int symmDims = (dims * (dims + 1)) / 2;
165 pa_data.SetSize(symmDims * nq * ne, Device::GetDeviceMemoryType());
167 MFEM_VERIFY(!VQ && !MQ,
168 "Only scalar coefficient supported for partial assembly for VectorDiffusionIntegrator");
178 coeff(0) = cQ->constant;
181 dynamic_cast<QuadratureFunctionCoefficient*>(Q))
184 MFEM_VERIFY(qFun.
Size() == ne*nq,
185 "Incompatible QuadratureFunction dimension \n");
188 "IntegrationRule used within integrator and in"
189 " QuadratureFunction appear to be different");
191 coeff.
MakeRef(const_cast<QuadratureFunction &>(qFun),0);
197 for (
int e = 0; e < ne; ++e)
200 for (
int q = 0; q < nq; ++q)
202 Co(q,e) = Q->Eval(T, ir->
IntPoint(q));
208 const Vector &j = geom->J;
210 if (dim == 1) { MFEM_ABORT(
"dim==1 not supported in PAVectorDiffusionSetup"); }
211 if (dim == 2 && sdim == 3)
213 constexpr
int DIM = 2;
214 constexpr
int SDIM = 3;
215 const int NQ = quad1D*quad1D;
220 const bool const_c = coeff.
Size() == 1;
221 const auto C = const_c ?
Reshape(coeff.
Read(), 1,1) :
226 for (
int q = 0; q < NQ; ++q)
228 const double wq = W[q];
229 const double J11 = J(q,0,0,e);
230 const double J21 = J(q,1,0,e);
231 const double J31 = J(q,2,0,e);
232 const double J12 = J(q,0,1,e);
233 const double J22 = J(q,1,1,e);
234 const double J32 = J(q,2,1,e);
235 const double E = J11*J11 + J21*J21 + J31*J31;
236 const double G = J12*J12 + J22*J22 + J32*J32;
237 const double F = J11*J12 + J21*J22 + J31*J32;
238 const double iw = 1.0 /
sqrt(E*G - F*F);
239 const double C1 = const_c ? C(0,0) : C(q,e);
240 const double alpha = wq * C1 * iw;
241 D(q,0,e) = alpha * G;
242 D(q,1,e) = -alpha * F;
243 D(q,2,e) = alpha * E;
249 PAVectorDiffusionSetup(dim, quad1D, ne, w, j, coeff, d);
254 template<
int T_D1D = 0,
int T_Q1D = 0,
int T_VDIM = 0>
static
255 void PAVectorDiffusionApply2D(
const int NE,
267 const int D1D = T_D1D ? T_D1D : d1d;
268 const int Q1D = T_Q1D ? T_Q1D : q1d;
269 const int VDIM = T_VDIM ? T_VDIM : vdim;
270 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
271 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
281 const int D1D = T_D1D ? T_D1D : d1d;
282 const int Q1D = T_Q1D ? T_Q1D : q1d;
283 const int VDIM = T_VDIM ? T_VDIM : vdim;
284 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
285 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
287 double grad[max_Q1D][max_Q1D][2];
288 for (
int c = 0; c < VDIM; c++)
290 for (
int qy = 0; qy < Q1D; ++qy)
292 for (
int qx = 0; qx < Q1D; ++qx)
294 grad[qy][qx][0] = 0.0;
295 grad[qy][qx][1] = 0.0;
298 for (
int dy = 0; dy < D1D; ++dy)
300 double gradX[max_Q1D][2];
301 for (
int qx = 0; qx < Q1D; ++qx)
306 for (
int dx = 0; dx < D1D; ++dx)
308 const double s = x(dx,dy,c,e);
309 for (
int qx = 0; qx < Q1D; ++qx)
311 gradX[qx][0] += s * B(qx,dx);
312 gradX[qx][1] += s * G(qx,dx);
315 for (
int qy = 0; qy < Q1D; ++qy)
317 const double wy = B(qy,dy);
318 const double wDy = G(qy,dy);
319 for (
int qx = 0; qx < Q1D; ++qx)
321 grad[qy][qx][0] += gradX[qx][1] * wy;
322 grad[qy][qx][1] += gradX[qx][0] * wDy;
327 for (
int qy = 0; qy < Q1D; ++qy)
329 for (
int qx = 0; qx < Q1D; ++qx)
331 const int q = qx + qy * Q1D;
332 const double O11 = D(q,0,e);
333 const double O12 = D(q,1,e);
334 const double O22 = D(q,2,e);
335 const double gradX = grad[qy][qx][0];
336 const double gradY = grad[qy][qx][1];
337 grad[qy][qx][0] = (O11 * gradX) + (O12 * gradY);
338 grad[qy][qx][1] = (O12 * gradX) + (O22 * gradY);
341 for (
int qy = 0; qy < Q1D; ++qy)
343 double gradX[max_D1D][2];
344 for (
int dx = 0; dx < D1D; ++dx)
349 for (
int qx = 0; qx < Q1D; ++qx)
351 const double gX = grad[qy][qx][0];
352 const double gY = grad[qy][qx][1];
353 for (
int dx = 0; dx < D1D; ++dx)
355 const double wx = Bt(dx,qx);
356 const double wDx = Gt(dx,qx);
357 gradX[dx][0] += gX * wDx;
358 gradX[dx][1] += gY * wx;
361 for (
int dy = 0; dy < D1D; ++dy)
363 const double wy = Bt(dy,qy);
364 const double wDy = Gt(dy,qy);
365 for (
int dx = 0; dx < D1D; ++dx)
367 y(dx,dy,c,e) += ((gradX[dx][0] * wy) + (gradX[dx][1] * wDy));
376 template<
const int T_D1D = 0,
377 const int T_Q1D = 0>
static
378 void PAVectorDiffusionApply3D(
const int NE,
379 const Array<double> &b,
380 const Array<double> &g,
381 const Array<double> &bt,
382 const Array<double> >,
386 int d1d = 0,
int q1d = 0)
388 const int D1D = T_D1D ? T_D1D : d1d;
389 const int Q1D = T_Q1D ? T_Q1D : q1d;
390 constexpr
int VDIM = 3;
391 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
392 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
393 auto B =
Reshape(b.Read(), Q1D, D1D);
394 auto G =
Reshape(g.Read(), Q1D, D1D);
395 auto Bt =
Reshape(bt.Read(), D1D, Q1D);
396 auto Gt =
Reshape(gt.Read(), D1D, Q1D);
397 auto op =
Reshape(op_.Read(), Q1D*Q1D*Q1D, 6, NE);
398 auto x =
Reshape(x_.Read(), D1D, D1D, D1D, VDIM, NE);
399 auto y =
Reshape(y_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
402 const int D1D = T_D1D ? T_D1D : d1d;
403 const int Q1D = T_Q1D ? T_Q1D : q1d;
404 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
405 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
406 for (
int c = 0; c < VDIM; ++ c)
408 double grad[max_Q1D][max_Q1D][max_Q1D][3];
409 for (
int qz = 0; qz < Q1D; ++qz)
411 for (
int qy = 0; qy < Q1D; ++qy)
413 for (
int qx = 0; qx < Q1D; ++qx)
415 grad[qz][qy][qx][0] = 0.0;
416 grad[qz][qy][qx][1] = 0.0;
417 grad[qz][qy][qx][2] = 0.0;
421 for (
int dz = 0; dz < D1D; ++dz)
423 double gradXY[max_Q1D][max_Q1D][3];
424 for (
int qy = 0; qy < Q1D; ++qy)
426 for (
int qx = 0; qx < Q1D; ++qx)
428 gradXY[qy][qx][0] = 0.0;
429 gradXY[qy][qx][1] = 0.0;
430 gradXY[qy][qx][2] = 0.0;
433 for (
int dy = 0; dy < D1D; ++dy)
435 double gradX[max_Q1D][2];
436 for (
int qx = 0; qx < Q1D; ++qx)
441 for (
int dx = 0; dx < D1D; ++dx)
443 const double s = x(dx,dy,dz,c,e);
444 for (
int qx = 0; qx < Q1D; ++qx)
446 gradX[qx][0] += s * B(qx,dx);
447 gradX[qx][1] += s * G(qx,dx);
450 for (
int qy = 0; qy < Q1D; ++qy)
452 const double wy = B(qy,dy);
453 const double wDy = G(qy,dy);
454 for (
int qx = 0; qx < Q1D; ++qx)
456 const double wx = gradX[qx][0];
457 const double wDx = gradX[qx][1];
458 gradXY[qy][qx][0] += wDx * wy;
459 gradXY[qy][qx][1] += wx * wDy;
460 gradXY[qy][qx][2] += wx * wy;
464 for (
int qz = 0; qz < Q1D; ++qz)
466 const double wz = B(qz,dz);
467 const double wDz = G(qz,dz);
468 for (
int qy = 0; qy < Q1D; ++qy)
470 for (
int qx = 0; qx < Q1D; ++qx)
472 grad[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
473 grad[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
474 grad[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
480 for (
int qz = 0; qz < Q1D; ++qz)
482 for (
int qy = 0; qy < Q1D; ++qy)
484 for (
int qx = 0; qx < Q1D; ++qx)
486 const int q = qx + (qy + qz * Q1D) * Q1D;
487 const double O11 = op(q,0,e);
488 const double O12 = op(q,1,e);
489 const double O13 = op(q,2,e);
490 const double O22 = op(q,3,e);
491 const double O23 = op(q,4,e);
492 const double O33 = op(q,5,e);
493 const double gradX = grad[qz][qy][qx][0];
494 const double gradY = grad[qz][qy][qx][1];
495 const double gradZ = grad[qz][qy][qx][2];
496 grad[qz][qy][qx][0] = (O11*gradX)+(O12*gradY)+(O13*gradZ);
497 grad[qz][qy][qx][1] = (O12*gradX)+(O22*gradY)+(O23*gradZ);
498 grad[qz][qy][qx][2] = (O13*gradX)+(O23*gradY)+(O33*gradZ);
502 for (
int qz = 0; qz < Q1D; ++qz)
504 double gradXY[max_D1D][max_D1D][3];
505 for (
int dy = 0; dy < D1D; ++dy)
507 for (
int dx = 0; dx < D1D; ++dx)
509 gradXY[dy][dx][0] = 0;
510 gradXY[dy][dx][1] = 0;
511 gradXY[dy][dx][2] = 0;
514 for (
int qy = 0; qy < Q1D; ++qy)
516 double gradX[max_D1D][3];
517 for (
int dx = 0; dx < D1D; ++dx)
523 for (
int qx = 0; qx < Q1D; ++qx)
525 const double gX = grad[qz][qy][qx][0];
526 const double gY = grad[qz][qy][qx][1];
527 const double gZ = grad[qz][qy][qx][2];
528 for (
int dx = 0; dx < D1D; ++dx)
530 const double wx = Bt(dx,qx);
531 const double wDx = Gt(dx,qx);
532 gradX[dx][0] += gX * wDx;
533 gradX[dx][1] += gY * wx;
534 gradX[dx][2] += gZ * wx;
537 for (
int dy = 0; dy < D1D; ++dy)
539 const double wy = Bt(dy,qy);
540 const double wDy = Gt(dy,qy);
541 for (
int dx = 0; dx < D1D; ++dx)
543 gradXY[dy][dx][0] += gradX[dx][0] * wy;
544 gradXY[dy][dx][1] += gradX[dx][1] * wDy;
545 gradXY[dy][dx][2] += gradX[dx][2] * wy;
549 for (
int dz = 0; dz < D1D; ++dz)
551 const double wz = Bt(dz,qz);
552 const double wDz = Gt(dz,qz);
553 for (
int dy = 0; dy < D1D; ++dy)
555 for (
int dx = 0; dx < D1D; ++dx)
558 ((gradXY[dy][dx][0] * wz) +
559 (gradXY[dy][dx][1] * wz) +
560 (gradXY[dy][dx][2] * wDz));
570 void VectorDiffusionIntegrator::AddMultPA(
const Vector &x,
Vector &y)
const
572 if (DeviceCanUseCeed())
574 ceedOp->AddMult(x, y);
578 const int D1D = dofs1D;
579 const int Q1D = quad1D;
584 const Vector &D = pa_data;
586 if (dim == 2 && sdim == 3)
588 switch ((dofs1D << 4 ) | quad1D)
590 case 0x22:
return PAVectorDiffusionApply2D<2,2,3>(ne,B,G,Bt,Gt,D,x,y);
591 case 0x33:
return PAVectorDiffusionApply2D<3,3,3>(ne,B,G,Bt,Gt,D,x,y);
592 case 0x44:
return PAVectorDiffusionApply2D<4,4,3>(ne,B,G,Bt,Gt,D,x,y);
593 case 0x55:
return PAVectorDiffusionApply2D<5,5,3>(ne,B,G,Bt,Gt,D,x,y);
595 return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim);
598 if (dim == 2 && sdim == 2)
599 {
return PAVectorDiffusionApply2D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D,sdim); }
601 if (dim == 3 && sdim == 3)
602 {
return PAVectorDiffusionApply3D(ne,B,G,Bt,Gt,D,x,y,D1D,Q1D); }
604 MFEM_ABORT(
"Unknown kernel.");
608 template<
int T_D1D = 0,
int T_Q1D = 0>
609 static void PAVectorDiffusionDiagonal2D(
const int NE,
617 const int D1D = T_D1D ? T_D1D : d1d;
618 const int Q1D = T_Q1D ? T_Q1D : q1d;
619 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
620 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
629 const int D1D = T_D1D ? T_D1D : d1d;
630 const int Q1D = T_Q1D ? T_Q1D : q1d;
631 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
632 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
634 double QD0[MQ1][MD1];
635 double QD1[MQ1][MD1];
636 double QD2[MQ1][MD1];
637 for (
int qx = 0; qx < Q1D; ++qx)
639 for (
int dy = 0; dy < D1D; ++dy)
644 for (
int qy = 0; qy < Q1D; ++qy)
646 const int q = qx + qy * Q1D;
647 const double D0 = D(q,0,e);
648 const double D1 = D(q,1,e);
649 const double D2 = D(q,2,e);
650 QD0[qx][dy] += B(qy, dy) * B(qy, dy) * D0;
651 QD1[qx][dy] += B(qy, dy) * G(qy, dy) * D1;
652 QD2[qx][dy] += G(qy, dy) * G(qy, dy) * D2;
656 for (
int dy = 0; dy < D1D; ++dy)
658 for (
int dx = 0; dx < D1D; ++dx)
661 for (
int qx = 0; qx < Q1D; ++qx)
663 temp += G(qx, dx) * G(qx, dx) * QD0[qx][dy];
664 temp += G(qx, dx) * B(qx, dx) * QD1[qx][dy];
665 temp += B(qx, dx) * G(qx, dx) * QD1[qx][dy];
666 temp += B(qx, dx) * B(qx, dx) * QD2[qx][dy];
668 Y(dx,dy,0,e) += temp;
669 Y(dx,dy,1,e) += temp;
675 template<
int T_D1D = 0,
int T_Q1D = 0>
676 static void PAVectorDiffusionDiagonal3D(
const int NE,
677 const Array<double> &b,
678 const Array<double> &g,
684 constexpr
int DIM = 3;
685 const int D1D = T_D1D ? T_D1D : d1d;
686 const int Q1D = T_Q1D ? T_Q1D : q1d;
687 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
688 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
689 MFEM_VERIFY(D1D <= MD1,
"");
690 MFEM_VERIFY(Q1D <= MQ1,
"");
691 auto B =
Reshape(b.Read(), Q1D, D1D);
692 auto G =
Reshape(g.Read(), Q1D, D1D);
693 auto Q =
Reshape(d.Read(), Q1D*Q1D*Q1D, 6, NE);
694 auto Y =
Reshape(y.ReadWrite(), D1D, D1D, D1D, 3, NE);
697 const int D1D = T_D1D ? T_D1D : d1d;
698 const int Q1D = T_Q1D ? T_Q1D : q1d;
699 constexpr
int MD1 = T_D1D ? T_D1D :
MAX_D1D;
700 constexpr
int MQ1 = T_Q1D ? T_Q1D :
MAX_Q1D;
701 double QQD[MQ1][MQ1][MD1];
702 double QDD[MQ1][MD1][MD1];
703 for (
int i = 0; i <
DIM; ++i)
705 for (
int j = 0; j <
DIM; ++j)
708 for (
int qx = 0; qx < Q1D; ++qx)
710 for (
int qy = 0; qy < Q1D; ++qy)
712 for (
int dz = 0; dz < D1D; ++dz)
714 QQD[qx][qy][dz] = 0.0;
715 for (
int qz = 0; qz < Q1D; ++qz)
717 const int q = qx + (qy + qz * Q1D) * Q1D;
718 const int k = j >= i ?
719 3 - (3-i)*(2-i)/2 + j:
720 3 - (3-j)*(2-j)/2 + i;
721 const double O = Q(q,k,e);
722 const double Bz = B(qz,dz);
723 const double Gz = G(qz,dz);
724 const double L = i==2 ? Gz : Bz;
725 const double R = j==2 ? Gz : Bz;
726 QQD[qx][qy][dz] += L * O * R;
732 for (
int qx = 0; qx < Q1D; ++qx)
734 for (
int dz = 0; dz < D1D; ++dz)
736 for (
int dy = 0; dy < D1D; ++dy)
738 QDD[qx][dy][dz] = 0.0;
739 for (
int qy = 0; qy < Q1D; ++qy)
741 const double By = B(qy,dy);
742 const double Gy = G(qy,dy);
743 const double L = i==1 ? Gy : By;
744 const double R = j==1 ? Gy : By;
745 QDD[qx][dy][dz] += L * QQD[qx][qy][dz] * R;
751 for (
int dz = 0; dz < D1D; ++dz)
753 for (
int dy = 0; dy < D1D; ++dy)
755 for (
int dx = 0; dx < D1D; ++dx)
758 for (
int qx = 0; qx < Q1D; ++qx)
760 const double Bx = B(qx,dx);
761 const double Gx = G(qx,dx);
762 const double L = i==0 ? Gx : Bx;
763 const double R = j==0 ? Gx : Bx;
764 temp += L * QDD[qx][dy][dz] * R;
766 Y(dx, dy, dz, 0, e) += temp;
767 Y(dx, dy, dz, 1, e) += temp;
768 Y(dx, dy, dz, 2, e) += temp;
777 static void PAVectorDiffusionAssembleDiagonal(
const int dim,
781 const Array<double> &B,
782 const Array<double> &G,
788 return PAVectorDiffusionDiagonal2D(NE, B, G, op, y, D1D, Q1D);
792 return PAVectorDiffusionDiagonal3D(NE, B, G, op, y, D1D, Q1D);
794 MFEM_ABORT(
"Dimension not implemented.");
797 void VectorDiffusionIntegrator::AssembleDiagonalPA(
Vector &diag)
799 if (DeviceCanUseCeed())
801 ceedOp->GetDiagonal(diag);
805 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.
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.
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
void SetSize(int s)
Resize the vector to size s.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
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.
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.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
Mesh * GetMesh() const
Returns the mesh.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
int SpaceDimension() const
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
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...
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Class representing a function through its values (scalar or vector) at quadrature points...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).