Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
finite_element.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Serac 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 serac {
21 
28 template <int q>
32 
34  SERAC_HOST_DEVICE double weight(int ix, int iy) const { return weights1D[ix] * weights1D[iy]; }
35 
37  SERAC_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>
45 };
46 
54 template <mfem::Geometry::Type g, int q>
56 
58 template <int q>
59 struct batched_jacobian<mfem::Geometry::CUBE, q> {
62 };
63 
65 template <int q>
66 struct batched_jacobian<mfem::Geometry::SQUARE, q> {
69 };
70 
72 template <int q>
73 struct batched_jacobian<mfem::Geometry::TRIANGLE, q> {
75  using type = tensor<double, 2, 2, q*(q + 1) / 2>;
76 };
77 
79 template <int q>
80 struct batched_jacobian<mfem::Geometry::TETRAHEDRON, q> {
82  using type = tensor<double, 3, 3, (q * (q + 1) * (q + 2)) / 6>;
83 };
84 
92 template <mfem::Geometry::Type g, int q>
94 
96 template <int q>
97 struct batched_position<mfem::Geometry::CUBE, q> {
100 };
101 
103 template <int q>
104 struct batched_position<mfem::Geometry::SQUARE, q> {
107 };
108 
110 template <int q>
111 struct batched_position<mfem::Geometry::TRIANGLE, q> {
113  using type = tensor<double, 2, q*(q + 1) / 2>;
114 };
115 
117 template <int q>
118 struct batched_position<mfem::Geometry::TETRAHEDRON, q> {
120  using type = tensor<double, 3, (q * (q + 1) * (q + 2)) / 6>;
121 };
122 
124 template <int q>
125 struct batched_position<mfem::Geometry::SEGMENT, q> {
128 };
129 
140 template <mfem::Geometry::Type g>
142 {
143  if (g == mfem::Geometry::CUBE) {
144  switch (q) {
145  case 1:
146  return 64;
147  case 2:
148  return 16;
149  case 3:
150  return 4;
151  default:
152  return 1;
153  }
154  }
155 
156  if (g == mfem::Geometry::SQUARE) {
157  switch (q) {
158  case 1:
159  return 128;
160  case 2:
161  return 32;
162  case 3:
163  return 16;
164  case 4:
165  return 8;
166  default:
167  return 1;
168  }
169  }
170 }
171 
181 enum class Family
182 {
183  QOI,
184  H1,
185  HCURL,
186  HDIV,
187  L2
188 };
189 
195 template <int p, int c = 1>
196 struct H1 {
197  static constexpr int order = p;
198  static constexpr int components = c;
199  static constexpr Family family = Family::H1;
200 };
201 
207 template <int p, int c = 1>
208 struct Hcurl {
209  static constexpr int order = p;
210  static constexpr int components = c;
211  static constexpr Family family = Family::HCURL;
212 };
213 
219 template <int p, int c = 1>
220 struct L2 {
221  static constexpr int order = p;
222  static constexpr int components = c;
223  static constexpr Family family = Family::L2;
224 };
225 
229 struct QOI {
230  static constexpr int order = 0;
231  static constexpr int components = 1;
232  static constexpr Family family = Family::QOI;
233 };
234 
240  int order;
242 
248  std::tuple<int, int, int> as_tuple() const { return std::tuple<int, int, int>(int(family), order, components); }
249 
253  bool operator<(FunctionSpace other) const { return this->as_tuple() < other.as_tuple(); }
254 };
255 
267 template <Family f, typename T, int q, int dim>
269 {
270  [[maybe_unused]] constexpr int VALUE = 0;
271  [[maybe_unused]] constexpr int DERIVATIVE = 1;
272 
273  for (int k = 0; k < q; k++) {
275  for (int row = 0; row < dim; row++) {
276  for (int col = 0; col < dim; col++) {
277  J[row][col] = jacobians(col, row, k);
278  }
279  }
280 
281  if constexpr (f == Family::H1 || f == Family::L2) {
282  // note: no transformation necessary for the values of H1-field
283  get<DERIVATIVE>(qf_input[k]) = dot(get<DERIVATIVE>(qf_input[k]), inv(J));
284  }
285 
286  if constexpr (f == Family::HCURL) {
287  get<VALUE>(qf_input[k]) = dot(get<VALUE>(qf_input[k]), inv(J));
288  get<DERIVATIVE>(qf_input[k]) = get<DERIVATIVE>(qf_input[k]) / det(J);
289  if constexpr (dim == 3) {
290  get<DERIVATIVE>(qf_input[k]) = dot(get<DERIVATIVE>(qf_input[k]), transpose(J));
291  }
292  }
293  }
294 }
295 
308 template <Family f, typename T, int q, int dim>
310 {
311  [[maybe_unused]] constexpr int SOURCE = 0;
312  [[maybe_unused]] constexpr int FLUX = 1;
313 
314  for (int k = 0; k < q; k++) {
316  for (int row = 0; row < dim; row++) {
317  for (int col = 0; col < dim; col++) {
318  J_T[row][col] = jacobians(row, col, k);
319  }
320  }
321 
322  auto dv = det(J_T);
323 
324  if constexpr (f == Family::H1 || f == Family::L2) {
325  get<SOURCE>(qf_output[k]) = get<SOURCE>(qf_output[k]) * dv;
326  get<FLUX>(qf_output[k]) = dot(get<FLUX>(qf_output[k]), inv(J_T)) * dv;
327  }
328 
329  // note: the flux term here is usually divided by detJ, but
330  // physical_to_parent also multiplies every quadrature-point value by det(J)
331  // so that part cancels out
332  if constexpr (f == Family::HCURL) {
333  get<SOURCE>(qf_output[k]) = dot(get<SOURCE>(qf_output[k]), inv(J_T)) * dv;
334  if constexpr (dim == 3) {
335  get<FLUX>(qf_output[k]) = dot(get<FLUX>(qf_output[k]), transpose(J_T));
336  }
337  }
338 
339  if constexpr (f == Family::QOI) {
340  qf_output[k] = qf_output[k] * dv;
341  }
342  }
343 }
344 
370 template <mfem::Geometry::Type g, typename family>
372 
373 #include "detail/segment_H1.inl"
374 #include "detail/segment_Hcurl.inl"
375 #include "detail/segment_L2.inl"
376 
377 #include "detail/triangle_H1.inl"
378 #include "detail/triangle_L2.inl"
379 
383 
384 #include "detail/tetrahedron_H1.inl"
385 #include "detail/tetrahedron_L2.inl"
386 
387 #include "detail/hexahedron_H1.inl"
389 #include "detail/hexahedron_L2.inl"
390 
391 #include "detail/qoi.inl"
392 
393 } // namespace serac
#define SERAC_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
Definition: accelerator.hpp:38
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: serac.cpp:36
constexpr SERAC_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
constexpr SERAC_HOST_DEVICE int elements_per_block(int q)
this function returns information about how many elements should be processed by a single thread bloc...
SERAC_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 SERAC_HOST_DEVICE auto det(const isotropic_tensor< T, m, m > &I)
compute the determinant of an isotropic tensor
constexpr SERAC_HOST_DEVICE auto inv(const isotropic_tensor< T, m, m > &I)
return the inverse of an isotropic tensor
constexpr SERAC_HOST_DEVICE auto transpose(const isotropic_tensor< T, m, m > &I)
return the transpose of an isotropic tensor
Family
Element conformity.
SERAC_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,...
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
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
Family family
either H1, Hcurl, L2
H1 elements of order p.
static constexpr int order
the polynomial order of the elements
static constexpr Family family
the family of the basis functions
static constexpr int components
the number of components at each node
H(curl) elements of order p.
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
Discontinuous elements of order p.
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
"Quantity of Interest" elements (i.e. elements with a single shape function, 1)
static constexpr int order
the polynomial order of the elements
static constexpr Family family
the family of the basis functions
static constexpr int components
the number of components at each node
a convenience class for generating information about tensor product integration rules from the underl...
SERAC_HOST_DEVICE double weight(int ix, int iy) const
return the quadrature weight for a quadrilateral
SERAC_HOST_DEVICE double weight(int ix, int iy, int iz) const
return the quadrature weight for a hexahedron
tensor< double, q > points1D
the abscissae of the underlying 1D quadrature rule
tensor< double, q > weights1D
the weights of the underlying 1D quadrature rule
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.