MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
lininteg_boundary_flux.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> static
20void BFLFEvalAssemble2D(const int nbe, const int d, const int q,
21 const int *markers, const real_t *b,
22 const real_t *weights, const Vector &coeff, real_t *y)
23{
24 const auto F = coeff.Read();
25 const auto M = Reshape(markers, nbe);
26 const auto B = Reshape(b, q, d);
27 const auto W = Reshape(weights, q);
28 const bool const_coeff = coeff.Size() == 1;
29 const auto C = const_coeff ? Reshape(F,1,1) : Reshape(F,q,nbe);
30 auto Y = Reshape(y, d, nbe);
31
32 mfem::forall(nbe, [=] MFEM_HOST_DEVICE (int e)
33 {
34 if (M(e) == 0) { return; } // ignore (in a lambda return acts as continue)
35
36 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
37 real_t QQ[Q];
38
39 for (int qx = 0; qx < q; ++qx)
40 {
41 const real_t coeff_val = const_coeff ? C(0,0) : C(qx,e);
42 QQ[qx] = W(qx) * coeff_val;
43 }
44 for (int dx = 0; dx < d; ++dx)
45 {
46 real_t u = 0;
47 for (int qx = 0; qx < q; ++qx) { u += QQ[qx] * B(qx,dx); }
48 Y(dx,e) += u;
49 }
50 });
51}
52
53template<int T_D1D = 0, int T_Q1D = 0> static
54void BFLFEvalAssemble3D(const int nbe, const int d, const int q,
55 const int *markers, const real_t *b,
56 const real_t *weights, const Vector &coeff, real_t *y)
57{
58 const auto F = coeff.Read();
59 const auto M = Reshape(markers, nbe);
60 const auto B = Reshape(b, q, d);
61 const auto W = Reshape(weights, q, q);
62 const bool const_coeff = coeff.Size() == 1;
63 const auto C = const_coeff ? Reshape(F,1,1,1) : Reshape(F,q,q,nbe);
64 auto Y = Reshape(y, d, d, nbe);
65
66 mfem::forall_2D(nbe, q, q, [=] MFEM_HOST_DEVICE (int e)
67 {
68 if (M(e) == 0) { return; } // ignore
69
70 constexpr int Q = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
71 constexpr int D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
72
73 MFEM_SHARED real_t sBt[Q*D];
74 MFEM_SHARED real_t sQQ[Q*Q];
75 MFEM_SHARED real_t sQD[Q*D];
76
77 const DeviceMatrix Bt(sBt, d, q);
78 kernels::internal::LoadB<D,Q>(d, q, B, sBt);
79
80 const DeviceMatrix QQ(sQQ, q, q);
81 const DeviceMatrix QD(sQD, q, d);
82
83 MFEM_FOREACH_THREAD(x,x,q)
84 {
85 MFEM_FOREACH_THREAD(y,y,q)
86 {
87 const real_t coeff_val = const_coeff ? C(0,0,0) : C(x,y,e);
88 QQ(y,x) = W(x,y) * coeff_val;
89 }
90 }
91 MFEM_SYNC_THREAD;
92 MFEM_FOREACH_THREAD(qy,y,q)
93 {
94 MFEM_FOREACH_THREAD(dx,x,d)
95 {
96 real_t u = 0.0;
97 for (int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
98 QD(qy,dx) = u;
99 }
100 }
101 MFEM_SYNC_THREAD;
102 MFEM_FOREACH_THREAD(dy,y,d)
103 {
104 MFEM_FOREACH_THREAD(dx,x,d)
105 {
106 real_t u = 0.0;
107 for (int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
108 Y(dx,dy,e) += u;
109 }
110 }
111 MFEM_SYNC_THREAD;
112 });
113}
114
115static void BFLFEvalAssemble(const FiniteElementSpace &fes,
116 const IntegrationRule &ir,
117 const Array<int> &markers,
118 const Vector &coeff,
119 Vector &y)
120{
121 Mesh &mesh = *fes.GetMesh();
122 const int dim = mesh.Dimension();
123 const FiniteElement &el = *fes.GetBE(0);
124 const DofToQuad &maps = el.GetDofToQuad(ir, DofToQuad::TENSOR);
125 const int d = maps.ndof, q = maps.nqpt;
126 auto ker = (dim == 2) ? BFLFEvalAssemble2D<> : BFLFEvalAssemble3D<>;
127
128 if (dim==2)
129 {
130 if (d==1 && q==1) { ker=BFLFEvalAssemble2D<1,1>; }
131 if (d==2 && q==2) { ker=BFLFEvalAssemble2D<2,2>; }
132 if (d==3 && q==3) { ker=BFLFEvalAssemble2D<3,3>; }
133 if (d==4 && q==4) { ker=BFLFEvalAssemble2D<4,4>; }
134 if (d==5 && q==5) { ker=BFLFEvalAssemble2D<5,5>; }
135 if (d==2 && q==3) { ker=BFLFEvalAssemble2D<2,3>; }
136 if (d==3 && q==4) { ker=BFLFEvalAssemble2D<3,4>; }
137 if (d==4 && q==5) { ker=BFLFEvalAssemble2D<4,5>; }
138 if (d==5 && q==6) { ker=BFLFEvalAssemble2D<5,6>; }
139 }
140
141 if (dim==3)
142 {
143 if (d==1 && q==1) { ker=BFLFEvalAssemble3D<1,1>; }
144 if (d==2 && q==2) { ker=BFLFEvalAssemble3D<2,2>; }
145 if (d==3 && q==3) { ker=BFLFEvalAssemble3D<3,3>; }
146 if (d==4 && q==4) { ker=BFLFEvalAssemble3D<4,4>; }
147 if (d==5 && q==5) { ker=BFLFEvalAssemble3D<5,5>; }
148 if (d==2 && q==3) { ker=BFLFEvalAssemble3D<2,3>; }
149 if (d==3 && q==4) { ker=BFLFEvalAssemble3D<3,4>; }
150 if (d==4 && q==5) { ker=BFLFEvalAssemble3D<4,5>; }
151 if (d==5 && q==6) { ker=BFLFEvalAssemble3D<5,6>; }
152 }
153
154 MFEM_VERIFY(ker, "No kernel ndof " << d << " nqpt " << q);
155
156 const int nbe = fes.GetMesh()->GetNFbyType(FaceType::Boundary);
157 const int *M = markers.Read();
158 const real_t *B = maps.B.Read();
159 const real_t *W = ir.GetWeights().Read();
160 real_t *Y = y.ReadWrite();
161 ker(nbe, d, q, M, B, W, coeff, Y);
162}
163
165 const FiniteElementSpace &fes,
166 const Array<int> &markers,
167 Vector &b)
168{
169 const FiniteElement &fe = *fes.GetBE(0);
170 const int qorder = oa * fe.GetOrder() + ob;
171 const Geometry::Type gtype = fe.GetGeomType();
172 const IntegrationRule &ir = IntRule ? *IntRule : IntRules.Get(gtype, qorder);
173 Mesh &mesh = *fes.GetMesh();
174
177 BFLFEvalAssemble(fes, ir, markers, coeff, b);
178}
179
180} // namespace mfem
Class to represent a coefficient evaluated at quadrature points.
@ 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
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
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
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