Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
functional_weak_form.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 
17 #include "smith/physics/mesh.hpp"
21 
22 namespace smith {
23 
24 template <int spatial_dim, typename OutputSpace, typename inputs = Parameters<>,
25  typename input_indices = std::make_integer_sequence<int, inputs::n>>
27 
35 template <int spatial_dim, typename OutputSpace, typename... InputSpaces, int... input_indices>
36 class FunctionalWeakForm<spatial_dim, OutputSpace, Parameters<InputSpaces...>,
37  std::integer_sequence<int, input_indices...>> : public WeakForm {
38  public:
39  using SpacesT = std::vector<const mfem::ParFiniteElementSpace*>;
40 
42 
51  FunctionalWeakForm(std::string physics_name, std::shared_ptr<Mesh> mesh,
52  const mfem::ParFiniteElementSpace& output_mfem_space, const SpacesT& input_mfem_spaces)
53  : WeakForm(physics_name), mesh_(mesh)
54  {
55  std::array<const mfem::ParFiniteElementSpace*, sizeof...(InputSpaces)> trial_spaces;
56  std::array<const mfem::ParFiniteElementSpace*, sizeof...(InputSpaces) + 1> vector_residual_trial_spaces{
57  &output_mfem_space};
58 
59  SLIC_ERROR_ROOT_IF(
60  sizeof...(InputSpaces) != input_mfem_spaces.size(),
61  axom::fmt::format("{} parameter spaces given in the template argument but {} input mfem spaces were supplied.",
62  sizeof...(InputSpaces), input_mfem_spaces.size()));
63 
64  if constexpr (sizeof...(InputSpaces) > 0) {
65  for_constexpr<sizeof...(InputSpaces)>([&](auto i) { trial_spaces[i] = input_mfem_spaces[i]; });
66  for_constexpr<sizeof...(InputSpaces)>(
67  [&](auto i) { vector_residual_trial_spaces[i + 1] = input_mfem_spaces[i]; });
68  }
69 
70  const auto& shape_disp_space = mesh->shapeDisplacementSpace();
71 
72  weak_form_ = std::make_unique<ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>>(
73  &shape_disp_space, &output_mfem_space, trial_spaces);
74 
75  v_dot_weak_form_residual_ =
76  std::make_unique<ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>>(
77  &shape_disp_space, vector_residual_trial_spaces);
78  }
79 
103  template <int... active_parameters, typename BodyIntegralType>
104  void addBodyIntegral(DependsOn<active_parameters...>, std::string body_name, BodyIntegralType integrand)
105  {
106  weak_form_->AddDomainIntegral(Dimension<spatial_dim>{}, DependsOn<active_parameters...>{}, integrand,
107  mesh_->domain(body_name));
108 
109  v_dot_weak_form_residual_->AddDomainIntegral(
110  Dimension<spatial_dim>{}, DependsOn<0, 1 + active_parameters...>{},
111  [integrand](double t, auto X, auto V, auto... inputs) {
112  auto orig_tuple = integrand(t, X, inputs...);
113  return smith::inner(get<VALUE>(V), get<VALUE>(orig_tuple)) +
114  smith::inner(get<DERIVATIVE>(V), get<DERIVATIVE>(orig_tuple));
115  },
116  mesh_->domain(body_name));
117  }
118 
120  template <typename BodyForceType>
121  void addBodyIntegral(std::string body_name, BodyForceType body_integral)
122  {
123  addBodyIntegral(DependsOn<>{}, body_name, body_integral);
124  }
125 
146  template <int... active_parameters, typename BodyLoadType>
147  void addBodySource(DependsOn<active_parameters...> depends_on, std::string body_name, BodyLoadType load_function)
148  {
149  addBodyIntegral(depends_on, body_name, [load_function](double t, auto X, auto... inputs) {
150  return smith::tuple{-load_function(t, get<VALUE>(X), get<VALUE>(inputs)...), smith::zero{}};
151  });
152  }
153 
155  template <int... active_parameters, typename BodyLoadType>
156  void addBodySource(std::string body_name, BodyLoadType load_function)
157  {
158  return addBodySource(DependsOn<>{}, body_name, load_function);
159  }
160 
184  template <int... active_parameters, typename BoundaryIntegrandType>
185  void addBoundaryIntegral(DependsOn<active_parameters...>, std::string boundary_name, BoundaryIntegrandType integrand)
186  {
187  weak_form_->AddBoundaryIntegral(Dimension<spatial_dim - 1>{}, DependsOn<active_parameters...>{}, integrand,
188  mesh_->domain(boundary_name));
189 
190  v_dot_weak_form_residual_->AddBoundaryIntegral(
191  Dimension<spatial_dim - 1>{}, DependsOn<0, 1 + active_parameters...>{},
192  [integrand](double t, auto X, auto V, auto... params) {
193  auto orig_surface_flux = integrand(t, X, params...);
194  return smith::inner(get<VALUE>(V), orig_surface_flux);
195  },
196  mesh_->domain(boundary_name));
197  }
198 
200  template <typename BoundaryIntegrandType>
201  void addBoundaryIntegral(std::string boundary_name, const BoundaryIntegrandType& integrand)
202  {
203  addBoundaryIntegral(DependsOn<>{}, boundary_name, integrand);
204  }
205 
227  template <int... active_parameters, typename BoundaryFluxType>
228  void addBoundaryFlux(DependsOn<active_parameters...> depends_on, std::string boundary_name,
229  BoundaryFluxType flux_function)
230  {
231  addBoundaryIntegral(depends_on, boundary_name, [flux_function](double t, auto X, auto... inputs) {
232  auto n = cross(get<DERIVATIVE>(X));
233  return -flux_function(t, get<VALUE>(X), normalize(n), get<VALUE>(inputs)...);
234  });
235  }
236 
238  template <typename BoundaryFluxType>
239  void addBoundaryFlux(std::string boundary_name, const BoundaryFluxType& integrand)
240  {
241  addBoundaryFlux(DependsOn<>{}, boundary_name, integrand);
242  }
243 
245  mfem::Vector residual(double time, double dt, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
246  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {},
247  int block_row = 0) const override
248  {
249  SLIC_ERROR_IF(block_row != 0, "Invalid block row and column requested in fieldJacobian for FunctionalResidual");
250  dt_ = dt;
251  auto ret = (*weak_form_)(time, *shape_disp, *fields[input_indices]...);
252  return ret;
253  }
254 
256  std::unique_ptr<mfem::HypreParMatrix> jacobian(
257  double time, double dt, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
258  const std::vector<double>& jacobian_weights,
259  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {}, int block_row = 0) const override
260  {
261  SLIC_ERROR_IF(block_row != 0, "Invalid block row and column requested in fieldJacobian for FunctionalResidual");
262 
263  dt_ = dt;
264 
265  std::unique_ptr<mfem::HypreParMatrix> J;
266 
267  auto addToJ = [&J](double factor, std::unique_ptr<mfem::HypreParMatrix> jac_contrib) {
268  if (J) {
269  SLIC_ERROR_IF(J->N() != jac_contrib->N(),
270  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
271  SLIC_ERROR_IF(J->M() != jac_contrib->M(),
272  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
273  J->Add(factor, *jac_contrib);
274  } else {
275  J.reset(jac_contrib.release());
276  if (factor != 1.0) (*J) *= factor;
277  }
278  };
279 
280  auto jacs =
281  jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time, shape_disp, fields);
282 
283  for (size_t input_col = 0; input_col < jacobian_weights.size(); ++input_col) {
284  if (jacobian_weights[input_col] != 0.0) {
285  auto K = smith::get<DERIVATIVE>(jacs[input_col](time, shape_disp, fields));
286  addToJ(jacobian_weights[input_col], assemble(K));
287  }
288  }
289 
290  return J;
291  }
292 
294  void jvp(double time, double dt, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
295  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields,
296  [[maybe_unused]] ConstFieldPtr v_shape_disp, const std::vector<ConstFieldPtr>& v_fields,
297  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& v_quad_fields,
298  const std::vector<DualFieldPtr>& jvp_reactions) const override
299  {
300  SLIC_ERROR_IF(v_fields.size() != fields.size(),
301  "Invalid number of field sensitivities relative to the number of fields");
302  SLIC_ERROR_IF(jvp_reactions.size() != 1, "FunctionalResidual nonlinear systems only supports 1 output residual");
303 
304  dt_ = dt;
305  auto jacs =
306  jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time, shape_disp, fields);
307 
308  *jvp_reactions[0] = 0.0;
309 
310  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
311  if (v_fields[input_col] != nullptr) {
312  auto K = smith::get<DERIVATIVE>(jacs[input_col](time, shape_disp, fields));
313  K.AddMult(*v_fields[input_col], *jvp_reactions[0]);
314  }
315  }
316  }
317 
319  void vjp(double time, double dt, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
320  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields,
321  const std::vector<ConstFieldPtr>& v_fields, DualFieldPtr vjp_shape_disp_sensitivity,
322  const std::vector<DualFieldPtr>& vjp_sensitivities,
323  [[maybe_unused]] const std::vector<QuadratureFieldPtr>& vjp_quad_field_sensitivities) const override
324  {
325  SLIC_ERROR_IF(vjp_sensitivities.size() != fields.size(),
326  "Invalid number of field sensitivities relative to the number of fields");
327  SLIC_ERROR_IF(v_fields.size() != 1, "FunctionalResidual nonlinear systems only supports 1 output residual");
328 
329  dt_ = dt;
330  auto vecJacs = vectorJacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time,
331  shape_disp, v_fields[0], fields);
332  {
333  auto shape_vjp = smith::get<DERIVATIVE>((*v_dot_weak_form_residual_)(DifferentiateWRT<0>{}, time, *shape_disp,
334  *v_fields[0], *fields[input_indices]...));
335  auto shape_vjp_vector = assemble(shape_vjp);
336  *vjp_shape_disp_sensitivity += *shape_vjp_vector;
337  }
338 
339  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
340  if (vjp_sensitivities[input_col] != nullptr) {
341  auto vec_jac = smith::get<DERIVATIVE>(vecJacs[input_col](time, shape_disp, v_fields[0], fields));
342  auto vec_jac_mfem_vector = assemble(vec_jac);
343  *vjp_sensitivities[input_col] += *vec_jac_mfem_vector;
344  }
345  }
346  }
347 
350  ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>& getShapeAwareResidual() { return *weak_form_; }
351 
355  ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>& getShapeAwareVectorTimesResidual()
356  {
357  return *v_dot_weak_form_residual_;
358  }
359 
360  protected:
362  template <int... i>
363  auto jacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp,
364  const std::vector<ConstFieldPtr>& fs) const
365  {
366  using JacFuncType = std::function<decltype((*weak_form_)(DifferentiateWRT<1>{}, time, *shape_disp, *fs[i]...))(
367  double, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
368  return std::array<JacFuncType, sizeof...(i)>{
369  [=](double _time, ConstFieldPtr _shape_disp, const std::vector<ConstFieldPtr>& _fs) {
370  return (*weak_form_)(DifferentiateWRT<i + 1>{}, _time, *_shape_disp, *_fs[i]...);
371  }...};
372  };
373 
375  template <int... i>
376  auto vectorJacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp, ConstFieldPtr v,
377  const std::vector<ConstFieldPtr>& fs) const
378  {
379  using GradFuncType =
380  std::function<decltype((*v_dot_weak_form_residual_)(DifferentiateWRT<1>{}, time, *shape_disp, *v, *fs[i]...))(
381  double, ConstFieldPtr, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
382  return std::array<GradFuncType, sizeof...(i)>{
383  [=](double _time, ConstFieldPtr _shape_disp, ConstFieldPtr _v, const std::vector<ConstFieldPtr>& _fs) {
384  return (*v_dot_weak_form_residual_)(DifferentiateWRT<i + 2>{}, _time, *_shape_disp, *_v, *_fs[i]...);
385  }...};
386  };
387 
389  mutable double dt_ = std::numeric_limits<double>::max();
390 
392  std::shared_ptr<Mesh> mesh_;
393 
395  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>> weak_form_;
396 
398  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>> v_dot_weak_form_residual_;
399 };
400 
403 inline std::vector<const mfem::ParFiniteElementSpace*> getSpaces(const std::vector<smith::FiniteElementState>& states)
404 {
405  std::vector<const mfem::ParFiniteElementSpace*> spaces;
406  for (auto& f : states) {
407  spaces.push_back(&f.space());
408  }
409  return spaces;
410 }
411 
412 } // namespace smith
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
mfem::Vector residual(double time, double dt, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}, int block_row=0) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void vjp(double time, double dt, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields, const std::vector< ConstFieldPtr > &v_fields, DualFieldPtr vjp_shape_disp_sensitivity, const std::vector< DualFieldPtr > &vjp_sensitivities, [[maybe_unused]] const std::vector< QuadratureFieldPtr > &vjp_quad_field_sensitivities) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
auto vectorJacobianFunctions(std::integer_sequence< int, i... >, double time, ConstFieldPtr shape_disp, ConstFieldPtr v, const std::vector< ConstFieldPtr > &fs) const
Utility to get array of jvp functions, one for each input field in fs.
void jvp(double time, double dt, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields, [[maybe_unused]] ConstFieldPtr v_shape_disp, const std::vector< ConstFieldPtr > &v_fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &v_quad_fields, const std::vector< DualFieldPtr > &jvp_reactions) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< ShapeAwareFunctional< ShapeDispSpace, double(OutputSpace, InputSpaces...)> > v_dot_weak_form_residual_
functional residual times and arbitrary vector v (same space as residual) evaluator,...
void addBodyIntegral(std::string body_name, BodyForceType body_integral)
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< ShapeAwareFunctional< ShapeDispSpace, OutputSpace(InputSpaces...)> > weak_form_
functional residual evaluator, shape aware
ShapeAwareFunctional< ShapeDispSpace, OutputSpace(InputSpaces...)> & getShapeAwareResidual()
Accessor to get a reference to the underlying ShapeAwareFunctional in case more direct access is need...
FunctionalWeakForm(std::string physics_name, std::shared_ptr< Mesh > mesh, const mfem::ParFiniteElementSpace &output_mfem_space, const SpacesT &input_mfem_spaces)
Construct a new FunctionalWeakForm object.
void addBodySource(std::string body_name, BodyLoadType load_function)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBoundaryIntegral(DependsOn< active_parameters... >, std::string boundary_name, BoundaryIntegrandType integrand)
Add a boundary integral term to the weak form.
auto jacobianFunctions(std::integer_sequence< int, i... >, double time, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fs) const
Utility to get array of jacobian functions, one for each input field in fs.
void addBoundaryIntegral(std::string boundary_name, const BoundaryIntegrandType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBoundaryFlux(std::string boundary_name, const BoundaryFluxType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBodySource(DependsOn< active_parameters... > depends_on, std::string body_name, BodyLoadType load_function)
Add a body source (body load) to the weak form.
void addBoundaryFlux(DependsOn< active_parameters... > depends_on, std::string boundary_name, BoundaryFluxType flux_function)
Add a boundary flux term to the weak form.
std::unique_ptr< mfem::HypreParMatrix > jacobian(double time, double dt, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, const std::vector< double > &jacobian_weights, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}, int block_row=0) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBodyIntegral(DependsOn< active_parameters... >, std::string body_name, BodyIntegralType integrand)
Add a body integral contribution to the residual.
ShapeAwareFunctional< ShapeDispSpace, double(OutputSpace, InputSpaces...)> & getShapeAwareVectorTimesResidual()
Accessor to get a reference to the underlying ShapeAwareFunctional vector-residual in case more direc...
Abstract WeakForm class.
Definition: weak_form.hpp:36
This contains a class that represents the dual of a finite element vector space, i....
This file contains the declaration of structure that manages the MFEM objects that make up the state ...
Smith mesh class which assists in constructing the appropriate parallel mfem meshes and registering a...
constexpr SMITH_HOST_DEVICE void for_constexpr(const lambda &f)
multidimensional loop tool that evaluates the lambda body inside the innermost loop.
Accelerator functionality.
Definition: smith.cpp:36
std::vector< const mfem::ParFiniteElementSpace * > getSpaces(const std::vector< smith::FiniteElementState > &states)
Helper function to construct vector of spaces from an existing vector of FiniteElementState.
SMITH_HOST_DEVICE auto max(dual< gradient_type > a, double b)
Implementation of max for dual numbers.
Definition: dual.hpp:229
constexpr SMITH_HOST_DEVICE auto inner(const dual< S > &A, const dual< T > &B)
Definition: dual.hpp:281
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:959
SMITH_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
Definition: tensor.hpp:1115
FiniteElementState const * ConstFieldPtr
using
Definition: field_types.hpp:36
Wrapper of smith::Functional for evaluating integrals and derivatives of quantities with shape displa...
Compile-time alias for a dimension.
Definition: geometry.hpp:17
H1 elements of order p.
a struct that is used in the physics modules to clarify which template arguments are user-controlled ...
Definition: common.hpp:21
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
Definition: tuple.hpp:28
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122
Specifies interface for evaluating weak form residuals and their gradients.