Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
integral.hpp
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 #pragma once
8 
9 #include <array>
10 #include <memory>
11 
12 #include "mfem.hpp"
13 
15 #include "smith/numerics/functional/geometric_factors.hpp"
16 #include "smith/numerics/functional/function_signature.hpp"
17 #include "smith/numerics/functional/domain_integral_kernels.hpp"
18 #include "smith/numerics/functional/boundary_integral_kernels.hpp"
19 #include "smith/numerics/functional/interior_face_integral_kernels.hpp"
20 #include "smith/numerics/functional/differentiate_wrt.hpp"
21 
22 namespace smith {
23 
25 struct Integral {
27  static constexpr std::size_t num_types = 3;
28 
38  Integral(const Domain& d, std::vector<uint32_t> trial_space_indices)
39  : domain_(d), active_trial_spaces_(trial_space_indices)
40  {
41  std::size_t num_trial_spaces = trial_space_indices.size();
42  evaluation_with_AD_.resize(num_trial_spaces);
43  jvp_.resize(num_trial_spaces);
44  element_gradient_.resize(num_trial_spaces);
45 
46  for (uint32_t i = 0; i < num_trial_spaces; i++) {
48  }
49  }
50 
66  void Mult(double t, const std::vector<mfem::BlockVector>& input_E, mfem::BlockVector& output_E,
67  uint32_t differentiation_index, bool update_state) const
68  {
69  output_E = 0.0;
70 
71  bool with_AD =
72  (functional_to_integral_index_.count(differentiation_index) > 0 && differentiation_index != NO_DIFFERENTIATION);
73  auto& kernels =
74  (with_AD) ? evaluation_with_AD_[functional_to_integral_index_.at(differentiation_index)] : evaluation_;
75  for (auto& [geometry, func] : kernels) {
76  std::vector<const double*> inputs(active_trial_spaces_.size());
77  for (std::size_t i = 0; i < active_trial_spaces_.size(); i++) {
78  inputs[i] = input_E[uint32_t(active_trial_spaces_[i])].GetBlock(geometry).Read();
79  }
80  func(t, inputs, output_E.GetBlock(geometry).ReadWrite(), update_state);
81  }
82  }
83 
94  void GradientMult(const mfem::BlockVector& dinput_E, mfem::BlockVector& doutput_E,
95  uint32_t differentiation_index) const
96  {
97  doutput_E = 0.0;
98 
99  // if this integral actually depends on the specified variable
100  if (functional_to_integral_index_.count(differentiation_index) > 0) {
101  for (auto& [geometry, func] : jvp_[functional_to_integral_index_.at(differentiation_index)]) {
102  func(dinput_E.GetBlock(geometry).Read(), doutput_E.GetBlock(geometry).ReadWrite());
103  }
104  }
105  }
106 
114  void ComputeElementGradients(std::map<mfem::Geometry::Type, ExecArray<double, 3, ExecutionSpace::CPU> >& K_e,
115  uint32_t differentiation_index) const
116  {
117  // if this integral actually depends on the specified variable
118  if (functional_to_integral_index_.count(differentiation_index) > 0) {
119  for (auto& [geometry, func] : element_gradient_[functional_to_integral_index_.at(differentiation_index)]) {
120  func(view(K_e[geometry]));
121  }
122  }
123  }
124 
131  bool DependsOn(uint32_t which) const
132  {
133  for (uint32_t i : active_trial_spaces_) {
134  if (which == i) return true;
135  }
136  return false;
137  }
138 
141 
143  using eval_func = std::function<void(double, const std::vector<const double*>&, double*, bool)>;
144 
146  std::map<mfem::Geometry::Type, eval_func> evaluation_;
147 
149  std::vector<std::map<mfem::Geometry::Type, eval_func> > evaluation_with_AD_;
150 
152  using jacobian_vector_product_func = std::function<void(const double*, double*)>;
153 
155  std::vector<std::map<mfem::Geometry::Type, jacobian_vector_product_func> > jvp_;
156 
159 
161  std::vector<std::map<mfem::Geometry::Type, grad_func> > element_gradient_;
162 
164  std::vector<uint32_t> active_trial_spaces_;
165 
180  std::map<uint32_t, uint32_t> functional_to_integral_index_;
181 
183  std::map<mfem::Geometry::Type, GeometricFactors> geometric_factors_;
184 };
185 
200 template <mfem::Geometry::Type geom, int Q, typename test, typename... trials, typename lambda_type,
201  typename qpt_data_type>
202 void generate_kernels(FunctionSignature<test(trials...)> s, Integral& integral, const lambda_type& qf,
203  std::shared_ptr<QuadratureData<qpt_data_type> > qdata)
204 {
205  integral.geometric_factors_[geom] = GeometricFactors(integral.domain_, Q, geom);
206  GeometricFactors& gf = integral.geometric_factors_[geom];
207  if (gf.num_elements == 0) return;
208 
209  const double* positions = gf.X.Read();
210  const double* jacobians = gf.J.Read();
211  const uint32_t num_elements = uint32_t(gf.num_elements);
212  const uint32_t qpts_per_element = num_quadrature_points(geom, Q);
213 
214  std::shared_ptr<zero> dummy_derivatives;
215  integral.evaluation_[geom] = domain_integral::evaluation_kernel<NO_DIFFERENTIATION, Q, geom>(
216  s, qf, positions, jacobians, qdata, dummy_derivatives, num_elements);
217 
218  constexpr std::size_t num_args = s.num_args;
219  [[maybe_unused]] static constexpr int dim = dimension_of(geom);
220  for_constexpr<num_args>([&](auto index) {
221  // allocate memory for the derivatives of the q-function at each quadrature point
222  //
223  // Note: ptrs' lifetime is managed in an unusual way! It is captured by-value in the
224  // action_of_gradient functor below to augment the reference count, and extend its lifetime to match
225  // that of the DomainIntegral that allocated it.
226  using derivative_type = decltype(domain_integral::get_derivative_type<index, dim, trials...>(qf, qpt_data_type{}));
227  auto ptr = accelerator::make_shared_array<ExecutionSpace::CPU, derivative_type>(num_elements * qpts_per_element);
228 
229  integral.evaluation_with_AD_[index][geom] =
230  domain_integral::evaluation_kernel<index, Q, geom>(s, qf, positions, jacobians, qdata, ptr, num_elements);
231 
232  integral.jvp_[index][geom] = domain_integral::jacobian_vector_product_kernel<index, Q, geom>(s, ptr, num_elements);
233  integral.element_gradient_[index][geom] =
234  domain_integral::element_gradient_kernel<index, Q, geom>(s, ptr, num_elements);
235  });
236 }
237 
252 template <typename s, int Q, int dim, typename lambda_type, typename qpt_data_type>
253 Integral MakeDomainIntegral(const Domain& domain, const lambda_type& qf,
254  std::shared_ptr<QuadratureData<qpt_data_type> > qdata,
255  std::vector<uint32_t> argument_indices)
256 {
257  FunctionSignature<s> signature;
258 
259  SLIC_ERROR_IF(domain.type_ != Domain::Type::Elements, "Error: trying to evaluate a domain integral over a boundary");
260 
261  Integral integral(domain, argument_indices);
262 
263  if constexpr (dim == 2) {
264  generate_kernels<mfem::Geometry::TRIANGLE, Q>(signature, integral, qf, qdata);
265  generate_kernels<mfem::Geometry::SQUARE, Q>(signature, integral, qf, qdata);
266  }
267 
268  if constexpr (dim == 3) {
269  generate_kernels<mfem::Geometry::TETRAHEDRON, Q>(signature, integral, qf, qdata);
270  generate_kernels<mfem::Geometry::CUBE, Q>(signature, integral, qf, qdata);
271  }
272 
273  return integral;
274 }
275 
289 template <mfem::Geometry::Type geom, int Q, typename test, typename... trials, typename lambda_type>
290 void generate_bdr_kernels(FunctionSignature<test(trials...)> s, Integral& integral, const lambda_type& qf)
291 {
292  integral.geometric_factors_[geom] = GeometricFactors(integral.domain_, Q, geom);
293  GeometricFactors& gf = integral.geometric_factors_[geom];
294  if (gf.num_elements == 0) return;
295 
296  const double* positions = gf.X.Read();
297  const double* jacobians = gf.J.Read();
298  const uint32_t num_elements = uint32_t(gf.num_elements);
299  const uint32_t qpts_per_element = num_quadrature_points(geom, Q);
300 
301  std::shared_ptr<zero> dummy_derivatives;
302  integral.evaluation_[geom] = boundary_integral::evaluation_kernel<NO_DIFFERENTIATION, Q, geom>(
303  s, qf, positions, jacobians, dummy_derivatives, num_elements);
304 
305  constexpr std::size_t num_args = s.num_args;
306  [[maybe_unused]] static constexpr int dim = dimension_of(geom);
307  for_constexpr<num_args>([&](auto index) {
308  // allocate memory for the derivatives of the q-function at each quadrature point
309  //
310  // Note: ptrs' lifetime is managed in an unusual way! It is captured by-value in the
311  // action_of_gradient functor below to augment the reference count, and extend its lifetime to match
312  // that of the boundaryIntegral that allocated it.
313  using derivative_type = decltype(boundary_integral::get_derivative_type<index, dim, trials...>(qf));
314  auto ptr = accelerator::make_shared_array<ExecutionSpace::CPU, derivative_type>(num_elements * qpts_per_element);
315 
316  integral.evaluation_with_AD_[index][geom] =
317  boundary_integral::evaluation_kernel<index, Q, geom>(s, qf, positions, jacobians, ptr, num_elements);
318 
319  integral.jvp_[index][geom] =
320  boundary_integral::jacobian_vector_product_kernel<index, Q, geom>(s, ptr, num_elements);
321  integral.element_gradient_[index][geom] =
322  boundary_integral::element_gradient_kernel<index, Q, geom>(s, ptr, num_elements);
323  });
324 }
325 
340 template <typename s, int Q, int dim, typename lambda_type>
341 Integral MakeBoundaryIntegral(const Domain& domain, const lambda_type& qf, std::vector<uint32_t> argument_indices)
342 {
343  FunctionSignature<s> signature;
344 
345  SLIC_ERROR_IF(domain.type_ != Domain::Type::BoundaryElements,
346  "Error: trying to evaluate a boundary integral over a non-boundary domain of integration");
347 
348  Integral integral(domain, argument_indices);
349 
350  if constexpr (dim == 1) {
351  generate_bdr_kernels<mfem::Geometry::SEGMENT, Q>(signature, integral, qf);
352  }
353 
354  if constexpr (dim == 2) {
355  generate_bdr_kernels<mfem::Geometry::TRIANGLE, Q>(signature, integral, qf);
356  generate_bdr_kernels<mfem::Geometry::SQUARE, Q>(signature, integral, qf);
357  }
358 
359  return integral;
360 }
361 
375 template <mfem::Geometry::Type geom, int Q, typename test, typename... trials, typename lambda_type>
376 void generate_interior_face_kernels(FunctionSignature<test(trials...)> s, Integral& integral, const lambda_type& qf)
377 {
378  integral.geometric_factors_[geom] = GeometricFactors(integral.domain_, Q, geom);
379  GeometricFactors& gf = integral.geometric_factors_[geom];
380  if (gf.num_elements == 0) return;
381 
382  const double* positions = gf.X.Read();
383  const double* jacobians = gf.J.Read();
384  const uint32_t num_elements = uint32_t(gf.num_elements);
385  const uint32_t qpts_per_element = num_quadrature_points(geom, Q);
386 
387  std::shared_ptr<zero> dummy_derivatives;
388  integral.evaluation_[geom] = interior_face_integral::evaluation_kernel<NO_DIFFERENTIATION, Q, geom>(
389  s, qf, positions, jacobians, dummy_derivatives, num_elements);
390 
391  constexpr std::size_t num_args = s.num_args;
392  [[maybe_unused]] static constexpr int dim = dimension_of(geom);
393  for_constexpr<num_args>([&](auto index) {
394  // allocate memory for the derivatives of the q-function at each quadrature point
395  //
396  // Note: ptrs' lifetime is managed in an unusual way! It is captured by-value in the
397  // action_of_gradient functor below to augment the reference count, and extend its lifetime to match
398  // that of the boundaryIntegral that allocated it.
399  using derivative_type = decltype(interior_face_integral::get_derivative_type<index, dim, trials...>(qf));
400  auto ptr = accelerator::make_shared_array<ExecutionSpace::CPU, derivative_type>(num_elements * qpts_per_element);
401 
402  integral.evaluation_with_AD_[index][geom] =
403  interior_face_integral::evaluation_kernel<index, Q, geom>(s, qf, positions, jacobians, ptr, num_elements);
404 
405  integral.jvp_[index][geom] =
406  interior_face_integral::jacobian_vector_product_kernel<index, Q, geom>(s, ptr, num_elements);
407  integral.element_gradient_[index][geom] =
408  interior_face_integral::element_gradient_kernel<index, Q, geom>(s, ptr, num_elements);
409  });
410 }
411 
426 template <typename s, int Q, int dim, typename lambda_type>
427 Integral MakeInteriorFaceIntegral(const Domain& domain, const lambda_type& qf, std::vector<uint32_t> argument_indices)
428 {
429  FunctionSignature<s> signature;
430 
431  SLIC_ERROR_IF(domain.type_ != Domain::Type::InteriorBoundaryElements,
432  "Error: trying to evaluate a boundary integral over a non-boundary domain of integration");
433 
434  Integral integral(domain, argument_indices);
435 
436  if constexpr (dim == 1) {
437  generate_interior_face_kernels<mfem::Geometry::SEGMENT, Q>(signature, integral, qf);
438  }
439 
440  if constexpr (dim == 2) {
441  generate_interior_face_kernels<mfem::Geometry::TRIANGLE, Q>(signature, integral, qf);
442  generate_interior_face_kernels<mfem::Geometry::SQUARE, Q>(signature, integral, qf);
443  }
444 
445  return integral;
446 }
447 
448 } // namespace smith
This file contains the interface used for initializing/terminating any hardware accelerator-related f...
Accelerator functionality.
Definition: smith.cpp:36
auto view(axom::Array< T, dim, space > &arr)
convenience function for creating a view of an axom::Array type
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
void generate_kernels(FunctionSignature< test(trials...)> s, Integral &integral, const lambda_type &qf, std::shared_ptr< QuadratureData< qpt_data_type > > qdata)
function to generate kernels held by an Integral object of type "Domain", with a specific element typ...
Definition: integral.hpp:202
constexpr int dimension_of(mfem::Geometry::Type g)
Returns the dimension of an element geometry.
Definition: geometry.hpp:55
axom::Array< T, dim, detail::execution_to_memory_v< space > > ExecArray
Alias for an Array corresponding to a particular ExecutionSpace.
void generate_bdr_kernels(FunctionSignature< test(trials...)> s, Integral &integral, const lambda_type &qf)
function to generate kernels held by an Integral object of type "BoundaryDomain", with a specific ele...
Definition: integral.hpp:290
Integral MakeInteriorFaceIntegral(const Domain &domain, const lambda_type &qf, std::vector< uint32_t > argument_indices)
function to generate kernels held by an Integral object of type "Boundary", for all element types
Definition: integral.hpp:427
Integral MakeBoundaryIntegral(const Domain &domain, const lambda_type &qf, std::vector< uint32_t > argument_indices)
function to generate kernels held by an Integral object of type "Boundary", for all element types
Definition: integral.hpp:341
void generate_interior_face_kernels(FunctionSignature< test(trials...)> s, Integral &integral, const lambda_type &qf)
function to generate kernels held by an Integral object of type "InteriorFaceDomain",...
Definition: integral.hpp:376
Integral MakeDomainIntegral(const Domain &domain, const lambda_type &qf, std::shared_ptr< QuadratureData< qpt_data_type > > qdata, std::vector< uint32_t > argument_indices)
function to generate kernels held by an Integral object of type "Domain", for all element types
Definition: integral.hpp:253
axom::ArrayView< T, dim, detail::execution_to_memory_v< space > > ExecArrayView
Alias for an ArrayView corresponding to a particular ExecutionSpace.
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
Type type_
whether the elements in this domain are on the boundary or not
Definition: domain.hpp:51
a class that computes and stores positions and jacobians at each quadrature point
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.
a class for representing a Integral calculations and their derivatives
Definition: integral.hpp:25
std::map< mfem::Geometry::Type, GeometricFactors > geometric_factors_
the spatial positions and jacobians (dx_dxi) for each element type and quadrature point
Definition: integral.hpp:183
std::function< void(double, const std::vector< const double * > &, double *, bool)> eval_func
signature of integral evaluation kernel
Definition: integral.hpp:143
void GradientMult(const mfem::BlockVector &dinput_E, mfem::BlockVector &doutput_E, uint32_t differentiation_index) const
evaluate the jacobian(with respect to some trial space)-vector product of this integral
Definition: integral.hpp:94
Integral(const Domain &d, std::vector< uint32_t > trial_space_indices)
Construct an "empty" Integral object, whose kernels are to be initialized later in one of the Make***...
Definition: integral.hpp:38
std::map< mfem::Geometry::Type, eval_func > evaluation_
kernels for integral evaluation over each type of element
Definition: integral.hpp:146
std::function< void(ExecArrayView< double, 3, ExecutionSpace::CPU >)> grad_func
signature of element gradient kernel
Definition: integral.hpp:158
std::vector< std::map< mfem::Geometry::Type, eval_func > > evaluation_with_AD_
kernels for integral evaluation + derivative w.r.t. specified argument over each type of element
Definition: integral.hpp:149
std::map< uint32_t, uint32_t > functional_to_integral_index_
a way of translating between the indices used by Functional and Integral to refer to the same trial s...
Definition: integral.hpp:180
std::vector< std::map< mfem::Geometry::Type, jacobian_vector_product_func > > jvp_
kernels for jacobian-vector product of integral calculation
Definition: integral.hpp:155
std::function< void(const double *, double *)> jacobian_vector_product_func
signature of element jvp kernel
Definition: integral.hpp:152
bool DependsOn(uint32_t which) const
returns whether or not this integral depends on argument which
Definition: integral.hpp:131
std::vector< std::map< mfem::Geometry::Type, grad_func > > element_gradient_
kernels for calculation of element jacobians
Definition: integral.hpp:161
void ComputeElementGradients(std::map< mfem::Geometry::Type, ExecArray< double, 3, ExecutionSpace::CPU > > &K_e, uint32_t differentiation_index) const
evaluate the jacobian (with respect to some trial space) of this integral
Definition: integral.hpp:114
std::vector< uint32_t > active_trial_spaces_
a list of the trial spaces that take part in this integrand
Definition: integral.hpp:164
static constexpr std::size_t num_types
the number of different kinds of integration domains
Definition: integral.hpp:27
Domain domain_
information about which elements to integrate over
Definition: integral.hpp:140
void Mult(double t, const std::vector< mfem::BlockVector > &input_E, mfem::BlockVector &output_E, uint32_t differentiation_index, bool update_state) const
evaluate the integral, optionally storing q-function derivatives with respect to a specific trial spa...
Definition: integral.hpp:66
A class for storing and access user-defined types at quadrature points.