Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
shape_aware_functional.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 
14 #pragma once
15 
18 
19 namespace smith {
20 
21 namespace detail {
22 
30 template <int dim, typename shape_type>
32  public:
39  : J_(Identity<dim>() + get<DERIVATIVE>(p)), detJ_(det(J_)), inv_J_(inv(J_)), inv_JT_(transpose(inv_J_))
40  {
41  }
42 
54  template <typename trial_type, typename space_type>
55  SMITH_HOST_DEVICE auto modify_trial_argument(space_type /* space */, const trial_type& u) const
56  {
57  if constexpr (space_type{}.family == Family::H1 || space_type{}.family == Family::L2) {
58  // Our updated spatial coordinate is x' = x + p. We want to calculate
59  // du/dx' = du/dx * dx/dx'
60  // = du/dx * (dx'/dx)^-1
61  // = du_dx * (J)^-1
62  auto trial_derivative = dot(get<DERIVATIVE>(u), inv_J_);
63 
64  return smith::tuple{get<VALUE>(u), trial_derivative};
65  }
66 
67  if constexpr (space_type{}.family == Family::HCURL) {
68  auto modified_val = dot(get<VALUE>(u), inv_J_);
69  if constexpr (dim == 2) {
70  auto modified_derivative = get<DERIVATIVE>(u) / detJ_;
71  return smith::tuple{modified_val, modified_derivative};
72  }
73  if constexpr (dim == 3) {
74  auto modified_derivative = dot(get<DERIVATIVE>(u), transpose(J_));
75  return smith::tuple{modified_val, modified_derivative};
76  }
77  }
78  }
79 
91  template <typename test_space, typename q_func_type>
92  SMITH_HOST_DEVICE auto modify_shape_aware_qf_return(test_space /*test*/, const q_func_type& v)
93  {
94  if constexpr (std::is_same_v<test_space, double>) {
95  return detJ_ * v;
96  } else {
97  if constexpr (test_space{}.family == Family::H1 || test_space{}.family == Family::L2) {
98  auto modified_flux = dot(get<FLUX>(v), inv_JT_) * detJ_;
99  auto modified_source = get<SOURCE>(v) * detJ_;
100 
101  return smith::tuple{modified_source, modified_flux};
102  }
103  if constexpr (test_space{}.family == Family::HCURL) {
104  auto modified_source = dot(get<SOURCE>(v), inv_JT_) * detJ_;
105  if constexpr (dim == 3) {
106  auto modified_flux = dot(get<FLUX>(v), J_);
107  return smith::tuple{modified_source, modified_flux};
108  } else {
109  return smith::tuple{modified_source, get<FLUX>(v)};
110  }
111  }
112  }
113  }
114 
115  private:
117  using jacobian_type = std::remove_reference_t<decltype(get<DERIVATIVE>(std::declval<shape_type>()))>;
118  using detJ_type = decltype(det(std::declval<jacobian_type>()));
119  using inv_J_type = decltype(inv(std::declval<jacobian_type>()));
120  using inv_JT_type = decltype(inv(transpose(std::declval<jacobian_type>())));
122 
124  jacobian_type J_;
125 
127  detJ_type detJ_;
128 
130  inv_J_type inv_J_;
131 
133  inv_JT_type inv_JT_;
134 };
135 
148 template <typename position_type, typename shape_type>
149 SMITH_HOST_DEVICE auto compute_boundary_area_correction(const position_type& X, const shape_type& shape)
150 {
151  auto x_prime = X + shape;
152 
153  auto n = cross(get<DERIVATIVE>(x_prime));
154 
155  // smith::Functional's boundary integrals multiply the q-function output by
156  // norm(cross(dX_dxi)) at that quadrature point, but if we impose a shape displacement
157  // then that weight needs to be corrected. The new weight should be
158  // norm(cross(dX_dxi + dp_dxi)), so we multiply by the ratio w_new / w_old
159  // to get
160  // q * area_correction * w_old
161  // = q * (w_new / w_old) * w_old
162  // = q * w_new
163 
164  auto area_correction = norm(n) / norm(cross(get<DERIVATIVE>(X)));
165 
166  return area_correction;
167 }
168 
193 template <typename lambda, typename coord_type, typename shape_type, typename space_types, typename trial_types,
194  typename correction_type, int... i>
195 SMITH_HOST_DEVICE auto apply_shape_aware_qf_helper(const lambda& qf, double t, const coord_type& position,
196  const shape_type& shape, const space_types& space_tuple,
197  const trial_types& arg_tuple, const correction_type& correction,
198  std::integer_sequence<int, i...>)
199 {
201  "Argument and finite element space tuples are not the same size.");
202 
203  auto x = smith::tuple{get<VALUE>(position) + get<VALUE>(shape),
204 
205  // x := X + u,
206  // so, dx/dxi = dX/dxi + du/dxi
207  // = dX/dxi + du/dX * dX/dxi
208  get<DERIVATIVE>(position) + get<DERIVATIVE>(shape) * get<DERIVATIVE>(position)};
209 
210  return qf(t, x, correction.modify_trial_argument(smith::get<i>(space_tuple), smith::get<i>(arg_tuple))...);
211 }
212 
239 template <typename lambda, typename coord_type, typename state_type, typename shape_type, typename space_types,
240  typename trial_types, typename correction_type, int... i>
241 SMITH_HOST_DEVICE auto apply_shape_aware_qf_helper_with_state(const lambda& qf, double t, const coord_type& position,
242  state_type& state, const shape_type& shape,
243  const space_types& space_tuple,
244  const trial_types& arg_tuple,
245  const correction_type& correction,
246  std::integer_sequence<int, i...>)
247 {
249  "Argument and finite element space tuples are not the same size.");
250 
251  auto x = smith::tuple{get<VALUE>(position) + get<VALUE>(shape),
252 
253  // x := X + u,
254  // so, dx/dxi = dX/dxi + du/dxi
255  // = dX/dxi + du/dX * dX/dxi
256  get<DERIVATIVE>(position) + get<DERIVATIVE>(shape) * get<DERIVATIVE>(position)};
257 
258  return qf(t, x, state, correction.modify_trial_argument(smith::get<i>(space_tuple), smith::get<i>(arg_tuple))...);
259 }
260 
261 } // namespace detail
262 
264 template <typename T1, typename T2, ExecutionSpace exec = smith::ExecutionSpace::CPU>
265 class ShapeAwareFunctional;
267 
285 template <typename test, typename shape, typename... trials, ExecutionSpace exec>
286 class ShapeAwareFunctional<shape, test(trials...), exec> {
288  static constexpr test test_space{};
289 
291  static constexpr shape shape_space{};
292 
294  static constexpr uint32_t num_trial_spaces = sizeof...(trials);
295 
296  public:
307  template <typename test_type = test, typename = std::enable_if_t<!std::is_same_v<double, test_type>>>
308  ShapeAwareFunctional(const mfem::ParFiniteElementSpace* shape_fes, const mfem::ParFiniteElementSpace* test_fes,
309  std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_fes)
310  {
311  static_assert(shape_space.family == Family::H1, "Only H1 spaces allowed for shape displacements");
312 
313  std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces + 1> prepended_spaces;
314 
315  prepended_spaces[0] = shape_fes;
316 
317  for (uint32_t i = 0; i < num_trial_spaces; ++i) {
318  prepended_spaces[1 + i] = trial_fes[i];
319  }
320 
321  functional_ = std::make_unique<Functional<test(shape, trials...), exec>>(test_fes, prepended_spaces);
322  }
323 
333  template <typename test_type = test, typename = std::enable_if_t<std::is_same_v<double, test_type>>>
334  ShapeAwareFunctional(const mfem::ParFiniteElementSpace* shape_fes,
335  std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces> trial_fes)
336  {
337  static_assert(shape_space.family == Family::H1, "Only H1 spaces allowed for shape displacements");
338 
339  std::array<const mfem::ParFiniteElementSpace*, num_trial_spaces + 1> prepended_spaces;
340 
341  prepended_spaces[0] = shape_fes;
342 
343  for (uint32_t i = 0; i < num_trial_spaces; ++i) {
344  prepended_spaces[1 + i] = trial_fes[i];
345  }
346 
347  functional_ = std::make_unique<Functional<double(shape, trials...), exec>>(prepended_spaces);
348  }
349 
354  template <typename Integrand, int dim, int... args>
355  struct ShapeAwareIntegrandWrapper {
357  ShapeAwareIntegrandWrapper(Integrand integrand) : integrand_(integrand) {}
359  Integrand integrand_;
372  template <typename PositionType, typename ShapeValueType, typename... QFuncArgs>
373  SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val,
374  QFuncArgs... qfunc_args) const
375  {
376  auto qfunc_tuple = make_tuple(qfunc_args...);
377  [[maybe_unused]] smith::tuple<trials...> local_trial_spaces{};
378  auto reduced_trial_space_tuple = make_tuple(get<args>(local_trial_spaces)...);
379 
380  detail::ShapeCorrection shape_correction(Dimension<dim>{}, shape_val);
381  // TODO(CUDA): When this is compiled to device code, the below make_integer_sequence will
382  // need to change to a camp integer sequence.
383  auto unmodified_qf_return = detail::apply_shape_aware_qf_helper(
384  integrand_, time, x, shape_val, reduced_trial_space_tuple, qfunc_tuple, shape_correction,
385  std::make_integer_sequence<int, sizeof...(qfunc_args)>{});
386  return shape_correction.modify_shape_aware_qf_return(test_space, unmodified_qf_return);
387  }
388  };
389 
394  template <typename Integrand, int dim, int... args>
395  struct ShapeAwareIntegrandWrapperWithState {
397  ShapeAwareIntegrandWrapperWithState(Integrand integrand) : integrand_(integrand) {}
399  Integrand integrand_;
414  template <typename PositionType, typename StateType, typename ShapeValueType, typename... QFuncArgs>
415  SMITH_HOST_DEVICE auto operator()(double time, PositionType x, StateType& state, ShapeValueType shape_val,
416  QFuncArgs... qfunc_args) const
417  {
418  auto qfunc_tuple = make_tuple(qfunc_args...);
419  [[maybe_unused]] smith::tuple<trials...> local_trial_spaces{};
420  auto reduced_trial_space_tuple = make_tuple(get<args>(local_trial_spaces)...);
421 
422  detail::ShapeCorrection shape_correction(Dimension<dim>{}, shape_val);
423  // TODO(CUDA): When this is compiled to device code, the below make_integer_sequence will
424  // need to change to a camp integer sequence.
425  auto unmodified_qf_return = detail::apply_shape_aware_qf_helper_with_state(
426  integrand_, time, x, state, shape_val, reduced_trial_space_tuple, qfunc_tuple, shape_correction,
427  std::make_integer_sequence<int, sizeof...(qfunc_args)>{});
428  return shape_correction.modify_shape_aware_qf_return(test_space, unmodified_qf_return);
429  }
430  };
431 
447  template <int dim, int... args, typename lambda, typename qpt_data_type = Nothing>
448  void AddDomainIntegral(Dimension<dim>, DependsOn<args...>, const lambda& integrand, Domain& domain,
449  std::shared_ptr<QuadratureData<qpt_data_type>> qdata = NoQData)
450  {
451  if constexpr (std::is_same_v<qpt_data_type, Nothing>) {
452  functional_->AddDomainIntegral(Dimension<dim>{}, DependsOn<0, (args + 1)...>{},
453  ShapeAwareIntegrandWrapper<lambda, dim, args...>(integrand), domain, qdata);
454  } else {
455  functional_->AddDomainIntegral(Dimension<dim>{}, DependsOn<0, (args + 1)...>{},
456  ShapeAwareIntegrandWrapperWithState<lambda, dim, args...>(integrand), domain,
457  qdata);
458  }
459  }
460 
465  template <typename test_type, typename Integrand, int dim, int... args>
466  struct ShapeAwareInteriorBoundaryIntegrandWrapper {
468  ShapeAwareInteriorBoundaryIntegrandWrapper(Integrand integrand) : integrand_(integrand) {}
470  Integrand integrand_;
483  template <typename PositionType, typename ShapeValueType, typename... QFuncArgs>
484  SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val,
485  QFuncArgs... qfunc_args) const
486  {
487  auto unmodified_qf_returns = integrand_(time, x + shape_val, qfunc_args...);
488  auto correction = detail::compute_boundary_area_correction(x, shape_val);
489  return smith::tuple{correction * smith::get<0>(unmodified_qf_returns),
490  correction * smith::get<1>(unmodified_qf_returns)};
491  }
492  };
493 
498  template <typename Integrand, int dim, int... args>
499  struct ShapeAwareInteriorBoundaryIntegrandWrapper<double, Integrand, dim, args...> {
501  ShapeAwareInteriorBoundaryIntegrandWrapper(Integrand integrand) : integrand_(integrand) {}
503  Integrand integrand_;
516  template <typename PositionType, typename ShapeValueType, typename... QFuncArgs>
517  SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val,
518  QFuncArgs... qfunc_args) const
519  {
520  auto unmodified_qf_returns = integrand_(time, x + shape_val, qfunc_args...);
521  auto correction = detail::compute_boundary_area_correction(x, shape_val);
522  return correction * unmodified_qf_returns;
523  }
524  };
525 
540  template <int dim, int... args, typename lambda>
541  void AddInteriorFaceIntegral(Dimension<dim>, DependsOn<args...>, const lambda& integrand, Domain& domain)
542  {
543  functional_->AddInteriorFaceIntegral(
544  Dimension<dim>{}, DependsOn<0, (args + 1)...>{},
545  ShapeAwareInteriorBoundaryIntegrandWrapper<test, lambda, dim, args...>(integrand), domain);
546  }
547 
552  template <typename Integrand, int dim, int... args>
553  struct ShapeAwareBoundaryIntegrandWrapper {
555  ShapeAwareBoundaryIntegrandWrapper(Integrand integrand) : integrand_(integrand) {}
557  Integrand integrand_;
570  template <typename PositionType, typename ShapeValueType, typename... QFuncArgs>
571  SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val,
572  QFuncArgs... qfunc_args) const
573  {
574  auto unmodified_qf_return = integrand_(time, x + shape_val, qfunc_args...);
575  return unmodified_qf_return * detail::compute_boundary_area_correction(x, shape_val);
576  }
577  };
578 
593  template <int dim, int... args, typename lambda>
594  void AddBoundaryIntegral(Dimension<dim>, DependsOn<args...>, const lambda& integrand, Domain& domain)
595  {
596  functional_->AddBoundaryIntegral(Dimension<dim>{}, DependsOn<0, (args + 1)...>{},
597  ShapeAwareBoundaryIntegrandWrapper<lambda, dim, args...>(integrand), domain);
598  }
599 
615  template <uint32_t wrt, typename... T>
616  auto operator()(DifferentiateWRT<wrt>, double t, const T&... args)
617  {
618  return (*functional_)(DifferentiateWRT<wrt>{}, t, args...);
619  }
620 
634  template <typename... T>
635  auto operator()(double t, const T&... args)
636  {
637  return (*functional_)(t, args...);
638  }
639 
648  void updateQdata(bool update_flag) { functional_->updateQdata(update_flag); }
649 
650  private:
652  std::unique_ptr<Functional<test(shape, trials...), exec>> functional_;
653 };
654 
655 } // namespace smith
This file contains the interface used for initializing/terminating any hardware accelerator-related f...
#define SMITH_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc or amdclang and does nothing on ...
Definition: accelerator.hpp:37
void AddDomainIntegral(Dimension< dim >, DependsOn< args... >, const lambda &integrand, Domain &domain, std::shared_ptr< QuadratureData< qpt_data_type >> qdata=NoQData)
Adds a domain integral term to the weak formulation of the PDE.
auto operator()(double t, const T &... args)
This function lets the user evaluate the smith::ShapeAwareFunctional with the given trial space value...
void AddInteriorFaceIntegral(Dimension< dim >, DependsOn< args... >, const lambda &integrand, Domain &domain)
Adds an interior face integral term to the weak formulation of the PDE.
ShapeAwareFunctional(const mfem::ParFiniteElementSpace *shape_fes, std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
Constructs a QOI functional using mfem::ParFiniteElementSpace objects corresponding to the trial spac...
void updateQdata(bool update_flag)
A flag to update the quadrature data for this operator following the computation.
auto operator()(DifferentiateWRT< wrt >, double t, const T &... args)
This function lets the user evaluate the smith::ShapeAwareFunctional with the given trial space value...
ShapeAwareFunctional(const mfem::ParFiniteElementSpace *shape_fes, const mfem::ParFiniteElementSpace *test_fes, std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
Constructs using mfem::ParFiniteElementSpace objects corresponding to the test/trial spaces.
void AddBoundaryIntegral(Dimension< dim >, DependsOn< args... >, const lambda &integrand, Domain &domain)
Adds a boundary integral term to the weak formulationx of the PDE.
Implementation of the quadrature-function-based functional enabling rapid development of FEM formulat...
Accelerator functionality.
Definition: smith.cpp:36
constexpr T & get(variant< T0, T1 > &v)
Returns the variant member of specified type.
Definition: variant.hpp:338
constexpr SMITH_HOST_DEVICE auto inv(const isotropic_tensor< T, m, m > &I)
return the inverse of an isotropic tensor
SMITH_HOST_DEVICE auto cross(const tensor< T, 3, 2 > &A)
compute the cross product of the columns of A: A(:,1) x A(:,2)
Definition: tensor.hpp:966
std::shared_ptr< QuadratureData< Nothing > > NoQData
a single instance of a QuadratureData container of Nothings, since they are all interchangeable
constexpr SMITH_HOST_DEVICE auto norm(const isotropic_tensor< T, m, m > &I)
compute the Frobenius norm (sqrt(tr(dot(transpose(I), I)))) of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto transpose(const isotropic_tensor< T, m, m > &I)
return the transpose of an isotropic tensor
ExecutionSpace
enum used for signalling whether or not to perform certain calculations on the CPU or GPU
Definition: accelerator.hpp:88
SMITH_HOST_DEVICE tuple< T... > make_tuple(const T &... args)
helper function for combining a list of values into a tuple
Definition: tuple.hpp:266
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
constexpr SMITH_HOST_DEVICE isotropic_tensor< double, m, m > Identity()
return the identity matrix of the specified size
SMITH_HOST_DEVICE auto compute_boundary_area_correction(const position_type &X, const shape_type &shape)
Compute the boundary area correction term for boundary integrals with a shape displacement field.
SMITH_HOST_DEVICE auto apply_shape_aware_qf_helper_with_state(const lambda &qf, double t, const coord_type &position, state_type &state, const shape_type &shape, const space_types &space_tuple, const trial_types &arg_tuple, const correction_type &correction, std::integer_sequence< int, i... >)
A helper function to modify all of the trial function input derivatives according to the given shape ...
SMITH_HOST_DEVICE auto apply_shape_aware_qf_helper(const lambda &qf, double t, const coord_type &position, const shape_type &shape, const space_types &space_tuple, const trial_types &arg_tuple, const correction_type &correction, std::integer_sequence< int, i... >)
A helper function to modify all of the trial function input derivatives according to the given shape ...
Compile-time alias for a dimension.
Definition: geometry.hpp:17
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
these classes are a little confusing. These two special types represent the similar (but different) c...
A class for storing and access user-defined types at quadrature points.
SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val, QFuncArgs... qfunc_args) const
integrand call
SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val, QFuncArgs... qfunc_args) const
integrand call
SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val, QFuncArgs... qfunc_args) const
integrand call
SMITH_HOST_DEVICE auto operator()(double time, PositionType x, StateType &state, ShapeValueType shape_val, QFuncArgs... qfunc_args) const
integrand call
SMITH_HOST_DEVICE auto operator()(double time, PositionType x, ShapeValueType shape_val, QFuncArgs... qfunc_args) const
integrand call
A helper struct that contains the appropriate parent-to-physical and physical-to-parent transformatio...
SMITH_HOST_DEVICE auto modify_trial_argument(space_type, const trial_type &u) const
Modify the trial argument using the correct physical to reference to shape-adjusted transformation fo...
SMITH_HOST_DEVICE ShapeCorrection(Dimension< dim >, shape_type p)
Construct a new Shape Correction object with the appropriate transformations for a shape field.
SMITH_HOST_DEVICE auto modify_shape_aware_qf_return(test_space, const q_func_type &v)
Modify the quadrature function return value using the correct shape-adjusted to reference transformat...
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
Definition: tuple.hpp:28