MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lininteg_domain.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 DLFEvalAssemble2D(const int vdim, const int ne, const int d, const int q,
21  const int map_type, const int *markers, const double *b,
22  const double *detj, const double *weights,
23  const Vector &coeff, double *y)
24 {
25  const auto F = coeff.Read();
26  const auto M = Reshape(markers, ne);
27  const auto B = Reshape(b, q, d);
28  const auto DETJ = Reshape(detj, q, q, ne);
29  const auto W = Reshape(weights, q, q);
30  const bool cst = coeff.Size() == vdim;
31  const auto C = cst ? Reshape(F,vdim,1,1,1) : Reshape(F,vdim,q,q,ne);
32  auto Y = Reshape(y, d,d, vdim, ne);
33 
34  MFEM_FORALL_2D(e, ne, q, q, 1,
35  {
36  if (M(e) == 0) { return; } // ignore
37 
38  constexpr int Q = T_Q1D ? T_Q1D : MAX_Q1D;
39  constexpr int D = T_D1D ? T_D1D : MAX_D1D;
40 
41  MFEM_SHARED double sBt[Q*D];
42  MFEM_SHARED double sQQ[Q*Q];
43  MFEM_SHARED double sQD[Q*D];
44 
45  const DeviceMatrix Bt(sBt, d, q);
46  kernels::internal::LoadB<D,Q>(d, q, B, sBt);
47 
48  const DeviceMatrix QQ(sQQ, q, q);
49  const DeviceMatrix QD(sQD, q, d);
50 
51  for (int c = 0; c < vdim; ++c)
52  {
53  const double cst_val = C(c,0,0,0);
54  MFEM_FOREACH_THREAD(x,x,q)
55  {
56  MFEM_FOREACH_THREAD(y,y,q)
57  {
58  const double detJ = (map_type == FiniteElement::VALUE) ? DETJ(x,y,e) : 1.0;
59  const double coeff_val = cst ? cst_val : C(c,x,y,e);
60  QQ(y,x) = W(x,y) * coeff_val * detJ;
61  }
62  }
63  MFEM_SYNC_THREAD;
64  MFEM_FOREACH_THREAD(qy,y,q)
65  {
66  MFEM_FOREACH_THREAD(dx,x,d)
67  {
68  double u = 0.0;
69  for (int qx = 0; qx < q; ++qx) { u += QQ(qy,qx) * Bt(dx,qx); }
70  QD(qy,dx) = u;
71  }
72  }
73  MFEM_SYNC_THREAD;
74  MFEM_FOREACH_THREAD(dy,y,d)
75  {
76  MFEM_FOREACH_THREAD(dx,x,d)
77  {
78  double u = 0.0;
79  for (int qy = 0; qy < q; ++qy) { u += QD(qy,dx) * Bt(dy,qy); }
80  Y(dx,dy,c,e) += u;
81  }
82  }
83  MFEM_SYNC_THREAD;
84  }
85  });
86 }
87 
88 template<int T_D1D = 0, int T_Q1D = 0> static
89 void DLFEvalAssemble3D(const int vdim, const int ne, const int d, const int q,
90  const int map_type, const int *markers, const double *b,
91  const double *detj, const double *weights,
92  const Vector &coeff, double *y)
93 {
94  const auto F = coeff.Read();
95  const auto M = Reshape(markers, ne);
96  const auto B = Reshape(b, q,d);
97  const auto DETJ = Reshape(detj, q, q, q, ne);
98  const auto W = Reshape(weights, q,q,q);
99  const bool cst_coeff = coeff.Size() == vdim;
100  const auto C = cst_coeff ? Reshape(F,vdim,1,1,1,1):Reshape(F,vdim,q,q,q,ne);
101 
102  auto Y = Reshape(y, d,d,d, vdim, ne);
103 
104  MFEM_FORALL_2D(e, ne, q, q, 1,
105  {
106  if (M(e) == 0) { return; } // ignore
107 
108  constexpr int Q = T_Q1D ? T_Q1D : MAX_Q1D;
109  constexpr int D = T_D1D ? T_D1D : MAX_D1D;
110 
111  double u[D];
112 
113  MFEM_SHARED double sBt[Q*D];
114  const DeviceMatrix Bt(sBt, d,q);
115  kernels::internal::LoadB<D,Q>(d,q,B,sBt);
116 
117  MFEM_SHARED double sQQQ[Q*Q*Q];
118  const DeviceCube QQQ(sQQQ, q,q,q);
119 
120  for (int c = 0; c < vdim; ++c)
121  {
122  const double cst_val = C(c,0,0,0,0);
123  MFEM_FOREACH_THREAD(x,x,q)
124  {
125  MFEM_FOREACH_THREAD(y,y,q)
126  {
127  for (int z = 0; z < q; ++z)
128  {
129  const double detJ = (map_type == FiniteElement::VALUE) ? DETJ(x,y,z,e) : 1.0;
130  const double coeff_val = cst_coeff ? cst_val : C(c,x,y,z,e);
131  QQQ(z,y,x) = W(x,y,z) * coeff_val * detJ;
132  }
133  }
134  }
135  MFEM_SYNC_THREAD;
136  MFEM_FOREACH_THREAD(qx,x,q)
137  {
138  MFEM_FOREACH_THREAD(qy,y,q)
139  {
140  for (int dz = 0; dz < d; ++dz) { u[dz] = 0.0; }
141  for (int qz = 0; qz < q; ++qz)
142  {
143  const double ZYX = QQQ(qz,qy,qx);
144  for (int dz = 0; dz < d; ++dz) { u[dz] += ZYX * Bt(dz,qz); }
145  }
146  for (int dz = 0; dz < d; ++dz) { QQQ(dz,qy,qx) = u[dz]; }
147  }
148  }
149  MFEM_SYNC_THREAD;
150  MFEM_FOREACH_THREAD(dz,y,d)
151  {
152  MFEM_FOREACH_THREAD(qx,x,q)
153  {
154  for (int dy = 0; dy < d; ++dy) { u[dy] = 0.0; }
155  for (int qy = 0; qy < q; ++qy)
156  {
157  const double zYX = QQQ(dz,qy,qx);
158  for (int dy = 0; dy < d; ++dy) { u[dy] += zYX * Bt(dy,qy); }
159  }
160  for (int dy = 0; dy < d; ++dy) { QQQ(dz,dy,qx) = u[dy]; }
161  }
162  }
163  MFEM_SYNC_THREAD;
164  MFEM_FOREACH_THREAD(dz,y,d)
165  {
166  MFEM_FOREACH_THREAD(dy,x,d)
167  {
168  for (int dx = 0; dx < d; ++dx) { u[dx] = 0.0; }
169  for (int qx = 0; qx < q; ++qx)
170  {
171  const double zyX = QQQ(dz,dy,qx);
172  for (int dx = 0; dx < d; ++dx) { u[dx] += zyX * Bt(dx,qx); }
173  }
174  for (int dx = 0; dx < d; ++dx) { Y(dx,dy,dz,c,e) += u[dx]; }
175  }
176  }
177  MFEM_SYNC_THREAD;
178  }
179  });
180 }
181 
182 static void DLFEvalAssemble(const FiniteElementSpace &fes,
183  const IntegrationRule *ir,
184  const Array<int> &markers,
185  const Vector &coeff,
186  Vector &y)
187 {
188  Mesh *mesh = fes.GetMesh();
189  const int dim = mesh->Dimension();
190  const FiniteElement &el = *fes.GetFE(0);
192  const DofToQuad &maps = el.GetDofToQuad(*ir, DofToQuad::TENSOR);
193  const int d = maps.ndof, q = maps.nqpt;
194  constexpr int flags = GeometricFactors::DETERMINANTS;
195  const GeometricFactors *geom = mesh->GetGeometricFactors(*ir, flags, mt);
196  const int map_type = fes.GetFE(0)->GetMapType();
197  decltype(&DLFEvalAssemble2D<>) ker =
198  dim == 2 ? DLFEvalAssemble2D<> : DLFEvalAssemble3D<>;
199 
200  if (dim==2)
201  {
202  if (d==1 && q==1) { ker=DLFEvalAssemble2D<1,1>; }
203  if (d==2 && q==2) { ker=DLFEvalAssemble2D<2,2>; }
204  if (d==3 && q==3) { ker=DLFEvalAssemble2D<3,3>; }
205  if (d==4 && q==4) { ker=DLFEvalAssemble2D<4,4>; }
206  if (d==5 && q==5) { ker=DLFEvalAssemble2D<5,5>; }
207  if (d==2 && q==3) { ker=DLFEvalAssemble2D<2,3>; }
208  if (d==3 && q==4) { ker=DLFEvalAssemble2D<3,4>; }
209  if (d==4 && q==5) { ker=DLFEvalAssemble2D<4,5>; }
210  if (d==5 && q==6) { ker=DLFEvalAssemble2D<5,6>; }
211  }
212 
213  if (dim==3)
214  {
215  if (d==1 && q==1) { ker=DLFEvalAssemble3D<1,1>; }
216  if (d==2 && q==2) { ker=DLFEvalAssemble3D<2,2>; }
217  if (d==3 && q==3) { ker=DLFEvalAssemble3D<3,3>; }
218  if (d==4 && q==4) { ker=DLFEvalAssemble3D<4,4>; }
219  if (d==5 && q==5) { ker=DLFEvalAssemble3D<5,5>; }
220  if (d==2 && q==3) { ker=DLFEvalAssemble3D<2,3>; }
221  if (d==3 && q==4) { ker=DLFEvalAssemble3D<3,4>; }
222  if (d==4 && q==5) { ker=DLFEvalAssemble3D<4,5>; }
223  if (d==5 && q==6) { ker=DLFEvalAssemble3D<5,6>; }
224  }
225 
226  MFEM_VERIFY(ker, "No kernel ndof " << d << " nqpt " << q);
227 
228  const int vdim = fes.GetVDim();
229  const int ne = fes.GetMesh()->GetNE();
230  const int *M = markers.Read();
231  const double *B = maps.B.Read();
232  const double *detJ = geom->detJ.Read();
233  const double *W = ir->GetWeights().Read();
234  double *Y = y.ReadWrite();
235  ker(vdim, ne, d, q, map_type, M, B, detJ, W, coeff, Y);
236 }
237 
239  const Array<int> &markers,
240  Vector &b)
241 {
242  const FiniteElement &fe = *fes.GetFE(0);
243  const int qorder = oa * fe.GetOrder() + ob;
244  const Geometry::Type gtype = fe.GetGeomType();
245  const IntegrationRule *ir = IntRule ? IntRule : &IntRules.Get(gtype, qorder);
246 
247  QuadratureSpace qs(*fes.GetMesh(), *ir);
249  DLFEvalAssemble(fes, ir, markers, coeff, b);
250 }
251 
253  const Array<int> &markers,
254  Vector &b)
255 {
256  const FiniteElement &fe = *fes.GetFE(0);
257  const int qorder = 2 * fe.GetOrder();
258  const Geometry::Type gtype = fe.GetGeomType();
259  const IntegrationRule *ir = IntRule ? IntRule : &IntRules.Get(gtype, qorder);
260 
261  QuadratureSpace qs(*fes.GetMesh(), *ir);
263  DLFEvalAssemble(fes, ir, markers, coeff, b);
264 }
265 
266 } // 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.
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
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
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.
DeviceTensor< 3, double > DeviceCube
Definition: dtensor.hpp:146
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
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
int dim
Definition: ex24.cpp:53
const int MAX_D1D
Definition: forall.hpp:28
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
Vector data type.
Definition: vector.hpp:60
Class representing the storage layout of a QuadratureFunction.
Definition: qspace.hpp:92
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
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