Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
geometric_factors.cpp
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 
7 #include "smith/numerics/functional/geometric_factors.hpp"
8 
9 #include <iostream>
10 #include <vector>
11 
13 
14 namespace smith {
15 
26 template <int Q, mfem::Geometry::Type geom, typename function_space>
27 void compute_geometric_factors(mfem::Vector& positions_q, mfem::Vector& jacobians_q, const mfem::Vector& positions_e,
28  const std::vector<int>& elements)
29 {
30  static constexpr TensorProductQuadratureRule<Q> rule{};
31 
32  constexpr int spatial_dim = function_space::components;
33  constexpr int geometry_dim = dimension_of(geom);
34  constexpr int qpts_per_elem = num_quadrature_points(geom, Q);
35 
36  using element_type = finite_element<geom, function_space>;
37  using position_type = tensor<double, spatial_dim, qpts_per_elem>;
39 
40  auto X_q = reinterpret_cast<position_type*>(positions_q.ReadWrite());
41  auto J_q = reinterpret_cast<jacobian_type*>(jacobians_q.ReadWrite());
42  auto X = reinterpret_cast<const typename element_type::dof_type*>(positions_e.Read());
43 
44  std::size_t num_elements = elements.size();
45 
46  // for each element in the domain
47  for (uint32_t e = 0; e < num_elements; e++) {
48  // load the positions for the nodes in this element
49  auto X_e = X[e];
50 
51  // calculate the values and derivatives (w.r.t. xi) of X at each quadrature point
52  auto quadrature_values = element_type::interpolate(X_e, rule);
53 
54  // mfem wants to store this data in a different layout, so we have to transpose it
55  for (int q = 0; q < qpts_per_elem; q++) {
56  auto [value, gradient] = quadrature_values[q];
57  for (int i = 0; i < spatial_dim; i++) {
58  X_q[e](i, q) = value[i];
59  if constexpr (std::is_same_v<decltype(value), decltype(gradient)>) {
60  J_q[e](0, i, q) = gradient[i];
61  }
62  if constexpr (!std::is_same_v<decltype(value), decltype(gradient)>) {
63  for (int j = 0; j < geometry_dim; j++) {
64  J_q[e](j, i, q) = gradient(i, j);
65  }
66  }
67  }
68  }
69  }
70 }
71 
72 GeometricFactors::GeometricFactors(const Domain& domain, int q, mfem::Geometry::Type geom)
73 {
74  const mfem::ParGridFunction* nodes = static_cast<const mfem::ParGridFunction*>(domain.mesh_.GetNodes());
75  mfem::ParFiniteElementSpace* fes = nodes->ParFESpace();
76  // const mfem::GridFunction* nodes = domain.mesh_.GetNodes();
77  // const mfem::FiniteElementSpace* fes = nodes->FESpace();
78 
79  const std::vector<int>& element_ids = domain.get_mfem_ids(geom);
80 
81  auto restriction = smith::ElementRestriction(fes, geom, element_ids);
82  mfem::Vector X_e(int(restriction.ESize()));
83  restriction.Gather(*nodes, X_e);
84 
85  // assumes all elements are the same order
86  int p = fes->GetElementOrder(0);
87 
88  int spatial_dim = domain.mesh_.SpaceDimension();
89  int geometry_dim = dimension_of(geom);
90  int qpts_per_elem = num_quadrature_points(geom, q);
91 
92  num_elements = element_ids.size();
93 
94  X = mfem::Vector(int(num_elements) * qpts_per_elem * spatial_dim);
95  J = mfem::Vector(int(num_elements) * qpts_per_elem * spatial_dim * geometry_dim);
96 
97 #define DISPATCH_KERNEL(GEOM, P, Q, BDR) \
98  if (geom == mfem::Geometry::GEOM && p == P && q == Q && (spatial_dim - geometry_dim) == BDR) { \
99  compute_geometric_factors<Q, mfem::Geometry::GEOM, H1<P, dimension_of(mfem::Geometry::GEOM) + BDR> >(X, J, X_e, \
100  element_ids); \
101  return; \
102  }
103 
104  DISPATCH_KERNEL(SEGMENT, 1, 1, 1);
105  DISPATCH_KERNEL(SEGMENT, 1, 2, 1);
106  DISPATCH_KERNEL(SEGMENT, 1, 3, 1);
107  DISPATCH_KERNEL(SEGMENT, 1, 4, 1);
108 
109  DISPATCH_KERNEL(SEGMENT, 2, 1, 1);
110  DISPATCH_KERNEL(SEGMENT, 2, 2, 1);
111  DISPATCH_KERNEL(SEGMENT, 2, 3, 1);
112  DISPATCH_KERNEL(SEGMENT, 2, 4, 1);
113 
114  DISPATCH_KERNEL(SEGMENT, 3, 1, 1);
115  DISPATCH_KERNEL(SEGMENT, 3, 2, 1);
116  DISPATCH_KERNEL(SEGMENT, 3, 3, 1);
117  DISPATCH_KERNEL(SEGMENT, 3, 4, 1);
118 
120 
121  DISPATCH_KERNEL(TRIANGLE, 1, 1, 0);
122  DISPATCH_KERNEL(TRIANGLE, 1, 2, 0);
123  DISPATCH_KERNEL(TRIANGLE, 1, 3, 0);
124  DISPATCH_KERNEL(TRIANGLE, 1, 4, 0);
125 
126  DISPATCH_KERNEL(TRIANGLE, 2, 1, 0);
127  DISPATCH_KERNEL(TRIANGLE, 2, 2, 0);
128  DISPATCH_KERNEL(TRIANGLE, 2, 3, 0);
129  DISPATCH_KERNEL(TRIANGLE, 2, 4, 0);
130 
131  DISPATCH_KERNEL(TRIANGLE, 3, 1, 0);
132  DISPATCH_KERNEL(TRIANGLE, 3, 2, 0);
133  DISPATCH_KERNEL(TRIANGLE, 3, 3, 0);
134  DISPATCH_KERNEL(TRIANGLE, 3, 4, 0);
135 
136  DISPATCH_KERNEL(TRIANGLE, 1, 1, 1);
137  DISPATCH_KERNEL(TRIANGLE, 1, 2, 1);
138  DISPATCH_KERNEL(TRIANGLE, 1, 3, 1);
139  DISPATCH_KERNEL(TRIANGLE, 1, 4, 1);
140 
141  DISPATCH_KERNEL(TRIANGLE, 2, 1, 1);
142  DISPATCH_KERNEL(TRIANGLE, 2, 2, 1);
143  DISPATCH_KERNEL(TRIANGLE, 2, 3, 1);
144  DISPATCH_KERNEL(TRIANGLE, 2, 4, 1);
145 
146  DISPATCH_KERNEL(TRIANGLE, 3, 1, 1);
147  DISPATCH_KERNEL(TRIANGLE, 3, 2, 1);
148  DISPATCH_KERNEL(TRIANGLE, 3, 3, 1);
149  DISPATCH_KERNEL(TRIANGLE, 3, 4, 1);
150 
152 
153  DISPATCH_KERNEL(SQUARE, 1, 1, 0);
154  DISPATCH_KERNEL(SQUARE, 1, 2, 0);
155  DISPATCH_KERNEL(SQUARE, 1, 3, 0);
156  DISPATCH_KERNEL(SQUARE, 1, 4, 0);
157 
158  DISPATCH_KERNEL(SQUARE, 2, 1, 0);
159  DISPATCH_KERNEL(SQUARE, 2, 2, 0);
160  DISPATCH_KERNEL(SQUARE, 2, 3, 0);
161  DISPATCH_KERNEL(SQUARE, 2, 4, 0);
162 
163  DISPATCH_KERNEL(SQUARE, 3, 1, 0);
164  DISPATCH_KERNEL(SQUARE, 3, 2, 0);
165  DISPATCH_KERNEL(SQUARE, 3, 3, 0);
166  DISPATCH_KERNEL(SQUARE, 3, 4, 0);
167 
168  DISPATCH_KERNEL(SQUARE, 1, 1, 1);
169  DISPATCH_KERNEL(SQUARE, 1, 2, 1);
170  DISPATCH_KERNEL(SQUARE, 1, 3, 1);
171  DISPATCH_KERNEL(SQUARE, 1, 4, 1);
172 
173  DISPATCH_KERNEL(SQUARE, 2, 1, 1);
174  DISPATCH_KERNEL(SQUARE, 2, 2, 1);
175  DISPATCH_KERNEL(SQUARE, 2, 3, 1);
176  DISPATCH_KERNEL(SQUARE, 2, 4, 1);
177 
178  DISPATCH_KERNEL(SQUARE, 3, 1, 1);
179  DISPATCH_KERNEL(SQUARE, 3, 2, 1);
180  DISPATCH_KERNEL(SQUARE, 3, 3, 1);
181  DISPATCH_KERNEL(SQUARE, 3, 4, 1);
182 
184 
185  DISPATCH_KERNEL(TETRAHEDRON, 1, 1, 0);
186  DISPATCH_KERNEL(TETRAHEDRON, 1, 2, 0);
187  DISPATCH_KERNEL(TETRAHEDRON, 1, 3, 0);
188  DISPATCH_KERNEL(TETRAHEDRON, 1, 4, 0);
189 
190  DISPATCH_KERNEL(TETRAHEDRON, 2, 1, 0);
191  DISPATCH_KERNEL(TETRAHEDRON, 2, 2, 0);
192  DISPATCH_KERNEL(TETRAHEDRON, 2, 3, 0);
193  DISPATCH_KERNEL(TETRAHEDRON, 2, 4, 0);
194 
195  DISPATCH_KERNEL(TETRAHEDRON, 3, 1, 0);
196  DISPATCH_KERNEL(TETRAHEDRON, 3, 2, 0);
197  DISPATCH_KERNEL(TETRAHEDRON, 3, 3, 0);
198  DISPATCH_KERNEL(TETRAHEDRON, 3, 4, 0);
199 
201 
202  DISPATCH_KERNEL(CUBE, 1, 1, 0);
203  DISPATCH_KERNEL(CUBE, 1, 2, 0);
204  DISPATCH_KERNEL(CUBE, 1, 3, 0);
205  DISPATCH_KERNEL(CUBE, 1, 4, 0);
206 
207  DISPATCH_KERNEL(CUBE, 2, 1, 0);
208  DISPATCH_KERNEL(CUBE, 2, 2, 0);
209  DISPATCH_KERNEL(CUBE, 2, 3, 0);
210  DISPATCH_KERNEL(CUBE, 2, 4, 0);
211 
212  DISPATCH_KERNEL(CUBE, 3, 1, 0);
213  DISPATCH_KERNEL(CUBE, 3, 2, 0);
214  DISPATCH_KERNEL(CUBE, 3, 3, 0);
215  DISPATCH_KERNEL(CUBE, 3, 4, 0);
216 
217 #undef DISPATCH_KERNEL
218 
219  std::cout << "should never be reached " << std::endl;
220 }
221 
222 } // namespace smith
This file contains helper traits and enumerations for classifying finite elements.
Accelerator functionality.
Definition: smith.cpp:36
void compute_geometric_factors(mfem::Vector &positions_q, mfem::Vector &jacobians_q, const mfem::Vector &positions_e, const std::vector< int > &elements)
a kernel to compute the positions and jacobians at each quadrature point (mfem calls this "geometric ...
constexpr int num_quadrature_points(mfem::Geometry::Type g, int Q)
return the number of quadrature points in a Gauss-Legendre rule with parameter "Q"
Definition: geometry.hpp:31
constexpr int dimension_of(mfem::Geometry::Type g)
Returns the dimension of an element geometry.
Definition: geometry.hpp:55
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
const std::vector< int > & get_mfem_ids(mfem::Geometry::Type geom) const
get elements by geometry type
Definition: domain.hpp:195
const mesh_t & mesh_
the underyling mesh for this domain
Definition: domain.hpp:45
std::size_t num_elements
the number of elements in the domain
mfem::Vector J
Jacobians of the element transformations at all quadrature points.
mfem::Vector X
Mapped (physical) coordinates of all quadrature points.
GeometricFactors()
default ctor, leaving this object uninitialized
a convenience class for generating information about tensor product integration rules from the underl...
Template prototype for finite element implementations.
Arbitrary-rank tensor class.
Definition: tensor.hpp:28