Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
finite_element.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
15 #include "tuple.hpp"
16 #include "tensor.hpp"
17 #include "geometry.hpp"
18 #include "polynomials.hpp"
19 
20 namespace smith {
21 
28 template <int q>
32 
34  SMITH_HOST_DEVICE double weight(int ix, int iy) const { return weights1D[ix] * weights1D[iy]; }
35 
37  SMITH_HOST_DEVICE double weight(int ix, int iy, int iz) const
38  {
39  return weights1D[ix] * weights1D[iy] * weights1D[iz];
40  }
41 };
42 
43 template <auto val>
44 struct CompileTimeValue {};
45 
53 template <mfem::Geometry::Type g, int q>
55 
57 template <int q>
58 struct batched_jacobian<mfem::Geometry::CUBE, q> {
61 };
62 
64 template <int q>
65 struct batched_jacobian<mfem::Geometry::SQUARE, q> {
68 };
69 
71 template <int q>
72 struct batched_jacobian<mfem::Geometry::TRIANGLE, q> {
74  using type = tensor<double, 2, 2, q*(q + 1) / 2>;
75 };
76 
78 template <int q>
79 struct batched_jacobian<mfem::Geometry::TETRAHEDRON, q> {
81  using type = tensor<double, 3, 3, (q * (q + 1) * (q + 2)) / 6>;
82 };
83 
91 template <mfem::Geometry::Type g, int q>
93 
95 template <int q>
96 struct batched_position<mfem::Geometry::CUBE, q> {
99 };
100 
102 template <int q>
103 struct batched_position<mfem::Geometry::SQUARE, q> {
106 };
107 
109 template <int q>
110 struct batched_position<mfem::Geometry::TRIANGLE, q> {
112  using type = tensor<double, 2, q*(q + 1) / 2>;
113 };
114 
116 template <int q>
117 struct batched_position<mfem::Geometry::TETRAHEDRON, q> {
119  using type = tensor<double, 3, (q * (q + 1) * (q + 2)) / 6>;
120 };
121 
123 template <int q>
124 struct batched_position<mfem::Geometry::SEGMENT, q> {
127 };
128 
139 template <mfem::Geometry::Type g>
141 {
142  if (g == mfem::Geometry::CUBE) {
143  switch (q) {
144  case 1:
145  return 64;
146  case 2:
147  return 16;
148  case 3:
149  return 4;
150  default:
151  return 1;
152  }
153  }
154 
155  if (g == mfem::Geometry::SQUARE) {
156  switch (q) {
157  case 1:
158  return 128;
159  case 2:
160  return 32;
161  case 3:
162  return 16;
163  case 4:
164  return 8;
165  default:
166  return 1;
167  }
168  }
169 }
170 
180 enum class Family
181 {
182  QOI,
183  H1,
184  HCURL,
185  HDIV,
186  L2
187 };
188 
194 template <int p, int c = 1>
195 struct H1 {
196  static constexpr int order = p;
197  static constexpr int components = c;
198  static constexpr Family family = Family::H1;
199 };
200 
206 template <int p, int c = 1>
207 struct Hcurl {
208  static constexpr int order = p;
209  static constexpr int components = c;
210  static constexpr Family family = Family::HCURL;
211 };
212 
218 template <int p, int c = 1>
219 struct L2 {
220  static constexpr int order = p;
221  static constexpr int components = c;
222  static constexpr Family family = Family::L2;
223 };
224 
228 struct QOI {
229  static constexpr int order = 0;
230  static constexpr int components = 1;
231  static constexpr Family family = Family::QOI;
232 };
233 
239  int order;
241 
247  std::tuple<int, int, int> as_tuple() const { return std::tuple<int, int, int>(int(family), order, components); }
248 
252  bool operator<(FunctionSpace other) const { return this->as_tuple() < other.as_tuple(); }
253 };
254 
266 template <Family f, typename T, int q, int dim>
268 {
269  [[maybe_unused]] constexpr int VALUE = 0;
270  [[maybe_unused]] constexpr int DERIVATIVE = 1;
271 
272  for (int k = 0; k < q; k++) {
274  for (int row = 0; row < dim; row++) {
275  for (int col = 0; col < dim; col++) {
276  J[row][col] = jacobians(col, row, k);
277  }
278  }
279 
280  if constexpr (f == Family::H1 || f == Family::L2) {
281  // note: no transformation necessary for the values of H1-field
282  get<DERIVATIVE>(qf_input[k]) = dot(get<DERIVATIVE>(qf_input[k]), inv(J));
283  }
284 
285  if constexpr (f == Family::HCURL) {
286  get<VALUE>(qf_input[k]) = dot(get<VALUE>(qf_input[k]), inv(J));
287  get<DERIVATIVE>(qf_input[k]) = get<DERIVATIVE>(qf_input[k]) / det(J);
288  if constexpr (dim == 3) {
289  get<DERIVATIVE>(qf_input[k]) = dot(get<DERIVATIVE>(qf_input[k]), transpose(J));
290  }
291  }
292  }
293 }
294 
307 template <Family f, typename T, int q, int dim>
309 {
310  [[maybe_unused]] constexpr int SOURCE = 0;
311  [[maybe_unused]] constexpr int FLUX = 1;
312 
313  for (int k = 0; k < q; k++) {
315  for (int row = 0; row < dim; row++) {
316  for (int col = 0; col < dim; col++) {
317  J_T[row][col] = jacobians(row, col, k);
318  }
319  }
320 
321  auto dv = det(J_T);
322 
323  if constexpr (f == Family::H1 || f == Family::L2) {
324  get<SOURCE>(qf_output[k]) = get<SOURCE>(qf_output[k]) * dv;
325  get<FLUX>(qf_output[k]) = dot(get<FLUX>(qf_output[k]), inv(J_T)) * dv;
326  }
327 
328  // note: the flux term here is usually divided by detJ, but
329  // physical_to_parent also multiplies every quadrature-point value by det(J)
330  // so that part cancels out
331  if constexpr (f == Family::HCURL) {
332  get<SOURCE>(qf_output[k]) = dot(get<SOURCE>(qf_output[k]), inv(J_T)) * dv;
333  if constexpr (dim == 3) {
334  get<FLUX>(qf_output[k]) = dot(get<FLUX>(qf_output[k]), transpose(J_T));
335  }
336  }
337 
338  if constexpr (f == Family::QOI) {
339  qf_output[k] = qf_output[k] * dv;
340  }
341  }
342 }
343 
369 template <mfem::Geometry::Type g, typename family>
371 
372 #include "detail/segment_H1.inl"
373 #include "detail/segment_Hcurl.inl"
374 #include "detail/segment_L2.inl"
375 
376 #include "detail/triangle_H1.inl"
377 #include "detail/triangle_L2.inl"
378 
382 
383 #include "detail/tetrahedron_H1.inl"
384 #include "detail/tetrahedron_L2.inl"
385 
386 #include "detail/hexahedron_H1.inl"
388 #include "detail/hexahedron_L2.inl"
389 
390 #include "detail/qoi.inl"
391 
392 } // namespace smith
#define SMITH_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc or amdclang and does nothing on ...
Definition: accelerator.hpp:37
Specialization of finite_element for H1 on hexahedron geometry.
Specialization of finite_element for Hcurl on hexahedron geometry.
Specialization of finite_element for L2 on hexahedron geometry.
Accelerator functionality.
Definition: smith.cpp:36
constexpr SMITH_HOST_DEVICE int elements_per_block(int q)
this function returns information about how many elements should be processed by a single thread bloc...
SMITH_HOST_DEVICE void parent_to_physical(tensor< T, q > &qf_input, const tensor< double, dim, dim, q > &jacobians)
transform information in the parent space (i.e. values and derivatives w.r.t {xi, eta,...
Family
Element conformity.
constexpr SMITH_HOST_DEVICE auto inv(const isotropic_tensor< T, m, m > &I)
return the inverse of an isotropic tensor
SMITH_HOST_DEVICE void physical_to_parent(tensor< T, q > &qf_output, const tensor< double, dim, dim, q > &jacobians)
transform information in the physical space (i.e. sources and fluxes w.r.t {x, y, z}) back to the par...
constexpr SMITH_HOST_DEVICE auto transpose(const isotropic_tensor< T, m, m > &I)
return the transpose of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto det(const isotropic_tensor< T, m, m > &I)
compute the determinant of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto dot(const isotropic_tensor< S, m, m > &I, const tensor< T, m, n... > &A)
dot product between an isotropic and (nonisotropic) tensor
Definitions of 1D quadrature weights and node locations and polynomial basis functions.
Specialization of finite_element for expressing quantities of interest on any geometry.
Specialization of finite_element for H1 on quadrilateral geometry.
Specialization of finite_element for Hcurl on quadrilateral geometry.
Specialization of finite_element for L2 on quadrilateral geometry.
Specialization of finite_element for H1 on segment geometry.
Specialization of finite_element for Hcurl on segment geometry.
Specialization of finite_element for L2 on segment geometry.
a small POD class for tracking function space metadata
std::tuple< int, int, int > as_tuple() const
return the data contained in this struct as a tuple
Family family
either H1, Hcurl, L2
bool operator<(FunctionSpace other) const
defines an ordering over FunctionSpaces, to enable use in containers like std::map
int order
polynomial order
int components
how many values are stored at each node
H1 elements of order p.
static constexpr Family family
the family of the basis functions
static constexpr int order
the polynomial order of the elements
static constexpr int components
the number of components at each node
H(curl) elements of order p.
static constexpr Family family
the family of the basis functions
static constexpr int components
the number of components at each node
static constexpr int order
the polynomial order of the elements
Discontinuous elements of order p.
static constexpr int components
the number of components at each node
static constexpr Family family
the family of the basis functions
static constexpr int order
the polynomial order of the elements
"Quantity of Interest" elements (i.e. elements with a single shape function, 1)
static constexpr int components
the number of components at each node
static constexpr int order
the polynomial order of the elements
static constexpr Family family
the family of the basis functions
a convenience class for generating information about tensor product integration rules from the underl...
tensor< double, q > weights1D
the weights of the underlying 1D quadrature rule
SMITH_HOST_DEVICE double weight(int ix, int iy) const
return the quadrature weight for a quadrilateral
tensor< double, q > points1D
the abscissae of the underlying 1D quadrature rule
SMITH_HOST_DEVICE double weight(int ix, int iy, int iz) const
return the quadrature weight for a hexahedron
this struct is used to look up mfem's memory layout of the quadrature point jacobian matrices
this struct is used to look up mfem's memory layout of the quadrature point position vectors
Template prototype for finite element implementations.
Implementation of the tensor class used by Functional.
Specialization of finite_element for H1 on tetrahedron geometry.
Specialization of finite_element for L2 on tetrahedron geometry.
Specialization of finite_element for H1 on triangle geometry.
Specialization of finite_element for L2 on triangle geometry.
Implements a std::tuple-like object that works in CUDA kernels.