19template<
int T_D1D = 0, 
int T_Q1D = 0> 
static 
   20void BFLFEvalAssemble2D(
const int nbe, 
const int d, 
const int q,
 
   21                        const int *markers, 
const real_t *
b,
 
   24   const auto F = coeff.Read();
 
   25   const auto M = 
Reshape(markers, nbe);
 
   27   const auto W = 
Reshape(weights, q);
 
   28   const bool const_coeff = coeff.Size() == 1;
 
   34      if (M(e) == 0) { 
return; } 
 
   36      constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
 
   39      for (
int qx = 0; qx < q; ++qx)
 
   41         const real_t coeff_val = const_coeff ? C(0,0) : C(qx,e);
 
   42         QQ[qx] = W(qx) * coeff_val;
 
   44      for (
int dx = 0; dx < d; ++dx)
 
   47         for (
int qx = 0; qx < q; ++qx) { 
u += QQ[qx] * B(qx,dx); }
 
   53template<
int T_D1D = 0, 
int T_Q1D = 0> 
static 
   54void BFLFEvalAssemble3D(
const int nbe, 
const int d, 
const int q,
 
   55                        const int *markers, 
const real_t *
b,
 
   58   const auto F = coeff.Read();
 
   59   const auto M = 
Reshape(markers, nbe);
 
   61   const auto W = 
Reshape(weights, q, q);
 
   62   const bool const_coeff = coeff.Size() == 1;
 
   63   const auto C = const_coeff ? 
Reshape(F,1,1,1) : 
Reshape(F,q,q,nbe);
 
   68      if (M(e) == 0) { 
return; } 
 
   70      constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
 
   71      constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
 
   73      MFEM_SHARED 
real_t sBt[Q*D];
 
   74      MFEM_SHARED 
real_t sQQ[Q*Q];
 
   75      MFEM_SHARED 
real_t sQD[Q*D];
 
   78      kernels::internal::LoadB<D,Q>(d, q, B, sBt);
 
   83      MFEM_FOREACH_THREAD(x,x,q)
 
   85         MFEM_FOREACH_THREAD(y,y,q)
 
   87            const real_t coeff_val = const_coeff ? C(0,0,0) : C(x,y,e);
 
   88            QQ(y,x) = W(x,y) * coeff_val;
 
   92      MFEM_FOREACH_THREAD(qy,y,q)
 
   94         MFEM_FOREACH_THREAD(dx,x,d)
 
   97            for (
int qx = 0; qx < q; ++qx) { 
u += QQ(qy,qx) * Bt(dx,qx); }
 
  102      MFEM_FOREACH_THREAD(dy,y,d)
 
  104         MFEM_FOREACH_THREAD(dx,x,d)
 
  107            for (
int qy = 0; qy < q; ++qy) { 
u += QD(qy,dx) * Bt(dy,qy); }
 
  115static void BFLFEvalAssemble(
const FiniteElementSpace &fes,
 
  116                             const IntegrationRule &ir,
 
  117                             const Array<int> &markers,
 
  121   Mesh &mesh = *fes.GetMesh();
 
  122   const int dim = mesh.Dimension();
 
  123   const FiniteElement &el = *fes.GetBE(0);
 
  125   const int d = maps.ndof, q = maps.nqpt;
 
  126   auto ker = (
dim == 2) ? BFLFEvalAssemble2D<> : BFLFEvalAssemble3D<>;
 
  130      if (d==1 && q==1) { ker=BFLFEvalAssemble2D<1,1>; }
 
  131      if (d==2 && q==2) { ker=BFLFEvalAssemble2D<2,2>; }
 
  132      if (d==3 && q==3) { ker=BFLFEvalAssemble2D<3,3>; }
 
  133      if (d==4 && q==4) { ker=BFLFEvalAssemble2D<4,4>; }
 
  134      if (d==5 && q==5) { ker=BFLFEvalAssemble2D<5,5>; }
 
  135      if (d==2 && q==3) { ker=BFLFEvalAssemble2D<2,3>; }
 
  136      if (d==3 && q==4) { ker=BFLFEvalAssemble2D<3,4>; }
 
  137      if (d==4 && q==5) { ker=BFLFEvalAssemble2D<4,5>; }
 
  138      if (d==5 && q==6) { ker=BFLFEvalAssemble2D<5,6>; }
 
  143      if (d==1 && q==1) { ker=BFLFEvalAssemble3D<1,1>; }
 
  144      if (d==2 && q==2) { ker=BFLFEvalAssemble3D<2,2>; }
 
  145      if (d==3 && q==3) { ker=BFLFEvalAssemble3D<3,3>; }
 
  146      if (d==4 && q==4) { ker=BFLFEvalAssemble3D<4,4>; }
 
  147      if (d==5 && q==5) { ker=BFLFEvalAssemble3D<5,5>; }
 
  148      if (d==2 && q==3) { ker=BFLFEvalAssemble3D<2,3>; }
 
  149      if (d==3 && q==4) { ker=BFLFEvalAssemble3D<3,4>; }
 
  150      if (d==4 && q==5) { ker=BFLFEvalAssemble3D<4,5>; }
 
  151      if (d==5 && q==6) { ker=BFLFEvalAssemble3D<5,6>; }
 
  154   MFEM_VERIFY(ker, 
"No kernel ndof " << d << 
" nqpt " << q);
 
  157   const int *M = markers.Read();
 
  158   const real_t *B = maps.B.Read();
 
  159   const real_t *W = ir.GetWeights().Read();
 
  160   real_t *Y = y.ReadWrite();
 
  161   ker(nbe, d, q, M, B, W, coeff, Y);
 
  170   const int qorder = oa * fe.
GetOrder() + ob;
 
  177   BFLFEvalAssemble(fes, ir, markers, coeff, 
b);
 
 
Class to represent a coefficient evaluated at quadrature points.
@ 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...
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
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
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)
void forall(int N, lambda &&body)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)