MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
det.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 "../quadinterpolator.hpp"
13 #include "../../general/forall.hpp"
14 #include "../../linalg/dtensor.hpp"
15 #include "../../fem/kernels.hpp"
16 #include "../../linalg/kernels.hpp"
17 
18 using namespace mfem;
19 
20 namespace mfem
21 {
22 
23 namespace internal
24 {
25 
26 namespace quadrature_interpolator
27 {
28 
29 static void Det1D(const int NE,
30  const double *g,
31  const double *x,
32  double *y,
33  const int d1d,
34  const int q1d)
35 {
36  const auto G = Reshape(g, q1d, d1d);
37  const auto X = Reshape(x, d1d, NE);
38 
39  auto Y = Reshape(y, q1d, NE);
40 
41  MFEM_FORALL(e, NE,
42  {
43  for (int q = 0; q < q1d; q++)
44  {
45  double u = 0.0;
46  for (int d = 0; d < d1d; d++)
47  {
48  u += G(q, d) * X(d, e);
49  }
50  Y(q, e) = u;
51  }
52  });
53 }
54 
55 template<int T_D1D = 0, int T_Q1D = 0, int MAX_D1D = 0, int MAX_Q1D = 0>
56 static void Det2D(const int NE,
57  const double *b,
58  const double *g,
59  const double *x,
60  double *y,
61  const int vdim = 1,
62  const int d1d = 0,
63  const int q1d = 0)
64 {
65  constexpr int DIM = 2;
66  static constexpr int NBZ = 1;
67 
68  const int D1D = T_D1D ? T_D1D : d1d;
69  const int Q1D = T_Q1D ? T_Q1D : q1d;
70 
71  const auto B = Reshape(b, Q1D, D1D);
72  const auto G = Reshape(g, Q1D, D1D);
73  const auto X = Reshape(x, D1D, D1D, DIM, NE);
74  auto Y = Reshape(y, Q1D, Q1D, NE);
75 
76  MFEM_FORALL_2D(e, NE, Q1D, Q1D, NBZ,
77  {
78  constexpr int MQ1 = T_Q1D ? T_Q1D : MAX_Q1D;
79  constexpr int MD1 = T_D1D ? T_D1D : MAX_D1D;
80  const int D1D = T_D1D ? T_D1D : d1d;
81  const int Q1D = T_Q1D ? T_Q1D : q1d;
82 
83  MFEM_SHARED double BG[2][MQ1*MD1];
84  MFEM_SHARED double XY[2][NBZ][MD1*MD1];
85  MFEM_SHARED double DQ[4][NBZ][MD1*MQ1];
86  MFEM_SHARED double QQ[4][NBZ][MQ1*MQ1];
87 
88  kernels::internal::LoadX<MD1,NBZ>(e,D1D,X,XY);
89  kernels::internal::LoadBG<MD1,MQ1>(D1D,Q1D,B,G,BG);
90 
91  kernels::internal::GradX<MD1,MQ1,NBZ>(D1D,Q1D,BG,XY,DQ);
92  kernels::internal::GradY<MD1,MQ1,NBZ>(D1D,Q1D,BG,DQ,QQ);
93 
94  MFEM_FOREACH_THREAD(qy,y,Q1D)
95  {
96  MFEM_FOREACH_THREAD(qx,x,Q1D)
97  {
98  double J[4];
99  kernels::internal::PullGrad<MQ1,NBZ>(Q1D,qx,qy,QQ,J);
100  Y(qx,qy,e) = kernels::Det<2>(J);
101  }
102  }
103  });
104 }
105 
106 template<int T_D1D = 0, int T_Q1D = 0, int MAX_D1D = 0, int MAX_Q1D = 0,
107  bool SMEM = true>
108 static void Det3D(const int NE,
109  const double *b,
110  const double *g,
111  const double *x,
112  double *y,
113  const int vdim = 1,
114  const int d1d = 0,
115  const int q1d = 0,
116  Vector *d_buff = nullptr) // used only with SMEM = false
117 {
118  constexpr int DIM = 3;
119  static constexpr int MQ1 = T_Q1D ? T_Q1D : MAX_Q1D;
120  static constexpr int MD1 = T_D1D ? T_D1D : MAX_D1D;
121  static constexpr int MDQ = MQ1 > MD1 ? MQ1 : MD1;
122  static constexpr int MSZ = MDQ * MDQ * MDQ * 9;
123  static constexpr int GRID = SMEM ? 0 : 128;
124 
125  const int D1D = T_D1D ? T_D1D : d1d;
126  const int Q1D = T_Q1D ? T_Q1D : q1d;
127 
128  const auto B = Reshape(b, Q1D, D1D);
129  const auto G = Reshape(g, Q1D, D1D);
130  const auto X = Reshape(x, D1D, D1D, D1D, DIM, NE);
131  auto Y = Reshape(y, Q1D, Q1D, Q1D, NE);
132 
133  double *GM = nullptr;
134  if (!SMEM)
135  {
136  d_buff->SetSize(2*MSZ*GRID);
137  GM = d_buff->Write();
138  }
139 
140  MFEM_FORALL_3D_GRID(e, NE, Q1D, Q1D, Q1D, GRID,
141  {
142  const int bid = MFEM_BLOCK_ID(x);
143  MFEM_SHARED double BG[2][MQ1*MD1];
144  MFEM_SHARED double SM0[SMEM?MSZ:1];
145  MFEM_SHARED double SM1[SMEM?MSZ:1];
146  double *lm0 = SMEM ? SM0 : GM + MSZ*bid;
147  double *lm1 = SMEM ? SM1 : GM + MSZ*(GRID+bid);
148  double (*DDD)[MD1*MD1*MD1] = (double (*)[MD1*MD1*MD1]) (lm0);
149  double (*DDQ)[MD1*MD1*MQ1] = (double (*)[MD1*MD1*MQ1]) (lm1);
150  double (*DQQ)[MD1*MQ1*MQ1] = (double (*)[MD1*MQ1*MQ1]) (lm0);
151  double (*QQQ)[MQ1*MQ1*MQ1] = (double (*)[MQ1*MQ1*MQ1]) (lm1);
152 
153  kernels::internal::LoadX<MD1>(e,D1D,X,DDD);
154  kernels::internal::LoadBG<MD1,MQ1>(D1D,Q1D,B,G,BG);
155 
156  kernels::internal::GradX<MD1,MQ1>(D1D,Q1D,BG,DDD,DDQ);
157  kernels::internal::GradY<MD1,MQ1>(D1D,Q1D,BG,DDQ,DQQ);
158  kernels::internal::GradZ<MD1,MQ1>(D1D,Q1D,BG,DQQ,QQQ);
159 
160  MFEM_FOREACH_THREAD(qz,z,Q1D)
161  {
162  MFEM_FOREACH_THREAD(qy,y,Q1D)
163  {
164  MFEM_FOREACH_THREAD(qx,x,Q1D)
165  {
166  double J[9];
167  kernels::internal::PullGrad<MQ1>(Q1D, qx,qy,qz, QQQ, J);
168  Y(qx,qy,qz,e) = kernels::Det<3>(J);
169  }
170  }
171  }
172  });
173 }
174 
175 // Tensor-product evaluation of quadrature point determinants: dispatch
176 // function.
177 void TensorDeterminants(const int NE,
178  const int vdim,
179  const DofToQuad &maps,
180  const Vector &e_vec,
181  Vector &q_det,
182  Vector &d_buff)
183 {
184  if (NE == 0) { return; }
185  const int dim = maps.FE->GetDim();
186  const int D1D = maps.ndof;
187  const int Q1D = maps.nqpt;
188  const double *B = maps.B.Read();
189  const double *G = maps.G.Read();
190  const double *X = e_vec.Read();
191  double *Y = q_det.Write();
192 
193  const int id = (vdim<<8) | (D1D<<4) | Q1D;
194 
195  if (dim == 1)
196  {
197  MFEM_VERIFY(D1D <= MAX_D1D, "Orders higher than " << MAX_D1D-1
198  << " are not supported!");
199  MFEM_VERIFY(Q1D <= MAX_Q1D, "Quadrature rules with more than "
200  << MAX_Q1D << " 1D points are not supported!");
201  Det1D(NE, G, X, Y, D1D, Q1D);
202  return;
203  }
204  if (dim == 2)
205  {
206  switch (id)
207  {
208  case 0x222: return Det2D<2,2>(NE,B,G,X,Y);
209  case 0x223: return Det2D<2,3>(NE,B,G,X,Y);
210  case 0x224: return Det2D<2,4>(NE,B,G,X,Y);
211  case 0x226: return Det2D<2,6>(NE,B,G,X,Y);
212  case 0x234: return Det2D<3,4>(NE,B,G,X,Y);
213  case 0x236: return Det2D<3,6>(NE,B,G,X,Y);
214  case 0x244: return Det2D<4,4>(NE,B,G,X,Y);
215  case 0x246: return Det2D<4,6>(NE,B,G,X,Y);
216  case 0x256: return Det2D<5,6>(NE,B,G,X,Y);
217  default:
218  {
219  constexpr int MD = MAX_D1D;
220  constexpr int MQ = MAX_Q1D;
221  MFEM_VERIFY(D1D <= MD, "Orders higher than " << MD-1
222  << " are not supported!");
223  MFEM_VERIFY(Q1D <= MQ, "Quadrature rules with more than "
224  << MQ << " 1D points are not supported!");
225  Det2D<0,0,MD,MQ>(NE,B,G,X,Y,vdim,D1D,Q1D);
226  return;
227  }
228  }
229  }
230  if (dim == 3)
231  {
232  switch (id)
233  {
234  case 0x324: return Det3D<2,4>(NE,B,G,X,Y);
235  case 0x333: return Det3D<3,3>(NE,B,G,X,Y);
236  case 0x335: return Det3D<3,5>(NE,B,G,X,Y);
237  case 0x336: return Det3D<3,6>(NE,B,G,X,Y);
238  default:
239  {
240  constexpr int MD = 6;
241  constexpr int MQ = 6;
242  // Highest orders that fit in shared memory
243  if (D1D <= MD && Q1D <= MQ)
244  { return Det3D<0,0,MD,MQ>(NE,B,G,X,Y,vdim,D1D,Q1D); }
245  // Last fall-back will use global memory
246  return Det3D<0,0,MAX_D1D,MAX_Q1D,false>(
247  NE,B,G,X,Y,vdim,D1D,Q1D,&d_buff);
248  }
249  }
250  }
251  MFEM_ABORT("Kernel " << std::hex << id << std::dec << " not supported yet");
252 }
253 
254 } // namespace quadrature_interpolator
255 
256 } // namespace internal
257 
258 } // namespace mfem
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:174
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:170
constexpr int DIM
const int MAX_Q1D
Definition: forall.hpp:29
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:457
double b
Definition: lissajous.cpp:42
class FiniteElement * FE
The FiniteElement that created and owns this object.
Definition: fe_base.hpp:141
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
MFEM_HOST_DEVICE double Det3D(DeviceMatrix &J)
Definition: lor_util.hpp:188
MFEM_HOST_DEVICE double Det2D(DeviceMatrix &J)
Definition: lor_util.hpp:183
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:185
int dim
Definition: ex24.cpp:53
const int MAX_D1D
Definition: forall.hpp:28
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:206
Vector data type.
Definition: vector.hpp:60
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
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