19template<
int T_D1D = 0, 
int T_Q1D = 0>
 
   20static void BLFEvalAssemble2D(
const int vdim, 
const int nbe, 
const int d,
 
   22                              const bool normals, 
const int *markers, 
const real_t *
b,
 
   24                              const Vector &coeff, 
real_t *y)
 
   26   const auto F = coeff.Read();
 
   27   const auto M = 
Reshape(markers, nbe);
 
   29   const auto detJ = 
Reshape(detj, q, nbe);
 
   30   const auto N = 
Reshape(n, q, 2, nbe);
 
   31   const auto W = 
Reshape(weights, q);
 
   32   const int cvdim = normals ? 2 : 1;
 
   33   const bool cst = coeff.Size() == cvdim;
 
   34   const auto C = cst ? 
Reshape(F,cvdim,1,1) : 
Reshape(F,cvdim,q,nbe);
 
   35   auto Y = 
Reshape(y, d, vdim, nbe);
 
   39      if (M(e) == 0) { 
return; } 
 
   41      constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
 
   44      for (
int c = 0; c < vdim; ++c)
 
   46         for (
int qx = 0; qx < q; ++qx)
 
   51               for (
int cd = 0; cd < 2; ++cd)
 
   53                  const real_t cval = cst ? C(cd,0,0) : C(cd,qx,e);
 
   54                  coeff_val += cval * N(qx, cd, e);
 
   59               coeff_val = cst ? C(0,0,0) : C(0,qx,e);
 
   61            QQ[qx] = W(qx) * coeff_val * detJ(qx,e);
 
   63         for (
int dx = 0; dx < d; ++dx)
 
   66            for (
int qx = 0; qx < q; ++qx) { 
u += QQ[qx] * B(qx,dx); }
 
   73template<
int T_D1D = 0, 
int T_Q1D = 0>
 
   74static void BLFEvalAssemble3D(
const int vdim, 
const int nbe, 
const int d,
 
   76                              const bool normals, 
const int *markers, 
const real_t *
b,
 
   78                              const Vector &coeff, 
real_t *y)
 
   80   const auto F = coeff.Read();
 
   81   const auto M = 
Reshape(markers, nbe);
 
   83   const auto detJ = 
Reshape(detj, q, q, nbe);
 
   84   const auto N = 
Reshape(n, q, q, 3, nbe);
 
   85   const auto W = 
Reshape(weights, q, q);
 
   86   const int cvdim = normals ? 3 : 1;
 
   87   const bool cst = coeff.Size() == cvdim;
 
   88   const auto C = cst ? 
Reshape(F,cvdim,1,1,1) : 
Reshape(F,cvdim,q,q,nbe);
 
   89   auto Y = 
Reshape(y, d, d, vdim, nbe);
 
   93      if (M(e) == 0) { 
return; } 
 
   95      constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
 
   96      constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
 
   98      MFEM_SHARED 
real_t sBt[Q*D];
 
   99      MFEM_SHARED 
real_t sQQ[Q*Q];
 
  100      MFEM_SHARED 
real_t sQD[Q*D];
 
  103      kernels::internal::LoadB<D,Q>(d, q, B, sBt);
 
  108      for (
int c = 0; c < vdim; ++c)
 
  110         MFEM_FOREACH_THREAD(x,x,q)
 
  112            MFEM_FOREACH_THREAD(y,y,q)
 
  117                  for (
int cd = 0; cd < 3; ++cd)
 
  119                     real_t cval = cst ? C(cd,0,0,0) : C(cd,x,y,e);
 
  120                     coeff_val += cval * N(x,y,cd,e);
 
  125                  coeff_val = cst ? C(0,0,0,0) : C(0,x,y,e);
 
  127               QQ(y,x) = W(x,y) * coeff_val * detJ(x,y,e);
 
  131         MFEM_FOREACH_THREAD(qy,y,q)
 
  133            MFEM_FOREACH_THREAD(dx,x,d)
 
  136               for (
int qx = 0; qx < q; ++qx) { 
u += QQ(qy,qx) * Bt(dx,qx); }
 
  141         MFEM_FOREACH_THREAD(dy,y,d)
 
  143            MFEM_FOREACH_THREAD(dx,x,d)
 
  146               for (
int qy = 0; qy < q; ++qy) { 
u += QD(qy,dx) * Bt(dy,qy); }
 
  155static void BLFEvalAssemble(
const FiniteElementSpace &fes,
 
  156                            const IntegrationRule &ir,
 
  157                            const Array<int> &markers,
 
  162   if (fes.GetNBE() == 0) { 
return; }
 
  163   Mesh &mesh = *fes.GetMesh();
 
  164   const int dim = mesh.Dimension();
 
  165   const FiniteElement &el = *fes.GetBE(0);
 
  168   const int d = maps.ndof, q = maps.nqpt;
 
  171   const FaceGeometricFactors *geom = mesh.GetFaceGeometricFactors(
 
  173   auto ker = (
dim == 2) ? BLFEvalAssemble2D<> : BLFEvalAssemble3D<>;
 
  177      if (d==1 && q==1) { ker=BLFEvalAssemble2D<1,1>; }
 
  178      if (d==2 && q==2) { ker=BLFEvalAssemble2D<2,2>; }
 
  179      if (d==3 && q==3) { ker=BLFEvalAssemble2D<3,3>; }
 
  180      if (d==4 && q==4) { ker=BLFEvalAssemble2D<4,4>; }
 
  181      if (d==5 && q==5) { ker=BLFEvalAssemble2D<5,5>; }
 
  182      if (d==2 && q==3) { ker=BLFEvalAssemble2D<2,3>; }
 
  183      if (d==3 && q==4) { ker=BLFEvalAssemble2D<3,4>; }
 
  184      if (d==4 && q==5) { ker=BLFEvalAssemble2D<4,5>; }
 
  185      if (d==5 && q==6) { ker=BLFEvalAssemble2D<5,6>; }
 
  190      if (d==1 && q==1) { ker=BLFEvalAssemble3D<1,1>; }
 
  191      if (d==2 && q==2) { ker=BLFEvalAssemble3D<2,2>; }
 
  192      if (d==3 && q==3) { ker=BLFEvalAssemble3D<3,3>; }
 
  193      if (d==4 && q==4) { ker=BLFEvalAssemble3D<4,4>; }
 
  194      if (d==5 && q==5) { ker=BLFEvalAssemble3D<5,5>; }
 
  195      if (d==2 && q==3) { ker=BLFEvalAssemble3D<2,3>; }
 
  196      if (d==3 && q==4) { ker=BLFEvalAssemble3D<3,4>; }
 
  197      if (d==4 && q==5) { ker=BLFEvalAssemble3D<4,5>; }
 
  198      if (d==5 && q==6) { ker=BLFEvalAssemble3D<5,6>; }
 
  201   MFEM_VERIFY(ker, 
"No kernel ndof " << d << 
" nqpt " << q);
 
  203   const int vdim = fes.GetVDim();
 
  205   const int *M = markers.Read();
 
  206   const real_t *B = maps.B.Read();
 
  207   const real_t *detJ = geom->detJ.Read();
 
  208   const real_t *n = geom->normal.Read();
 
  209   const real_t *W = ir.GetWeights().Read();
 
  210   real_t *Y = y.ReadWrite();
 
  211   ker(vdim, nbe, d, q, normals, M, B, detJ, n, W, coeff, Y);
 
  218   if (fes.
GetNBE() == 0) { 
return; }
 
  220   const int qorder = oa * fe.
GetOrder() + ob;
 
  227   BLFEvalAssemble(fes, ir, markers, coeff, 
false, 
b);
 
 
  234   if (fes.
GetNBE() == 0) { 
return; }
 
  236   const int qorder = oa * fe.
GetOrder() + ob;
 
  243   BLFEvalAssemble(fes, ir, markers, coeff, 
true, 
b);
 
 
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
Class to represent a coefficient evaluated at quadrature points.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Class representing the storage layout of a FaceQuadratureFunction.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
int GetNBE() const
Returns number of boundary elements in the mesh.
Mesh * GetMesh() const
Returns the mesh.
Abstract class for all finite elements.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
const IntegrationRule * IntRule
real_t u(const Vector &xvec)
DeviceTensor< 2, real_t > DeviceMatrix
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
@ COMPRESSED
Enable all above compressions.
void forall_2D(int N, int X, int Y, lambda &&body)
MemoryType
Memory types supported by MFEM.
void forall(int N, lambda &&body)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)