13 #include "../fem/kernels.hpp" 14 #include "../general/forall.hpp" 19 template<
int T_D1D = 0,
int T_Q1D = 0>
static 20 void DLFEvalAssemble2D(
const int vdim,
const int ne,
const int d,
const int q,
21 const int map_type,
const int *markers,
const double *
b,
22 const double *detj,
const double *weights,
23 const Vector &coeff,
double *y)
25 const auto F = coeff.Read();
26 const auto M =
Reshape(markers, ne);
28 const auto DETJ =
Reshape(detj, q, q, ne);
29 const auto W =
Reshape(weights, q, q);
30 const bool cst = coeff.Size() == vdim;
31 const auto C = cst ?
Reshape(F,vdim,1,1,1) :
Reshape(F,vdim,q,q,ne);
32 auto Y =
Reshape(y, d,d, vdim, ne);
34 MFEM_FORALL_2D(e, ne, q, q, 1,
36 if (M(e) == 0) {
return; }
38 constexpr
int Q = T_Q1D ? T_Q1D :
MAX_Q1D;
39 constexpr
int D = T_D1D ? T_D1D :
MAX_D1D;
41 MFEM_SHARED
double sBt[Q*D];
42 MFEM_SHARED
double sQQ[Q*Q];
43 MFEM_SHARED
double sQD[Q*D];
46 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
51 for (
int c = 0; c < vdim; ++c)
53 const double cst_val = C(c,0,0,0);
54 MFEM_FOREACH_THREAD(x,x,q)
56 MFEM_FOREACH_THREAD(y,y,q)
59 const double coeff_val = cst ? cst_val : C(c,x,y,e);
60 QQ(y,x) = W(x,y) * coeff_val * detJ;
64 MFEM_FOREACH_THREAD(qy,y,q)
66 MFEM_FOREACH_THREAD(dx,x,d)
69 for (
int qx = 0; qx < q; ++qx) {
u += QQ(qy,qx) * Bt(dx,qx); }
74 MFEM_FOREACH_THREAD(dy,y,d)
76 MFEM_FOREACH_THREAD(dx,x,d)
79 for (
int qy = 0; qy < q; ++qy) {
u += QD(qy,dx) * Bt(dy,qy); }
88 template<
int T_D1D = 0,
int T_Q1D = 0>
static 89 void DLFEvalAssemble3D(
const int vdim,
const int ne,
const int d,
const int q,
90 const int map_type,
const int *markers,
const double *
b,
91 const double *detj,
const double *weights,
92 const Vector &coeff,
double *y)
94 const auto F = coeff.Read();
95 const auto M =
Reshape(markers, ne);
97 const auto DETJ =
Reshape(detj, q, q, q, ne);
98 const auto W =
Reshape(weights, q,q,q);
99 const bool cst_coeff = coeff.Size() == vdim;
100 const auto C = cst_coeff ?
Reshape(F,vdim,1,1,1,1):
Reshape(F,vdim,q,q,q,ne);
102 auto Y =
Reshape(y, d,d,d, vdim, ne);
104 MFEM_FORALL_2D(e, ne, q, q, 1,
106 if (M(e) == 0) {
return; }
108 constexpr
int Q = T_Q1D ? T_Q1D :
MAX_Q1D;
109 constexpr
int D = T_D1D ? T_D1D :
MAX_D1D;
113 MFEM_SHARED
double sBt[Q*D];
115 kernels::internal::LoadB<D,Q>(d,q,B,sBt);
117 MFEM_SHARED
double sQQQ[Q*Q*Q];
120 for (
int c = 0; c < vdim; ++c)
122 const double cst_val = C(c,0,0,0,0);
123 MFEM_FOREACH_THREAD(x,x,q)
125 MFEM_FOREACH_THREAD(y,y,q)
127 for (
int z = 0; z < q; ++z)
130 const double coeff_val = cst_coeff ? cst_val : C(c,x,y,z,e);
131 QQQ(z,y,x) = W(x,y,z) * coeff_val * detJ;
136 MFEM_FOREACH_THREAD(qx,x,q)
138 MFEM_FOREACH_THREAD(qy,y,q)
140 for (
int dz = 0; dz < d; ++dz) {
u[dz] = 0.0; }
141 for (
int qz = 0; qz < q; ++qz)
143 const double ZYX = QQQ(qz,qy,qx);
144 for (
int dz = 0; dz < d; ++dz) {
u[dz] += ZYX * Bt(dz,qz); }
146 for (
int dz = 0; dz < d; ++dz) { QQQ(dz,qy,qx) =
u[dz]; }
150 MFEM_FOREACH_THREAD(dz,y,d)
152 MFEM_FOREACH_THREAD(qx,x,q)
154 for (
int dy = 0; dy < d; ++dy) {
u[dy] = 0.0; }
155 for (
int qy = 0; qy < q; ++qy)
157 const double zYX = QQQ(dz,qy,qx);
158 for (
int dy = 0; dy < d; ++dy) {
u[dy] += zYX * Bt(dy,qy); }
160 for (
int dy = 0; dy < d; ++dy) { QQQ(dz,dy,qx) =
u[dy]; }
164 MFEM_FOREACH_THREAD(dz,y,d)
166 MFEM_FOREACH_THREAD(dy,x,d)
168 for (
int dx = 0; dx < d; ++dx) {
u[dx] = 0.0; }
169 for (
int qx = 0; qx < q; ++qx)
171 const double zyX = QQQ(dz,dy,qx);
172 for (
int dx = 0; dx < d; ++dx) {
u[dx] += zyX * Bt(dx,qx); }
174 for (
int dx = 0; dx < d; ++dx) { Y(dx,dy,dz,c,e) +=
u[dx]; }
182 static void DLFEvalAssemble(
const FiniteElementSpace &fes,
183 const IntegrationRule *ir,
184 const Array<int> &markers,
188 Mesh *mesh = fes.GetMesh();
189 const int dim = mesh->Dimension();
190 const FiniteElement &el = *fes.GetFE(0);
193 const int d = maps.ndof, q = maps.nqpt;
195 const GeometricFactors *geom = mesh->GetGeometricFactors(*ir, flags, mt);
196 const int map_type = fes.GetFE(0)->GetMapType();
197 decltype(&DLFEvalAssemble2D<>) ker =
198 dim == 2 ? DLFEvalAssemble2D<> : DLFEvalAssemble3D<>;
202 if (d==1 && q==1) { ker=DLFEvalAssemble2D<1,1>; }
203 if (d==2 && q==2) { ker=DLFEvalAssemble2D<2,2>; }
204 if (d==3 && q==3) { ker=DLFEvalAssemble2D<3,3>; }
205 if (d==4 && q==4) { ker=DLFEvalAssemble2D<4,4>; }
206 if (d==5 && q==5) { ker=DLFEvalAssemble2D<5,5>; }
207 if (d==2 && q==3) { ker=DLFEvalAssemble2D<2,3>; }
208 if (d==3 && q==4) { ker=DLFEvalAssemble2D<3,4>; }
209 if (d==4 && q==5) { ker=DLFEvalAssemble2D<4,5>; }
210 if (d==5 && q==6) { ker=DLFEvalAssemble2D<5,6>; }
215 if (d==1 && q==1) { ker=DLFEvalAssemble3D<1,1>; }
216 if (d==2 && q==2) { ker=DLFEvalAssemble3D<2,2>; }
217 if (d==3 && q==3) { ker=DLFEvalAssemble3D<3,3>; }
218 if (d==4 && q==4) { ker=DLFEvalAssemble3D<4,4>; }
219 if (d==5 && q==5) { ker=DLFEvalAssemble3D<5,5>; }
220 if (d==2 && q==3) { ker=DLFEvalAssemble3D<2,3>; }
221 if (d==3 && q==4) { ker=DLFEvalAssemble3D<3,4>; }
222 if (d==4 && q==5) { ker=DLFEvalAssemble3D<4,5>; }
223 if (d==5 && q==6) { ker=DLFEvalAssemble3D<5,6>; }
226 MFEM_VERIFY(ker,
"No kernel ndof " << d <<
" nqpt " << q);
228 const int vdim = fes.GetVDim();
229 const int ne = fes.GetMesh()->GetNE();
230 const int *M = markers.Read();
231 const double *B = maps.B.Read();
232 const double *detJ = geom->detJ.Read();
233 const double *W = ir->GetWeights().Read();
234 double *Y = y.ReadWrite();
235 ker(vdim, ne, d, q, map_type, M, B, detJ, W, coeff, Y);
243 const int qorder = oa * fe.
GetOrder() + ob;
249 DLFEvalAssemble(fes, ir, markers, coeff,
b);
257 const int qorder = 2 * fe.
GetOrder();
263 DLFEvalAssemble(fes, ir, markers, coeff,
b);
Abstract class for all finite elements.
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.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Class to represent a coefficient evaluated at quadrature points.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
DeviceTensor< 2, double > DeviceMatrix
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Enable all above compressions.
Mesh * GetMesh() const
Returns the mesh.
DeviceTensor< 3, double > DeviceCube
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
MemoryType
Memory types supported by MFEM.
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Class representing the storage layout of a QuadratureFunction.
double u(const Vector &xvec)
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...