13 #include "../general/annotation.hpp" 14 #include "../general/forall.hpp" 15 #include "../linalg/dtensor.hpp" 16 #include "../linalg/kernels.hpp" 22 static void GetSigns(
const FiniteElementSpace &fes,
const FaceType type,
25 const Mesh &mesh = *fes.GetMesh();
26 const int dim = mesh.SpaceDimension();
29 for (
int f = 0;
f < mesh.GetNumFacesWithGhost(); ++
f)
31 Mesh::FaceInformation face = mesh.GetFaceInformation(
f);
32 face_id = face.element[0].local_face_id;
33 if (face.IsNonconformingCoarse())
39 else if ( face.IsOfFaceType(type) )
43 if (face_id==2 || face_id==3)
54 if (face_id==0 || face_id==3 || face_id==4)
71 : type(type_), nf(fes.GetNFbyType(type)), signs(nf),
79 GetSigns(*
fespace, type, signs);
85 MFEM_VERIFY(sfe != NULL,
"Only scalar finite elements are supported");
86 MFEM_VERIFY(tfe != NULL &&
89 "Only Gauss-Lobatto and Bernstein basis are supported in " 90 "FaceQuadratureInterpolator.");
93 template<const
int T_VDIM, const
int T_ND1D, const
int T_NQ1D>
105 const int eval_flags)
107 const int nd1d = maps.
ndof;
108 const int nq1d = maps.
nqpt;
109 const int ND1D = T_ND1D ? T_ND1D : nd1d;
110 const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
111 const int VDIM = T_VDIM ? T_VDIM : vdim;
114 MFEM_VERIFY(VDIM == 2 || !(eval_flags &
DETERMINANTS),
"");
118 auto sign = signs.
Read();
128 "Derivatives on the faces are not yet supported.");
132 const int ND1D = T_ND1D ? T_ND1D : nd1d;
133 const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
134 const int VDIM = T_VDIM ? T_VDIM : vdim;
135 constexpr
int max_ND1D = T_ND1D ? T_ND1D :
MAX_ND1D;
136 constexpr
int max_VDIM = T_VDIM ? T_VDIM :
MAX_VDIM2D;
137 double r_F[max_ND1D][max_VDIM];
138 for (
int d = 0; d < ND1D; d++)
140 for (
int c = 0; c < VDIM; c++)
142 r_F[d][c] = F(d,c,
f);
145 for (
int q = 0; q < NQ1D; ++q)
150 for (
int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
151 for (
int d = 0; d < ND1D; ++d)
153 const double b = B(q,d);
154 for (
int c = 0; c < VDIM; c++) { ed[c] +=
b*r_F[d][c]; }
156 for (
int c = 0; c < VDIM; c++)
167 for (
int i = 0; i < VDIM; i++) { D[i] = 0.0; }
168 for (
int d = 0; d < ND1D; ++d)
170 const double w = G(q,d);
171 for (
int c = 0; c < VDIM; c++)
173 double s_e = r_F[d][c];
181 const double norm = sqrt(D[0]*D[0]+D[1]*D[1]);
188 const double s = sign[
f] ? -1.0 : 1.0;
191 n(0,q,
f) =
s*D[1]/norm;
192 n(1,q,
f) = -
s*D[0]/norm;
196 n(q,0,
f) =
s*D[1]/norm;
197 n(q,1,
f) = -
s*D[0]/norm;
206 template<const
int T_VDIM, const
int T_ND1D, const
int T_NQ1D>
218 const int eval_flags)
220 const int nd1d = maps.
ndof;
221 const int nq1d = maps.
nqpt;
222 const int ND1D = T_ND1D ? T_ND1D : nd1d;
223 const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
224 const int VDIM = T_VDIM ? T_VDIM : vdim;
227 MFEM_VERIFY(VDIM == 3 || !(eval_flags &
DETERMINANTS),
"");
230 auto F =
Reshape(e_vec.
Read(), ND1D, ND1D, VDIM, NF);
231 auto sign = signs.
Read();
241 "Derivatives on the faces are not yet supported.");
244 constexpr
int max_ND1D = T_ND1D ? T_ND1D :
MAX_ND1D;
245 constexpr
int max_NQ1D = T_NQ1D ? T_NQ1D :
MAX_NQ1D;
246 constexpr
int max_VDIM = T_VDIM ? T_VDIM :
MAX_VDIM3D;
247 double r_F[max_ND1D][max_ND1D][max_VDIM];
248 for (
int d1 = 0; d1 < ND1D; d1++)
250 for (
int d2 = 0; d2 < ND1D; d2++)
252 for (
int c = 0; c < VDIM; c++)
254 r_F[d1][d2][c] = F(d1,d2,c,
f);
260 double Bu[max_NQ1D][max_ND1D][max_VDIM];
261 for (
int d2 = 0; d2 < ND1D; ++d2)
263 for (
int q = 0; q < NQ1D; ++q)
265 for (
int c = 0; c < VDIM; c++) { Bu[q][d2][c] = 0.0; }
266 for (
int d1 = 0; d1 < ND1D; ++d1)
268 const double b = B(q,d1);
269 for (
int c = 0; c < VDIM; c++)
271 Bu[q][d2][c] +=
b*r_F[d1][d2][c];
276 double BBu[max_NQ1D][max_NQ1D][max_VDIM];
277 for (
int q2 = 0; q2 < NQ1D; ++q2)
279 for (
int q1 = 0; q1 < NQ1D; ++q1)
281 for (
int c = 0; c < VDIM; c++) { BBu[q2][q1][c] = 0.0; }
282 for (
int d2 = 0; d2 < ND1D; ++d2)
284 const double b = B(q2,d2);
285 for (
int c = 0; c < VDIM; c++)
287 BBu[q2][q1][c] +=
b*Bu[q1][d2][c];
290 for (
int c = 0; c < VDIM; c++)
292 const double v = BBu[q2][q1][c];
304 double Gu[max_NQ1D][max_ND1D][max_VDIM];
305 double Bu[max_NQ1D][max_ND1D][max_VDIM];
306 for (
int d2 = 0; d2 < ND1D; ++d2)
308 for (
int q = 0; q < NQ1D; ++q)
310 for (
int c = 0; c < VDIM; c++)
315 for (
int d1 = 0; d1 < ND1D; ++d1)
317 const double b = B(q,d1);
318 const double g = G(q,d1);
319 for (
int c = 0; c < VDIM; c++)
321 const double u = r_F[d1][d2][c];
328 double BGu[max_NQ1D][max_NQ1D][max_VDIM];
329 double GBu[max_NQ1D][max_NQ1D][max_VDIM];
330 for (
int q2 = 0; q2 < NQ1D; ++q2)
332 for (
int q1 = 0; q1 < NQ1D; ++q1)
334 for (
int c = 0; c < VDIM; c++)
336 BGu[q2][q1][c] = 0.0;
337 GBu[q2][q1][c] = 0.0;
339 for (
int d2 = 0; d2 < ND1D; ++d2)
341 const double b = B(q2,d2);
342 const double g = G(q2,d2);
343 for (
int c = 0; c < VDIM; c++)
345 BGu[q2][q1][c] +=
b*Gu[q1][d2][c];
346 GBu[q2][q1][c] += g*Bu[q1][d2][c];
351 if (VDIM == 3 && ((eval_flags &
NORMALS) ||
355 for (
int q2 = 0; q2 < NQ1D; ++q2)
357 for (
int q1 = 0; q1 < NQ1D; ++q1)
359 const double s = sign[
f] ? -1.0 : 1.0;
360 n[0] =
s*( BGu[q2][q1][1]*GBu[q2][q1][2]-GBu[q2][q1][1]*
362 n[1] =
s*(-BGu[q2][q1][0]*GBu[q2][q1][2]+GBu[q2][q1][0]*
364 n[2] =
s*( BGu[q2][q1][0]*GBu[q2][q1][1]-GBu[q2][q1][0]*
366 const double norm = sqrt(n[0]*n[0]+n[1]*n[1]+n[2]*n[2]);
372 nor(0,q1,q2,
f) = n[0]/norm;
373 nor(1,q1,q2,
f) = n[1]/norm;
374 nor(2,q1,q2,
f) = n[2]/norm;
378 nor(q1,q2,0,
f) = n[0]/norm;
379 nor(q1,q2,1,
f) = n[1]/norm;
380 nor(q1,q2,2,
f) = n[2]/norm;
390 template<const
int T_VDIM, const
int T_ND1D, const
int T_NQ1D>
402 const int eval_flags)
404 MFEM_PERF_SCOPE(
"FaceQuadInterpolator::SmemEval3D");
405 const int nd1d = maps.
ndof;
406 const int nq1d = maps.
nqpt;
407 const int ND1D = T_ND1D ? T_ND1D : nd1d;
408 const int NQ1D = T_NQ1D ? T_NQ1D : nq1d;
409 const int VDIM = T_VDIM ? T_VDIM : vdim;
412 MFEM_VERIFY(VDIM == 3 || !(eval_flags &
DETERMINANTS),
"");
415 auto F =
Reshape(e_vec.
Read(), ND1D, ND1D, VDIM, NF);
416 auto sign = signs.
Read();
426 "Derivatives on the faces are not yet supported.");
428 MFEM_FORALL_3D(
f, NF, NQ1D, NQ1D, VDIM,
430 constexpr
int max_ND1D = T_ND1D ? T_ND1D :
MAX_ND1D;
431 constexpr
int max_NQ1D = T_NQ1D ? T_NQ1D :
MAX_NQ1D;
432 constexpr
int max_VDIM = T_VDIM ? T_VDIM :
MAX_VDIM3D;
434 MFEM_SHARED
double sm1[max_NQ1D*max_NQ1D*max_VDIM];
435 MFEM_SHARED
double sm2[max_NQ1D*max_ND1D*max_VDIM];
437 auto s_F = (double(*)[max_ND1D][max_VDIM])sm1;
438 MFEM_FOREACH_THREAD(d1,x,ND1D)
440 MFEM_FOREACH_THREAD(d2,y,ND1D)
442 MFEM_FOREACH_THREAD(c,z,VDIM)
444 s_F[d1][d2][c] = F(d1,d2,c,
f);
452 auto Bu = (double (*)[max_ND1D][max_VDIM])sm2;
453 MFEM_FOREACH_THREAD(d2,x,ND1D)
455 MFEM_FOREACH_THREAD(q1,y,NQ1D)
457 MFEM_FOREACH_THREAD(c,z,VDIM)
460 for (
int d1 = 0; d1 < ND1D; ++d1)
462 thrdBu += B(q1,d1)*s_F[d1][d2][c];
464 Bu[q1][d2][c] = thrdBu;
470 MFEM_FOREACH_THREAD(q2,x,NQ1D)
472 MFEM_FOREACH_THREAD(q1,y,NQ1D)
474 MFEM_FOREACH_THREAD(c,z,VDIM)
477 for (
int d2 = 0; d2 < ND1D; ++d2)
479 v += B(q2,d2)*Bu[q1][d2][c];
493 auto Gu = (double (*)[max_ND1D][max_VDIM])sm2;
494 MFEM_SHARED
double Bu[max_NQ1D][max_ND1D][max_VDIM];
495 MFEM_FOREACH_THREAD(d2,x,ND1D)
497 MFEM_FOREACH_THREAD(q1,y,NQ1D)
499 MFEM_FOREACH_THREAD(c,z,VDIM)
503 for (
int d1 = 0; d1 < ND1D; ++d1)
505 const double u = s_F[d1][d2][c];
506 thrdBu += B(q1,d1)*
u;
507 thrdGu += G(q1,d1)*
u;
509 Gu[q1][d2][c] = thrdGu;
510 Bu[q1][d2][c] = thrdBu;
516 auto BGu = (double (*)[max_NQ1D][max_VDIM])sm1;
517 MFEM_SHARED
double GBu[max_NQ1D][max_NQ1D][max_VDIM];
518 MFEM_FOREACH_THREAD(q2,x,NQ1D)
520 MFEM_FOREACH_THREAD(q1,y,NQ1D)
522 MFEM_FOREACH_THREAD(c,z,VDIM)
524 double thrdBGu = 0.0;
525 double thrdGBu = 0.0;
526 for (
int d2 = 0; d2 < ND1D; ++d2)
528 thrdBGu += B(q2,d2)*Gu[q1][d2][c];
529 thrdGBu += G(q2,d2)*Bu[q1][d2][c];
531 BGu[q2][q1][c] = thrdBGu;
532 GBu[q2][q1][c] = thrdGBu;
538 if (VDIM == 3 && ((eval_flags &
NORMALS) ||
542 MFEM_FOREACH_THREAD(q2,x,NQ1D)
544 MFEM_FOREACH_THREAD(q1,y,NQ1D)
546 if (MFEM_THREAD_ID(z) == 0)
548 const double s = sign[
f] ? -1.0 : 1.0;
549 n[0] =
s*( BGu[q2][q1][1]*GBu[q2][q1][2]-GBu[q2][q1][1]*
551 n[1] =
s*(-BGu[q2][q1][0]*GBu[q2][q1][2]+GBu[q2][q1][0]*
553 n[2] =
s*( BGu[q2][q1][0]*GBu[q2][q1][1]-GBu[q2][q1][0]*
556 const double norm = sqrt(n[0]*n[0]+n[1]*n[1]+n[2]*n[2]);
564 nor(0,q1,q2,
f) = n[0]/norm;
565 nor(1,q1,q2,
f) = n[1]/norm;
566 nor(2,q1,q2,
f) = n[2]/norm;
570 nor(q1,q2,0,
f) = n[0]/norm;
571 nor(q1,q2,1,
f) = n[1]/norm;
572 nor(q1,q2,2,
f) = n[2]/norm;
584 const Vector &e_vec,
unsigned eval_flags,
587 if (nf == 0) {
return; }
594 const int nd1d = maps.
ndof;
595 const int nq1d = maps.
nqpt;
607 const int eval_flags) = NULL;
612 switch (10*nd1d + nq1d)
615 case 11: eval_func = &Eval2D<1,1,1>;
break;
616 case 12: eval_func = &Eval2D<1,1,2>;
break;
618 case 22: eval_func = &Eval2D<1,2,2>;
break;
619 case 23: eval_func = &Eval2D<1,2,3>;
break;
621 case 33: eval_func = &Eval2D<1,3,3>;
break;
622 case 34: eval_func = &Eval2D<1,3,4>;
break;
624 case 44: eval_func = &Eval2D<1,4,4>;
break;
625 case 45: eval_func = &Eval2D<1,4,5>;
break;
626 case 46: eval_func = &Eval2D<1,4,6>;
break;
628 case 55: eval_func = &Eval2D<1,5,5>;
break;
629 case 56: eval_func = &Eval2D<1,5,6>;
break;
630 case 57: eval_func = &Eval2D<1,5,7>;
break;
631 case 58: eval_func = &Eval2D<1,5,8>;
break;
633 if (nq1d >= 10 || !eval_func)
635 eval_func = &Eval2D<1>;
640 switch (10*nd1d + nq1d)
643 case 11: eval_func = &SmemEval3D<1,1,1>;
break;
644 case 12: eval_func = &SmemEval3D<1,1,2>;
break;
646 case 22: eval_func = &SmemEval3D<1,2,2>;
break;
647 case 23: eval_func = &SmemEval3D<1,2,3>;
break;
648 case 24: eval_func = &SmemEval3D<1,2,4>;
break;
650 case 33: eval_func = &SmemEval3D<1,3,3>;
break;
651 case 34: eval_func = &SmemEval3D<1,3,4>;
break;
653 case 44: eval_func = &SmemEval3D<1,4,4>;
break;
654 case 45: eval_func = &SmemEval3D<1,4,5>;
break;
655 case 46: eval_func = &SmemEval3D<1,4,6>;
break;
657 case 55: eval_func = &SmemEval3D<1,5,5>;
break;
658 case 56: eval_func = &SmemEval3D<1,5,6>;
break;
660 if (nq1d >= 10 || !eval_func)
662 eval_func = &Eval3D<1>;
666 else if (vdim ==
dim)
670 switch (10*nd1d + nq1d)
673 case 22: eval_func = &Eval2D<2,2,2>;
break;
674 case 23: eval_func = &Eval2D<2,2,3>;
break;
676 case 33: eval_func = &Eval2D<2,3,3>;
break;
677 case 34: eval_func = &Eval2D<2,3,4>;
break;
679 case 44: eval_func = &Eval2D<2,4,4>;
break;
680 case 45: eval_func = &Eval2D<2,4,5>;
break;
681 case 46: eval_func = &Eval2D<2,4,6>;
break;
683 case 55: eval_func = &Eval2D<2,5,5>;
break;
684 case 56: eval_func = &Eval2D<2,5,6>;
break;
685 case 57: eval_func = &Eval2D<2,5,7>;
break;
686 case 58: eval_func = &Eval2D<2,5,8>;
break;
688 if (nq1d >= 10 || !eval_func)
690 eval_func = &Eval2D<2>;
695 switch (10*nd1d + nq1d)
698 case 22: eval_func = &SmemEval3D<3,2,2>;
break;
699 case 23: eval_func = &SmemEval3D<3,2,3>;
break;
700 case 24: eval_func = &SmemEval3D<3,2,4>;
break;
702 case 33: eval_func = &SmemEval3D<3,3,3>;
break;
703 case 34: eval_func = &SmemEval3D<3,3,4>;
break;
705 case 44: eval_func = &SmemEval3D<3,4,4>;
break;
706 case 45: eval_func = &SmemEval3D<3,4,5>;
break;
707 case 46: eval_func = &SmemEval3D<3,4,6>;
break;
709 case 55: eval_func = &SmemEval3D<3,5,5>;
break;
710 case 56: eval_func = &SmemEval3D<3,5,6>;
break;
712 if (nq1d >= 10 || !eval_func)
714 eval_func = &Eval3D<3>;
720 eval_func(nf, vdim,
q_layout, maps, signs, e_vec,
721 q_val, q_der, q_det, q_nor, eval_flags);
725 MFEM_ABORT(
"case not supported yet");
732 Vector q_der, q_det, q_nor;
733 Mult(e_vec,
VALUES, q_val, q_der, q_det, q_nor);
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Abstract class for all finite elements.
static void Eval2D(const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags)
Template compute kernel for 2D.
Class for an integration rule - an Array of IntegrationPoint.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
FaceQuadratureInterpolator(const FiniteElementSpace &fes, const IntegrationRule &ir, FaceType type)
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
const IntegrationRule * IntRule
Not owned.
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...
Class for finite elements with basis functions that return scalar values.
double f(const Vector &xvec)
static const int MAX_ND1D
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void Mult(const Vector &e_vec, unsigned eval_flags, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor) const
Interpolate the E-vector e_vec to quadrature points.
int GetNE() const
Returns number of elements in the mesh.
int GetVDim() const
Returns vector dimension.
Mesh * GetMesh() const
Returns the mesh.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
static void SmemEval3D(const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags)
void Values(const Vector &e_vec, Vector &q_val) const
Interpolate the values of the E-vector e_vec at quadrature points.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
static const int MAX_VDIM3D
Evaluate the values at quadrature points.
Array< double > B
Basis functions evaluated at quadrature points.
QVectorLayout
Type describing possible layouts for Q-vectors.
static void Eval3D(const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags)
Template compute kernel for 3D.
static const int MAX_NQ1D
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
Evaluate the derivatives at quadrature points.
void f_vec(const Vector &xvec, Vector &f)
double u(const Vector &xvec)
Assuming the derivative at quadrature points form a matrix, this flag can be used to compute and stor...
const FiniteElementSpace * fespace
Not owned.
static const int MAX_VDIM2D
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
QVectorLayout q_layout
Output Q-vector layout.
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)