13 #include "../fem/kernels.hpp"
14 #include "../general/forall.hpp"
19 template<
int T_D1D = 0,
int T_Q1D = 0>
static
20 void BLFEvalAssemble2D(
const int vdim,
const int nbe,
const int d,
const int q,
21 const bool normals,
const int *markers,
const double *
b,
22 const double *detj,
const double *n,
const double *weights,
23 const Vector &coeff,
double *y)
25 const auto F = coeff.Read();
26 const auto M =
Reshape(markers, nbe);
27 const auto B =
Reshape(b, q, d);
28 const auto detJ =
Reshape(detj, q, nbe);
29 const auto N =
Reshape(n, q, 2, nbe);
30 const auto W =
Reshape(weights, q);
31 const int cvdim = normals ? 2 : 1;
32 const bool cst = coeff.Size() == cvdim;
33 const auto C = cst ?
Reshape(F,cvdim,1,1) : Reshape(F,cvdim,q,nbe);
34 auto Y =
Reshape(y, d, vdim, nbe);
38 if (M(e) == 0) {
return; }
40 constexpr
int Q = T_Q1D ? T_Q1D :
MAX_Q1D;
43 for (
int c = 0; c < vdim; ++c)
45 for (
int qx = 0; qx < q; ++qx)
47 double coeff_val = 0.0;
50 for (
int cd = 0; cd < 2; ++cd)
52 const double cval = cst ? C(cd,0,0) : C(cd,qx,e);
53 coeff_val += cval * N(qx, cd, e);
58 coeff_val = cst ? C(0,0,0) : C(0,qx,e);
60 QQ[qx] = W(qx) * coeff_val * detJ(qx,e);
62 for (
int dx = 0; dx < d; ++dx)
65 for (
int qx = 0; qx < q; ++qx) { u += QQ[qx] * B(qx,dx); }
72 template<
int T_D1D = 0,
int T_Q1D = 0>
static
73 void BLFEvalAssemble3D(
const int vdim,
const int nbe,
const int d,
const int q,
74 const bool normals,
const int *markers,
const double *b,
75 const double *detj,
const double *n,
const double *weights,
76 const Vector &coeff,
double *y)
78 const auto F = coeff.Read();
79 const auto M =
Reshape(markers, nbe);
80 const auto B =
Reshape(b, q, d);
81 const auto detJ =
Reshape(detj, q, q, nbe);
82 const auto N =
Reshape(n, q, q, 3, nbe);
83 const auto W =
Reshape(weights, q, q);
84 const int cvdim = normals ? 3 : 1;
85 const bool cst = coeff.Size() == cvdim;
86 const auto C = cst ?
Reshape(F,cvdim,1,1,1) : Reshape(F,cvdim,q,q,nbe);
87 auto Y =
Reshape(y, d, d, vdim, nbe);
89 MFEM_FORALL_2D(e, nbe, q, q, 1,
91 if (M(e) == 0) {
return; }
93 constexpr
int Q = T_Q1D ? T_Q1D :
MAX_Q1D;
94 constexpr
int D = T_D1D ? T_D1D :
MAX_D1D;
96 MFEM_SHARED
double sBt[Q*D];
97 MFEM_SHARED
double sQQ[Q*Q];
98 MFEM_SHARED
double sQD[Q*D];
101 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
106 for (
int c = 0; c < vdim; ++c)
108 MFEM_FOREACH_THREAD(x,x,q)
110 MFEM_FOREACH_THREAD(y,y,q)
112 double coeff_val = 0.0;
115 for (
int cd = 0; cd < 3; ++cd)
117 double cval = cst ? C(cd,0,0,0) : C(cd,x,y,e);
118 coeff_val += cval * N(x,y,cd,e);
123 coeff_val = cst ? C(0,0,0,0) : C(0,x,y,e);
125 QQ(y,x) = W(x,y) * coeff_val * detJ(x,y,e);
129 MFEM_FOREACH_THREAD(qy,y,q)
131 MFEM_FOREACH_THREAD(dx,x,d)
134 for (
int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
139 MFEM_FOREACH_THREAD(dy,y,d)
141 MFEM_FOREACH_THREAD(dx,x,d)
144 for (
int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
153 static void BLFEvalAssemble(
const FiniteElementSpace &fes,
154 const IntegrationRule &ir,
155 const Array<int> &markers,
160 Mesh &mesh = *fes.GetMesh();
161 const int dim = mesh.Dimension();
162 const FiniteElement &el = *fes.GetBE(0);
165 const int d = maps.ndof, q = maps.nqpt;
168 const FaceGeometricFactors *geom = mesh.GetFaceGeometricFactors(
170 auto ker = (dim == 2) ? BLFEvalAssemble2D<> : BLFEvalAssemble3D<>;
174 if (d==1 && q==1) { ker=BLFEvalAssemble2D<1,1>; }
175 if (d==2 && q==2) { ker=BLFEvalAssemble2D<2,2>; }
176 if (d==3 && q==3) { ker=BLFEvalAssemble2D<3,3>; }
177 if (d==4 && q==4) { ker=BLFEvalAssemble2D<4,4>; }
178 if (d==5 && q==5) { ker=BLFEvalAssemble2D<5,5>; }
179 if (d==2 && q==3) { ker=BLFEvalAssemble2D<2,3>; }
180 if (d==3 && q==4) { ker=BLFEvalAssemble2D<3,4>; }
181 if (d==4 && q==5) { ker=BLFEvalAssemble2D<4,5>; }
182 if (d==5 && q==6) { ker=BLFEvalAssemble2D<5,6>; }
187 if (d==1 && q==1) { ker=BLFEvalAssemble3D<1,1>; }
188 if (d==2 && q==2) { ker=BLFEvalAssemble3D<2,2>; }
189 if (d==3 && q==3) { ker=BLFEvalAssemble3D<3,3>; }
190 if (d==4 && q==4) { ker=BLFEvalAssemble3D<4,4>; }
191 if (d==5 && q==5) { ker=BLFEvalAssemble3D<5,5>; }
192 if (d==2 && q==3) { ker=BLFEvalAssemble3D<2,3>; }
193 if (d==3 && q==4) { ker=BLFEvalAssemble3D<3,4>; }
194 if (d==4 && q==5) { ker=BLFEvalAssemble3D<4,5>; }
195 if (d==5 && q==6) { ker=BLFEvalAssemble3D<5,6>; }
198 MFEM_VERIFY(ker,
"No kernel ndof " << d <<
" nqpt " << q);
200 const int vdim = fes.GetVDim();
202 const int *M = markers.Read();
203 const double *B = maps.B.Read();
204 const double *detJ = geom->detJ.Read();
205 const double *n = geom->normal.Read();
206 const double *W = ir.GetWeights().Read();
207 double *Y = y.ReadWrite();
208 ker(vdim, nbe, d, q, normals, M, B, detJ, n, W, coeff, Y);
216 const int qorder = oa * fe.
GetOrder() + ob;
223 BLFEvalAssemble(fes, ir, markers, coeff,
false, b);
231 const int qorder = oa * fe.
GetOrder() + ob;
238 BLFEvalAssemble(fes, ir, markers, coeff,
true, 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...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Class to represent a coefficient evaluated at quadrature points.
DeviceTensor< 2, double > DeviceMatrix
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Mesh * GetMesh() const
Returns the mesh.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Enable all above compressions.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
MemoryType
Memory types supported by MFEM.
Class representing the storage layout of a FaceQuadratureFunction.
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
double u(const Vector &xvec)
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.