MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
lininteg_boundary.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 BLFEvalAssemble2D(const int vdim, const int nbe, const int d,
21 const int q,
22 const bool normals, const int *markers, const real_t *b,
23 const real_t *detj, const real_t *n, const real_t *weights,
24 const Vector &coeff, real_t *y)
25{
26 const auto F = coeff.Read();
27 const auto M = Reshape(markers, nbe);
28 const auto B = Reshape(b, q, d);
29 const auto detJ = Reshape(detj, q, nbe);
30 const auto N = Reshape(n, q, 2, nbe);
31 const auto W = Reshape(weights, q);
32 const int cvdim = normals ? 2 : 1;
33 const bool cst = coeff.Size() == cvdim;
34 const auto C = cst ? Reshape(F,cvdim,1,1) : Reshape(F,cvdim,q,nbe);
35 auto Y = Reshape(y, d, vdim, nbe);
36
37 mfem::forall(nbe, [=] MFEM_HOST_DEVICE (int e)
38 {
39 if (M(e) == 0) { return; } // ignore
40
41 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
42 real_t QQ[Q];
43
44 for (int c = 0; c < vdim; ++c)
45 {
46 for (int qx = 0; qx < q; ++qx)
47 {
48 real_t coeff_val = 0.0;
49 if (normals)
50 {
51 for (int cd = 0; cd < 2; ++cd)
52 {
53 const real_t cval = cst ? C(cd,0,0) : C(cd,qx,e);
54 coeff_val += cval * N(qx, cd, e);
55 }
56 }
57 else
58 {
59 coeff_val = cst ? C(0,0,0) : C(0,qx,e);
60 }
61 QQ[qx] = W(qx) * coeff_val * detJ(qx,e);
62 }
63 for (int dx = 0; dx < d; ++dx)
64 {
65 real_t u = 0;
66 for (int qx = 0; qx < q; ++qx) { u += QQ[qx] * B(qx,dx); }
67 Y(dx,c,e) += u;
68 }
69 }
70 });
71}
72
73template<int T_D1D = 0, int T_Q1D = 0>
74static void BLFEvalAssemble3D(const int vdim, const int nbe, const int d,
75 const int q,
76 const bool normals, const int *markers, const real_t *b,
77 const real_t *detj, const real_t *n, const real_t *weights,
78 const Vector &coeff, real_t *y)
79{
80 const auto F = coeff.Read();
81 const auto M = Reshape(markers, nbe);
82 const auto B = Reshape(b, q, d);
83 const auto detJ = Reshape(detj, q, q, nbe);
84 const auto N = Reshape(n, q, q, 3, nbe);
85 const auto W = Reshape(weights, q, q);
86 const int cvdim = normals ? 3 : 1;
87 const bool cst = coeff.Size() == cvdim;
88 const auto C = cst ? Reshape(F,cvdim,1,1,1) : Reshape(F,cvdim,q,q,nbe);
89 auto Y = Reshape(y, d, d, vdim, nbe);
90
91 mfem::forall_2D(nbe, q, q, [=] MFEM_HOST_DEVICE (int e)
92 {
93 if (M(e) == 0) { return; } // ignore
94
95 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
96 constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
97
98 MFEM_SHARED real_t sBt[Q*D];
99 MFEM_SHARED real_t sQQ[Q*Q];
100 MFEM_SHARED real_t sQD[Q*D];
101
102 const DeviceMatrix Bt(sBt, d, q);
103 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
104
105 const DeviceMatrix QQ(sQQ, q, q);
106 const DeviceMatrix QD(sQD, q, d);
107
108 for (int c = 0; c < vdim; ++c)
109 {
110 MFEM_FOREACH_THREAD(x,x,q)
111 {
112 MFEM_FOREACH_THREAD(y,y,q)
113 {
114 real_t coeff_val = 0.0;
115 if (normals)
116 {
117 for (int cd = 0; cd < 3; ++cd)
118 {
119 real_t cval = cst ? C(cd,0,0,0) : C(cd,x,y,e);
120 coeff_val += cval * N(x,y,cd,e);
121 }
122 }
123 else
124 {
125 coeff_val = cst ? C(0,0,0,0) : C(0,x,y,e);
126 }
127 QQ(y,x) = W(x,y) * coeff_val * detJ(x,y,e);
128 }
129 }
130 MFEM_SYNC_THREAD;
131 MFEM_FOREACH_THREAD(qy,y,q)
132 {
133 MFEM_FOREACH_THREAD(dx,x,d)
134 {
135 real_t u = 0.0;
136 for (int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
137 QD(qy,dx) = u;
138 }
139 }
140 MFEM_SYNC_THREAD;
141 MFEM_FOREACH_THREAD(dy,y,d)
142 {
143 MFEM_FOREACH_THREAD(dx,x,d)
144 {
145 real_t u = 0.0;
146 for (int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
147 Y(dx,dy,c,e) += u;
148 }
149 }
150 MFEM_SYNC_THREAD;
151 }
152 });
153}
154
155static void BLFEvalAssemble(const FiniteElementSpace &fes,
156 const IntegrationRule &ir,
157 const Array<int> &markers,
158 const Vector &coeff,
159 const bool normals,
160 Vector &y)
161{
162 if (fes.GetNBE() == 0) { return; }
163 Mesh &mesh = *fes.GetMesh();
164 const int dim = mesh.Dimension();
165 const FiniteElement &el = *fes.GetBE(0);
167 const DofToQuad &maps = el.GetDofToQuad(ir, DofToQuad::TENSOR);
168 const int d = maps.ndof, q = maps.nqpt;
170 if (normals) { flags |= FaceGeometricFactors::NORMALS; }
171 const FaceGeometricFactors *geom = mesh.GetFaceGeometricFactors(
172 ir, flags, FaceType::Boundary, mt);
173 auto ker = (dim == 2) ? BLFEvalAssemble2D<> : BLFEvalAssemble3D<>;
174
175 if (dim==2)
176 {
177 if (d==1 && q==1) { ker=BLFEvalAssemble2D<1,1>; }
178 if (d==2 && q==2) { ker=BLFEvalAssemble2D<2,2>; }
179 if (d==3 && q==3) { ker=BLFEvalAssemble2D<3,3>; }
180 if (d==4 && q==4) { ker=BLFEvalAssemble2D<4,4>; }
181 if (d==5 && q==5) { ker=BLFEvalAssemble2D<5,5>; }
182 if (d==2 && q==3) { ker=BLFEvalAssemble2D<2,3>; }
183 if (d==3 && q==4) { ker=BLFEvalAssemble2D<3,4>; }
184 if (d==4 && q==5) { ker=BLFEvalAssemble2D<4,5>; }
185 if (d==5 && q==6) { ker=BLFEvalAssemble2D<5,6>; }
186 }
187
188 if (dim==3)
189 {
190 if (d==1 && q==1) { ker=BLFEvalAssemble3D<1,1>; }
191 if (d==2 && q==2) { ker=BLFEvalAssemble3D<2,2>; }
192 if (d==3 && q==3) { ker=BLFEvalAssemble3D<3,3>; }
193 if (d==4 && q==4) { ker=BLFEvalAssemble3D<4,4>; }
194 if (d==5 && q==5) { ker=BLFEvalAssemble3D<5,5>; }
195 if (d==2 && q==3) { ker=BLFEvalAssemble3D<2,3>; }
196 if (d==3 && q==4) { ker=BLFEvalAssemble3D<3,4>; }
197 if (d==4 && q==5) { ker=BLFEvalAssemble3D<4,5>; }
198 if (d==5 && q==6) { ker=BLFEvalAssemble3D<5,6>; }
199 }
200
201 MFEM_VERIFY(ker, "No kernel ndof " << d << " nqpt " << q);
202
203 const int vdim = fes.GetVDim();
204 const int nbe = fes.GetMesh()->GetNFbyType(FaceType::Boundary);
205 const int *M = markers.Read();
206 const real_t *B = maps.B.Read();
207 const real_t *detJ = geom->detJ.Read();
208 const real_t *n = geom->normal.Read();
209 const real_t *W = ir.GetWeights().Read();
210 real_t *Y = y.ReadWrite();
211 ker(vdim, nbe, d, q, normals, M, B, detJ, n, W, coeff, Y);
212}
213
215 const Array<int> &markers,
216 Vector &b)
217{
218 if (fes.GetNBE() == 0) { return; }
219 const FiniteElement &fe = *fes.GetBE(0);
220 const int qorder = oa * fe.GetOrder() + ob;
221 const Geometry::Type gtype = fe.GetGeomType();
222 const IntegrationRule &ir = IntRule ? *IntRule : IntRules.Get(gtype, qorder);
223 Mesh &mesh = *fes.GetMesh();
224
227 BLFEvalAssemble(fes, ir, markers, coeff, false, b);
228}
229
231 const Array<int> &markers,
232 Vector &b)
233{
234 if (fes.GetNBE() == 0) { return; }
235 const FiniteElement &fe = *fes.GetBE(0);
236 const int qorder = oa * fe.GetOrder() + ob;
237 const Geometry::Type gtype = fe.GetGeomType();
238 const IntegrationRule &ir = IntRule ? *IntRule : IntRules.Get(gtype, qorder);
239 Mesh &mesh = *fes.GetMesh();
240
243 BLFEvalAssemble(fes, ir, markers, coeff, true, b);
244}
245
246} // namespace mfem
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
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
Class representing the storage layout of a FaceQuadratureFunction.
Definition qspace.hpp:170
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
Definition fespace.cpp:3881
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition fespace.hpp:900
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
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
Mesh data type.
Definition mesh.hpp:64
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.
void forall(int N, lambda &&body)
Definition forall.hpp:753
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492