MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lininteg_boundary.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 
12 #include "fem.hpp"
13 #include "../fem/kernels.hpp"
14 #include "../general/forall.hpp"
15 
16 namespace mfem
17 {
18 
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)
24 {
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);
35 
36  MFEM_FORALL(e, nbe,
37  {
38  if (M(e) == 0) { return; } // ignore
39 
40  constexpr int Q = T_Q1D ? T_Q1D : MAX_Q1D;
41  double QQ[Q];
42 
43  for (int c = 0; c < vdim; ++c)
44  {
45  for (int qx = 0; qx < q; ++qx)
46  {
47  double coeff_val = 0.0;
48  if (normals)
49  {
50  for (int cd = 0; cd < 2; ++cd)
51  {
52  const double cval = cst ? C(cd,0,0) : C(cd,qx,e);
53  coeff_val += cval * N(qx, cd, e);
54  }
55  }
56  else
57  {
58  coeff_val = cst ? C(0,0,0) : C(0,qx,e);
59  }
60  QQ[qx] = W(qx) * coeff_val * detJ(qx,e);
61  }
62  for (int dx = 0; dx < d; ++dx)
63  {
64  double u = 0;
65  for (int qx = 0; qx < q; ++qx) { u += QQ[qx] * B(qx,dx); }
66  Y(dx,c,e) += u;
67  }
68  }
69  });
70 }
71 
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)
77 {
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);
88 
89  MFEM_FORALL_2D(e, nbe, q, q, 1,
90  {
91  if (M(e) == 0) { return; } // ignore
92 
93  constexpr int Q = T_Q1D ? T_Q1D : MAX_Q1D;
94  constexpr int D = T_D1D ? T_D1D : MAX_D1D;
95 
96  MFEM_SHARED double sBt[Q*D];
97  MFEM_SHARED double sQQ[Q*Q];
98  MFEM_SHARED double sQD[Q*D];
99 
100  const DeviceMatrix Bt(sBt, d, q);
101  kernels::internal::LoadB<D,Q>(d, q, B, sBt);
102 
103  const DeviceMatrix QQ(sQQ, q, q);
104  const DeviceMatrix QD(sQD, q, d);
105 
106  for (int c = 0; c < vdim; ++c)
107  {
108  MFEM_FOREACH_THREAD(x,x,q)
109  {
110  MFEM_FOREACH_THREAD(y,y,q)
111  {
112  double coeff_val = 0.0;
113  if (normals)
114  {
115  for (int cd = 0; cd < 3; ++cd)
116  {
117  double cval = cst ? C(cd,0,0,0) : C(cd,x,y,e);
118  coeff_val += cval * N(x,y,cd,e);
119  }
120  }
121  else
122  {
123  coeff_val = cst ? C(0,0,0,0) : C(0,x,y,e);
124  }
125  QQ(y,x) = W(x,y) * coeff_val * detJ(x,y,e);
126  }
127  }
128  MFEM_SYNC_THREAD;
129  MFEM_FOREACH_THREAD(qy,y,q)
130  {
131  MFEM_FOREACH_THREAD(dx,x,d)
132  {
133  double u = 0.0;
134  for (int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
135  QD(qy,dx) = u;
136  }
137  }
138  MFEM_SYNC_THREAD;
139  MFEM_FOREACH_THREAD(dy,y,d)
140  {
141  MFEM_FOREACH_THREAD(dx,x,d)
142  {
143  double u = 0.0;
144  for (int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
145  Y(dx,dy,c,e) += u;
146  }
147  }
148  MFEM_SYNC_THREAD;
149  }
150  });
151 }
152 
153 static void BLFEvalAssemble(const FiniteElementSpace &fes,
154  const IntegrationRule &ir,
155  const Array<int> &markers,
156  const Vector &coeff,
157  const bool normals,
158  Vector &y)
159 {
160  Mesh &mesh = *fes.GetMesh();
161  const int dim = mesh.Dimension();
162  const FiniteElement &el = *fes.GetBE(0);
164  const DofToQuad &maps = el.GetDofToQuad(ir, DofToQuad::TENSOR);
165  const int d = maps.ndof, q = maps.nqpt;
167  if (normals) { flags |= FaceGeometricFactors::NORMALS; }
168  const FaceGeometricFactors *geom = mesh.GetFaceGeometricFactors(
169  ir, flags, FaceType::Boundary, mt);
170  auto ker = (dim == 2) ? BLFEvalAssemble2D<> : BLFEvalAssemble3D<>;
171 
172  if (dim==2)
173  {
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>; }
183  }
184 
185  if (dim==3)
186  {
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>; }
196  }
197 
198  MFEM_VERIFY(ker, "No kernel ndof " << d << " nqpt " << q);
199 
200  const int vdim = fes.GetVDim();
201  const int nbe = fes.GetMesh()->GetNFbyType(FaceType::Boundary);
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);
209 }
210 
212  const Array<int> &markers,
213  Vector &b)
214 {
215  const FiniteElement &fe = *fes.GetBE(0);
216  const int qorder = oa * fe.GetOrder() + ob;
217  const Geometry::Type gtype = fe.GetGeomType();
218  const IntegrationRule &ir = IntRule ? *IntRule : IntRules.Get(gtype, qorder);
219  Mesh &mesh = *fes.GetMesh();
220 
223  BLFEvalAssemble(fes, ir, markers, coeff, false, b);
224 }
225 
227  const Array<int> &markers,
228  Vector &b)
229 {
230  const FiniteElement &fe = *fes.GetBE(0);
231  const int qorder = oa * fe.GetOrder() + ob;
232  const Geometry::Type gtype = fe.GetGeomType();
233  const IntegrationRule &ir = IntRule ? *IntRule : IntRules.Get(gtype, qorder);
234  Mesh &mesh = *fes.GetMesh();
235 
238  BLFEvalAssemble(fes, ir, markers, coeff, true, b);
239 }
240 
241 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:162
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
const IntegrationRule * IntRule
Definition: lininteg.hpp:27
Class to represent a coefficient evaluated at quadrature points.
DeviceTensor< 2, double > DeviceMatrix
Definition: dtensor.hpp:143
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
const int MAX_Q1D
Definition: forall.hpp:29
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
double b
Definition: lissajous.cpp:42
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
Enable all above compressions.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
Class representing the storage layout of a FaceQuadratureFunction.
Definition: qspace.hpp:138
int dim
Definition: ex24.cpp:53
const int MAX_D1D
Definition: forall.hpp:28
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Vector data type.
Definition: vector.hpp:60
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3100
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
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