19template<
int T_D1D = 0,
int T_Q1D = 0>
20static void DLFEvalAssemble2D(
const int vdim,
const int ne,
const int d,
22 const int map_type,
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, ne);
29 const auto DETJ =
Reshape(detj, q, q, ne);
30 const auto W =
Reshape(weights, q, q);
31 const bool cst = coeff.Size() == vdim;
32 const auto C = cst ?
Reshape(F,vdim,1,1,1) :
Reshape(F,vdim,q,q,ne);
33 auto Y =
Reshape(y, d,d, vdim, ne);
37 if (M(e) == 0) {
return; }
39 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
40 constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
42 MFEM_SHARED
real_t sBt[Q*D];
43 MFEM_SHARED
real_t sQQ[Q*Q];
44 MFEM_SHARED
real_t sQD[Q*D];
47 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
52 for (
int c = 0; c < vdim; ++c)
54 const real_t cst_val = C(c,0,0,0);
55 MFEM_FOREACH_THREAD(x,x,q)
57 MFEM_FOREACH_THREAD(y,y,q)
60 const real_t coeff_val = cst ? cst_val : C(c,x,y,e);
61 QQ(y,x) = W(x,y) * coeff_val * detJ;
65 MFEM_FOREACH_THREAD(qy,y,q)
67 MFEM_FOREACH_THREAD(dx,x,d)
70 for (
int qx = 0; qx < q; ++qx) {
u += QQ(qy,qx) * Bt(dx,qx); }
75 MFEM_FOREACH_THREAD(dy,y,d)
77 MFEM_FOREACH_THREAD(dx,x,d)
80 for (
int qy = 0; qy < q; ++qy) {
u += QD(qy,dx) * Bt(dy,qy); }
89template<
int T_D1D = 0,
int T_Q1D = 0>
90static void DLFEvalAssemble3D(
const int vdim,
const int ne,
const int d,
92 const int map_type,
const int *markers,
const real_t *
b,
94 const Vector &coeff,
real_t *y)
96 const auto F = coeff.Read();
97 const auto M =
Reshape(markers, ne);
99 const auto DETJ =
Reshape(detj, q, q, q, ne);
100 const auto W =
Reshape(weights, q,q,q);
101 const bool cst_coeff = coeff.Size() == vdim;
102 const auto C = cst_coeff ?
Reshape(F,vdim,1,1,1,1):
Reshape(F,vdim,q,q,q,ne);
104 auto Y =
Reshape(y, d,d,d, vdim, ne);
108 if (M(e) == 0) {
return; }
110 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
111 constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
112 constexpr int MQD = (Q >= D) ? Q : D;
116 MFEM_SHARED
real_t sBt[Q*D];
118 kernels::internal::LoadB<D,Q>(d,q,B,sBt);
120 MFEM_SHARED
real_t sQQQ[MQD*MQD*MQD];
123 for (
int c = 0; c < vdim; ++c)
125 const real_t cst_val = C(c,0,0,0,0);
126 MFEM_FOREACH_THREAD(x,x,q)
128 MFEM_FOREACH_THREAD(y,y,q)
130 for (
int z = 0; z < q; ++z)
133 const real_t coeff_val = cst_coeff ? cst_val : C(c,x,y,z,e);
134 QQQ(z,y,x) = W(x,y,z) * coeff_val * detJ;
139 MFEM_FOREACH_THREAD(qx,x,q)
141 MFEM_FOREACH_THREAD(qy,y,q)
143 for (
int dz = 0; dz < d; ++dz) {
u[dz] = 0.0; }
144 for (
int qz = 0; qz < q; ++qz)
146 const real_t ZYX = QQQ(qz,qy,qx);
147 for (
int dz = 0; dz < d; ++dz) {
u[dz] += ZYX * Bt(dz,qz); }
149 for (
int dz = 0; dz < d; ++dz) { QQQ(dz,qy,qx) =
u[dz]; }
153 MFEM_FOREACH_THREAD(dz,y,d)
155 MFEM_FOREACH_THREAD(qx,x,q)
157 for (
int dy = 0; dy < d; ++dy) {
u[dy] = 0.0; }
158 for (
int qy = 0; qy < q; ++qy)
160 const real_t zYX = QQQ(dz,qy,qx);
161 for (
int dy = 0; dy < d; ++dy) {
u[dy] += zYX * Bt(dy,qy); }
163 for (
int dy = 0; dy < d; ++dy) { QQQ(dz,dy,qx) =
u[dy]; }
167 MFEM_FOREACH_THREAD(dz,y,d)
169 MFEM_FOREACH_THREAD(dy,x,d)
171 for (
int dx = 0; dx < d; ++dx) {
u[dx] = 0.0; }
172 for (
int qx = 0; qx < q; ++qx)
174 const real_t zyX = QQQ(dz,dy,qx);
175 for (
int dx = 0; dx < d; ++dx) {
u[dx] += zyX * Bt(dx,qx); }
177 for (
int dx = 0; dx < d; ++dx) { Y(dx,dy,dz,c,e) +=
u[dx]; }
185static void DLFEvalAssemble(
const FiniteElementSpace &fes,
186 const IntegrationRule *ir,
187 const Array<int> &markers,
191 Mesh *mesh = fes.GetMesh();
192 const int dim = mesh->Dimension();
193 const FiniteElement &el = *fes.GetTypicalFE();
196 const int d = maps.ndof, q = maps.nqpt;
198 const GeometricFactors *geom = mesh->GetGeometricFactors(*ir, flags, mt);
199 const int map_type = fes.GetTypicalFE()->GetMapType();
200 decltype(&DLFEvalAssemble2D<>) ker =
201 dim == 2 ? DLFEvalAssemble2D<> : DLFEvalAssemble3D<>;
205 if (d==1 && q==1) { ker=DLFEvalAssemble2D<1,1>; }
206 if (d==2 && q==2) { ker=DLFEvalAssemble2D<2,2>; }
207 if (d==3 && q==3) { ker=DLFEvalAssemble2D<3,3>; }
208 if (d==4 && q==4) { ker=DLFEvalAssemble2D<4,4>; }
209 if (d==5 && q==5) { ker=DLFEvalAssemble2D<5,5>; }
210 if (d==2 && q==3) { ker=DLFEvalAssemble2D<2,3>; }
211 if (d==3 && q==4) { ker=DLFEvalAssemble2D<3,4>; }
212 if (d==4 && q==5) { ker=DLFEvalAssemble2D<4,5>; }
213 if (d==5 && q==6) { ker=DLFEvalAssemble2D<5,6>; }
218 if (d==1 && q==1) { ker=DLFEvalAssemble3D<1,1>; }
219 if (d==2 && q==2) { ker=DLFEvalAssemble3D<2,2>; }
220 if (d==3 && q==3) { ker=DLFEvalAssemble3D<3,3>; }
221 if (d==4 && q==4) { ker=DLFEvalAssemble3D<4,4>; }
222 if (d==5 && q==5) { ker=DLFEvalAssemble3D<5,5>; }
223 if (d==2 && q==3) { ker=DLFEvalAssemble3D<2,3>; }
224 if (d==3 && q==4) { ker=DLFEvalAssemble3D<3,4>; }
225 if (d==4 && q==5) { ker=DLFEvalAssemble3D<4,5>; }
226 if (d==5 && q==6) { ker=DLFEvalAssemble3D<5,6>; }
229 MFEM_VERIFY(ker,
"No kernel ndof " << d <<
" nqpt " << q);
231 const int vdim = fes.GetVDim();
232 const int ne = fes.GetMesh()->GetNE();
233 const int *M = markers.Read();
234 const real_t *B = maps.B.Read();
235 const real_t *detJ = geom->detJ.Read();
236 const real_t *W = ir->GetWeights().Read();
237 real_t *Y = y.ReadWrite();
238 ker(vdim, ne, d, q, map_type, M, B, detJ, W, coeff, Y);
246 const int qorder = oa * fe.
GetOrder() + ob;
252 DLFEvalAssemble(fes, ir, markers, coeff,
b);
260 const int qorder = 2 * fe.
GetOrder();
266 DLFEvalAssemble(fes, ir, markers, coeff,
b);
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...
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Mesh * GetMesh() const
Returns the mesh.
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
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
Class representing the storage layout of a QuadratureFunction.
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)
MemoryType
Memory types supported by MFEM.
DeviceTensor< 3, real_t > DeviceCube
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)