MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
lininteg_domain.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
13#include "../../fem/kernels.hpp"
14#include "../fem.hpp"
15
16namespace mfem
17{
18
19template<int T_D1D = 0, int T_Q1D = 0>
20static void DLFEvalAssemble2D(const int vdim, const int ne, const int d,
21 const int q,
22 const int map_type, const int *markers, const real_t *b,
23 const real_t *detj, const real_t *weights,
24 const Vector &coeff, real_t *y)
25{
26 const auto F = coeff.Read();
27 const auto M = Reshape(markers, ne);
28 const auto B = Reshape(b, q, d);
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);
34
35 mfem::forall_2D(ne, q, q, [=] MFEM_HOST_DEVICE (int e)
36 {
37 if (M(e) == 0) { return; } // ignore
38
39 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
40 constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
41
42 MFEM_SHARED real_t sBt[Q*D];
43 MFEM_SHARED real_t sQQ[Q*Q];
44 MFEM_SHARED real_t sQD[Q*D];
45
46 const DeviceMatrix Bt(sBt, d, q);
47 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
48
49 const DeviceMatrix QQ(sQQ, q, q);
50 const DeviceMatrix QD(sQD, q, d);
51
52 for (int c = 0; c < vdim; ++c)
53 {
54 const real_t cst_val = C(c,0,0,0);
55 MFEM_FOREACH_THREAD(x,x,q)
56 {
57 MFEM_FOREACH_THREAD(y,y,q)
58 {
59 const real_t detJ = (map_type == FiniteElement::VALUE) ? DETJ(x,y,e) : 1.0;
60 const real_t coeff_val = cst ? cst_val : C(c,x,y,e);
61 QQ(y,x) = W(x,y) * coeff_val * detJ;
62 }
63 }
64 MFEM_SYNC_THREAD;
65 MFEM_FOREACH_THREAD(qy,y,q)
66 {
67 MFEM_FOREACH_THREAD(dx,x,d)
68 {
69 real_t u = 0.0;
70 for (int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
71 QD(qy,dx) = u;
72 }
73 }
74 MFEM_SYNC_THREAD;
75 MFEM_FOREACH_THREAD(dy,y,d)
76 {
77 MFEM_FOREACH_THREAD(dx,x,d)
78 {
79 real_t u = 0.0;
80 for (int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
81 Y(dx,dy,c,e) += u;
82 }
83 }
84 MFEM_SYNC_THREAD;
85 }
86 });
87}
88
89template<int T_D1D = 0, int T_Q1D = 0>
90static void DLFEvalAssemble3D(const int vdim, const int ne, const int d,
91 const int q,
92 const int map_type, const int *markers, const real_t *b,
93 const real_t *detj, const real_t *weights,
94 const Vector &coeff, real_t *y)
95{
96 const auto F = coeff.Read();
97 const auto M = Reshape(markers, ne);
98 const auto B = Reshape(b, q,d);
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);
103
104 auto Y = Reshape(y, d,d,d, vdim, ne);
105
106 mfem::forall_2D(ne, q, q, [=] MFEM_HOST_DEVICE (int e)
107 {
108 if (M(e) == 0) { return; } // ignore
109
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;
113
114 real_t u[D];
115
116 MFEM_SHARED real_t sBt[Q*D];
117 const DeviceMatrix Bt(sBt, d,q);
118 kernels::internal::LoadB<D,Q>(d,q,B,sBt);
119
120 MFEM_SHARED real_t sQQQ[MQD*MQD*MQD];
121 const DeviceCube QQQ(sQQQ, MQD, MQD, MQD);
122
123 for (int c = 0; c < vdim; ++c)
124 {
125 const real_t cst_val = C(c,0,0,0,0);
126 MFEM_FOREACH_THREAD(x,x,q)
127 {
128 MFEM_FOREACH_THREAD(y,y,q)
129 {
130 for (int z = 0; z < q; ++z)
131 {
132 const real_t detJ = (map_type == FiniteElement::VALUE) ? DETJ(x,y,z,e) : 1.0;
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;
135 }
136 }
137 }
138 MFEM_SYNC_THREAD;
139 MFEM_FOREACH_THREAD(qx,x,q)
140 {
141 MFEM_FOREACH_THREAD(qy,y,q)
142 {
143 for (int dz = 0; dz < d; ++dz) { u[dz] = 0.0; }
144 for (int qz = 0; qz < q; ++qz)
145 {
146 const real_t ZYX = QQQ(qz,qy,qx);
147 for (int dz = 0; dz < d; ++dz) { u[dz] += ZYX * Bt(dz,qz); }
148 }
149 for (int dz = 0; dz < d; ++dz) { QQQ(dz,qy,qx) = u[dz]; }
150 }
151 }
152 MFEM_SYNC_THREAD;
153 MFEM_FOREACH_THREAD(dz,y,d)
154 {
155 MFEM_FOREACH_THREAD(qx,x,q)
156 {
157 for (int dy = 0; dy < d; ++dy) { u[dy] = 0.0; }
158 for (int qy = 0; qy < q; ++qy)
159 {
160 const real_t zYX = QQQ(dz,qy,qx);
161 for (int dy = 0; dy < d; ++dy) { u[dy] += zYX * Bt(dy,qy); }
162 }
163 for (int dy = 0; dy < d; ++dy) { QQQ(dz,dy,qx) = u[dy]; }
164 }
165 }
166 MFEM_SYNC_THREAD;
167 MFEM_FOREACH_THREAD(dz,y,d)
168 {
169 MFEM_FOREACH_THREAD(dy,x,d)
170 {
171 for (int dx = 0; dx < d; ++dx) { u[dx] = 0.0; }
172 for (int qx = 0; qx < q; ++qx)
173 {
174 const real_t zyX = QQQ(dz,dy,qx);
175 for (int dx = 0; dx < d; ++dx) { u[dx] += zyX * Bt(dx,qx); }
176 }
177 for (int dx = 0; dx < d; ++dx) { Y(dx,dy,dz,c,e) += u[dx]; }
178 }
179 }
180 MFEM_SYNC_THREAD;
181 }
182 });
183}
184
185static void DLFEvalAssemble(const FiniteElementSpace &fes,
186 const IntegrationRule *ir,
187 const Array<int> &markers,
188 const Vector &coeff,
189 Vector &y)
190{
191 Mesh *mesh = fes.GetMesh();
192 const int dim = mesh->Dimension();
193 const FiniteElement &el = *fes.GetTypicalFE();
195 const DofToQuad &maps = el.GetDofToQuad(*ir, DofToQuad::TENSOR);
196 const int d = maps.ndof, q = maps.nqpt;
197 constexpr int flags = GeometricFactors::DETERMINANTS;
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<>;
202
203 if (dim==2)
204 {
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>; }
214 }
215
216 if (dim==3)
217 {
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>; }
227 }
228
229 MFEM_VERIFY(ker, "No kernel ndof " << d << " nqpt " << q);
230
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);
239}
240
242 const Array<int> &markers,
243 Vector &b)
244{
245 const FiniteElement &fe = *fes.GetTypicalFE();
246 const int qorder = oa * fe.GetOrder() + ob;
247 const Geometry::Type gtype = fe.GetGeomType();
248 const IntegrationRule *ir = IntRule ? IntRule : &IntRules.Get(gtype, qorder);
249
250 QuadratureSpace qs(*fes.GetMesh(), *ir);
252 DLFEvalAssemble(fes, ir, markers, coeff, b);
253}
254
256 const Array<int> &markers,
257 Vector &b)
258{
259 const FiniteElement &fe = *fes.GetTypicalFE();
260 const int qorder = 2 * fe.GetOrder();
261 const Geometry::Type gtype = fe.GetGeomType();
262 const IntegrationRule *ir = IntRule ? IntRule : &IntRules.Get(gtype, qorder);
263
264 QuadratureSpace qs(*fes.GetMesh(), *ir);
266 DLFEvalAssemble(fes, ir, markers, coeff, b);
267}
268
269} // namespace mfem
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...
Definition device.hpp:274
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:165
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...
Definition fespace.hpp:244
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
Abstract class for all finite elements.
Definition fe_base.hpp:244
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:338
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:331
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
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.
Definition qspace.hpp:120
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
Vector data type.
Definition vector.hpp:82
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
DeviceTensor< 2, real_t > DeviceMatrix
Definition dtensor.hpp:143
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
@ COMPRESSED
Enable all above compressions.
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:762
float real_t
Definition config.hpp:43
MemoryType
Memory types supported by MFEM.
DeviceTensor< 3, real_t > DeviceCube
Definition dtensor.hpp:146
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492