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  std::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 time, auto X, auto V, auto... inputs) {
112  auto orig_tuple = integrand(time, 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 
224  template <int... active_parameters, typename InteriorIntegrandType>
226  InteriorIntegrandType integrand)
227  {
228  weak_form_->AddInteriorFaceIntegral(Dimension<spatial_dim - 1>{}, DependsOn<active_parameters...>{}, integrand,
229  mesh_->domain(interior_name));
230 
231  v_dot_weak_form_residual_->AddInteriorFaceIntegral(
232  Dimension<spatial_dim - 1>{}, DependsOn<0, 1 + active_parameters...>{},
233  [integrand](double t, auto X, auto V, auto... params) {
234  auto [V1, V2] = V;
235  auto orig_surface_flux = integrand(t, X, params...);
236  auto [flux_pos, flux_neg] = orig_surface_flux;
237  return smith::inner(V1, flux_pos) + smith::inner(V2, flux_neg);
238  },
239  mesh_->domain(interior_name));
240  }
241 
243  template <typename InteriorIntegrandType>
244  void addInteriorBoundaryIntegral(std::string interior_name, const InteriorIntegrandType& integrand)
245  {
246  addInteriorBoundaryIntegral(DependsOn<>{}, interior_name, integrand);
247  }
248 
270  template <int... active_parameters, typename BoundaryFluxType>
271  void addBoundaryFlux(DependsOn<active_parameters...> depends_on, std::string boundary_name,
272  BoundaryFluxType flux_function)
273  {
274  addBoundaryIntegral(depends_on, boundary_name, [flux_function](double t, auto X, auto... inputs) {
275  auto n = cross(get<DERIVATIVE>(X));
276  return -flux_function(t, get<VALUE>(X), normalize(n), get<VALUE>(inputs)...);
277  });
278  }
279 
281  template <typename BoundaryFluxType>
282  void addBoundaryFlux(std::string boundary_name, const BoundaryFluxType& integrand)
283  {
284  addBoundaryFlux(DependsOn<>{}, boundary_name, integrand);
285  }
286 
288  mfem::Vector residual(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
289  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {}) const override
290  {
291  dt_ = time_info.dt();
292  cycle_ = time_info.cycle();
293  auto ret = (*weak_form_)(time_info.time(), *shape_disp, *fields[input_indices]...);
294  return ret;
295  }
296 
298  std::unique_ptr<mfem::HypreParMatrix> jacobian(
299  TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
300  const std::vector<double>& jacobian_weights,
301  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields = {}) const override
302  {
303  dt_ = time_info.dt();
304  cycle_ = time_info.cycle();
305 
306  std::unique_ptr<mfem::HypreParMatrix> J;
307 
308  auto addToJ = [&J](double factor, std::unique_ptr<mfem::HypreParMatrix> jac_contrib) {
309  if (J) {
310  SLIC_ERROR_IF(J->N() != jac_contrib->N(),
311  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
312  SLIC_ERROR_IF(J->M() != jac_contrib->M(),
313  "Multiple nonzero jacobian weights are being used on inconsistently sized input arguments.");
314  J->Add(factor, *jac_contrib);
315  } else {
316  J.reset(jac_contrib.release());
317  if (factor != 1.0) (*J) *= factor;
318  }
319  };
320 
321  auto jacs = jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time_info.time(),
322  shape_disp, fields);
323 
324  for (size_t input_col = 0; input_col < jacobian_weights.size(); ++input_col) {
325  if (jacobian_weights[input_col] != 0.0) {
326  auto K = smith::get<DERIVATIVE>(jacs[input_col](time_info.time(), shape_disp, fields));
327  addToJ(jacobian_weights[input_col], assemble(K));
328  }
329  }
330 
331  return J;
332  }
333 
335  void jvp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
336  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields,
337  [[maybe_unused]] ConstFieldPtr v_shape_disp, const std::vector<ConstFieldPtr>& v_fields,
338  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& v_quad_fields,
339  DualFieldPtr jvp_reaction) const override
340  {
341  SLIC_ERROR_IF(v_fields.size() != fields.size(),
342  "Invalid number of field sensitivities relative to the number of fields");
343 
344  dt_ = time_info.dt();
345  cycle_ = time_info.cycle();
346 
347  auto jacs = jacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{}, time_info.time(),
348  shape_disp, fields);
349 
350  *jvp_reaction = 0.0;
351 
352  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
353  if (v_fields[input_col] != nullptr) {
354  auto K = smith::get<DERIVATIVE>(jacs[input_col](time_info.time(), shape_disp, fields));
355  K.AddMult(*v_fields[input_col], *jvp_reaction);
356  }
357  }
358  }
359 
361  void vjp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector<ConstFieldPtr>& fields,
362  [[maybe_unused]] const std::vector<ConstQuadratureFieldPtr>& quad_fields, ConstFieldPtr v_field,
363  DualFieldPtr vjp_shape_disp_sensitivity, const std::vector<DualFieldPtr>& vjp_sensitivities,
364  [[maybe_unused]] const std::vector<QuadratureFieldPtr>& vjp_quad_field_sensitivities) const override
365  {
366  SLIC_ERROR_IF(vjp_sensitivities.size() != fields.size(),
367  "Invalid number of field sensitivities relative to the number of fields");
368 
369  dt_ = time_info.dt();
370  cycle_ = time_info.cycle();
371  auto vecJacs = vectorJacobianFunctions(std::make_integer_sequence<int, sizeof...(input_indices)>{},
372  time_info.time(), shape_disp, v_field, fields);
373  {
374  auto shape_vjp = smith::get<DERIVATIVE>((*v_dot_weak_form_residual_)(
375  DifferentiateWRT<0>{}, time_info.time(), *shape_disp, *v_field, *fields[input_indices]...));
376  auto shape_vjp_vector = assemble(shape_vjp);
377  *vjp_shape_disp_sensitivity += *shape_vjp_vector;
378  }
379 
380  for (size_t input_col = 0; input_col < fields.size(); ++input_col) {
381  if (vjp_sensitivities[input_col] != nullptr) {
382  auto vec_jac = smith::get<DERIVATIVE>(vecJacs[input_col](time_info.time(), shape_disp, v_field, fields));
383  auto vec_jac_mfem_vector = assemble(vec_jac);
384  *vjp_sensitivities[input_col] += *vec_jac_mfem_vector;
385  }
386  }
387  }
388 
391  ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>& getShapeAwareResidual() { return *weak_form_; }
392 
396  ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>& getShapeAwareVectorTimesResidual()
397  {
398  return *v_dot_weak_form_residual_;
399  }
400 
401  protected:
403  template <int... i>
404  auto jacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp,
405  const std::vector<ConstFieldPtr>& fs) const
406  {
407  using JacFuncType = std::function<decltype((*weak_form_)(DifferentiateWRT<1>{}, time, *shape_disp, *fs[i]...))(
408  double, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
409  return std::array<JacFuncType, sizeof...(i)>{
410  [this](double _time, ConstFieldPtr _shape_disp, const std::vector<ConstFieldPtr>& _fs) {
411  return (*weak_form_)(DifferentiateWRT<i + 1>{}, _time, *_shape_disp, *_fs[i]...);
412  }...};
413  }
414 
416  template <int... i>
417  auto vectorJacobianFunctions(std::integer_sequence<int, i...>, double time, ConstFieldPtr shape_disp, ConstFieldPtr v,
418  const std::vector<ConstFieldPtr>& fs) const
419  {
420  using GradFuncType =
421  std::function<decltype((*v_dot_weak_form_residual_)(DifferentiateWRT<1>{}, time, *shape_disp, *v, *fs[i]...))(
422  double, ConstFieldPtr, ConstFieldPtr, const std::vector<ConstFieldPtr>&)>;
423  return std::array<GradFuncType, sizeof...(i)>{
424  [this](double _time, ConstFieldPtr _shape_disp, ConstFieldPtr _v, const std::vector<ConstFieldPtr>& _fs) {
425  return (*v_dot_weak_form_residual_)(DifferentiateWRT<i + 2>{}, _time, *_shape_disp, *_v, *_fs[i]...);
426  }...};
427  }
428 
430  mutable double dt_ = std::numeric_limits<double>::max();
431 
433  mutable size_t cycle_ = 0;
434 
436  std::shared_ptr<Mesh> mesh_;
437 
439  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, OutputSpace(InputSpaces...)>> weak_form_;
440 
442  std::unique_ptr<ShapeAwareFunctional<ShapeDispSpace, double(OutputSpace, InputSpaces...)>> v_dot_weak_form_residual_;
443 };
444 
447 inline std::vector<const mfem::ParFiniteElementSpace*> getSpaces(const std::vector<smith::FiniteElementState>& states)
448 {
449  std::vector<const mfem::ParFiniteElementSpace*> spaces;
450  for (auto& f : states) {
451  spaces.push_back(&f.space());
452  }
453  return spaces;
454 }
455 
456 } // namespace smith
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
void jvp(TimeInfo time_info, 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, DualFieldPtr jvp_reaction) 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.
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
void addInteriorBoundaryIntegral(std::string interior_name, const InteriorIntegrandType &integrand)
This is an overloaded member function, provided for convenience. It differs from the above function o...
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...
std::unique_ptr< mfem::HypreParMatrix > jacobian(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, const std::vector< double > &jacobian_weights, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}) const override
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.
mfem::Vector residual(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields={}) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addInteriorBoundaryIntegral(DependsOn< active_parameters... >, std::string interior_name, InteriorIntegrandType integrand)
Add a interior boundary integral term to the weak form.
void vjp(TimeInfo time_info, ConstFieldPtr shape_disp, const std::vector< ConstFieldPtr > &fields, [[maybe_unused]] const std::vector< ConstQuadratureFieldPtr > &quad_fields, ConstFieldPtr v_field, 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...
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.
std::vector< const mfem::ParFiniteElementSpace * > spaces(const std::vector< FieldState > &states, const std::vector< FieldState > &params={})
Get the spaces from the primal fields of a vector of field states.
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
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
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:1122
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:45
struct storing time and timestep information
Definition: common.hpp:18
double dt() const
accessor for dt
Definition: common.hpp:29
size_t cycle() const
accessor for cycle
Definition: common.hpp:32
double time() const
accessor for the current time
Definition: common.hpp:26
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.