Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
solid_mechanics.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Serac Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
15 #include <cstddef>
16 #include <array>
17 #include <functional>
18 #include <memory>
19 #include <optional>
20 #include <string>
21 #include <unordered_map>
22 #include <utility>
23 #include <vector>
24 
25 #include "mfem.hpp"
26 
27 #include "serac/serac_config.hpp"
28 #include "serac/physics/common.hpp"
31 #include "serac/physics/boundary_conditions/components.hpp"
32 #include "serac/numerics/odes.hpp"
36 #include "serac/physics/mesh.hpp"
42 #include "serac/numerics/functional/differentiate_wrt.hpp"
44 #include "serac/numerics/functional/geometry.hpp"
48 #include "serac/numerics/petsc_solvers.hpp"
54 
55 namespace serac {
56 namespace solid_mechanics {
57 
58 namespace detail {
62 void adjoint_integrate(double dt_n, double dt_np1, mfem::HypreParMatrix* m_mat, mfem::HypreParMatrix* k_mat,
63  mfem::HypreParVector& disp_adjoint_load_vector, mfem::HypreParVector& velo_adjoint_load_vector,
64  mfem::HypreParVector& accel_adjoint_load_vector, mfem::HypreParVector& adjoint_displacement_,
65  mfem::HypreParVector& implicit_sensitivity_displacement_start_of_step_,
66  mfem::HypreParVector& implicit_sensitivity_velocity_start_of_step_,
67  mfem::HypreParVector& adjoint_essential_, BoundaryConditionManager& bcs_,
68  mfem::Solver& lin_solver);
69 } // namespace detail
70 
77  .preconditioner = serac::ordering == mfem::Ordering::byVDIM
80  .relative_tol = 1.0e-6,
81  .absolute_tol = 1.0e-16,
82  .max_iterations = 500,
83  .print_level = 0};
84 
86 #ifdef MFEM_USE_STRUMPACK
88 #else
90 #endif
91 
98  .relative_tol = 1.0e-4,
99  .absolute_tol = 1.0e-8,
100  .max_iterations = 10,
101  .print_level = 1};
102 
105 
109 
110 } // namespace solid_mechanics
111 
112 template <int order, int dim, typename parameters = Parameters<>,
113  typename parameter_indices = std::make_integer_sequence<int, parameters::n>>
115 
126 template <int order, int dim, typename... parameter_space, int... parameter_indices>
127 class SolidMechanics<order, dim, Parameters<parameter_space...>, std::integer_sequence<int, parameter_indices...>>
128  : public BasePhysics {
129  public:
131  static constexpr int VALUE = 0, DERIVATIVE = 1;
132  static constexpr int SHAPE = 0;
133  static constexpr auto I = Identity<dim>();
135 
138  static constexpr auto NUM_STATE_VARS = 2;
139 
142  template <typename T>
143  using qdata_type = std::shared_ptr<QuadratureData<T>>;
144 
163  SolidMechanics(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts,
164  const serac::TimesteppingOptions timestepping_opts, const std::string& physics_name,
165  std::shared_ptr<serac::Mesh> serac_mesh, std::vector<std::string> parameter_names = {}, int cycle = 0,
166  double time = 0.0, bool checkpoint_to_disk = false, bool use_warm_start = true)
167  : SolidMechanics(std::make_unique<EquationSolver>(nonlinear_opts, lin_opts, serac_mesh->getComm()),
168  timestepping_opts, physics_name, serac_mesh, parameter_names, cycle, time, checkpoint_to_disk,
169  use_warm_start)
170  {
171  }
172 
190  SolidMechanics(std::unique_ptr<serac::EquationSolver> solver, const serac::TimesteppingOptions timestepping_opts,
191  const std::string& physics_name, std::shared_ptr<serac::Mesh> serac_mesh,
192  std::vector<std::string> parameter_names = {}, int cycle = 0, double time = 0.0,
193  bool checkpoint_to_disk = false, bool use_warm_start = true)
194  : BasePhysics(physics_name, serac_mesh, cycle, time, checkpoint_to_disk),
195  displacement_(
196  StateManager::newState(H1<order, dim>{}, detail::addPrefix(physics_name, "displacement"), mesh_->tag())),
197  velocity_(StateManager::newState(H1<order, dim>{}, detail::addPrefix(physics_name, "velocity"), mesh_->tag())),
198  acceleration_(
199  StateManager::newState(H1<order, dim>{}, detail::addPrefix(physics_name, "acceleration"), mesh_->tag())),
200  adjoint_displacement_(StateManager::newState(
201  H1<order, dim>{}, detail::addPrefix(physics_name, "adjoint_displacement"), mesh_->tag())),
202  displacement_adjoint_load_(displacement_.space(), detail::addPrefix(physics_name, "displacement_adjoint_load")),
203  velocity_adjoint_load_(displacement_.space(), detail::addPrefix(physics_name, "velocity_adjoint_load")),
204  acceleration_adjoint_load_(displacement_.space(), detail::addPrefix(physics_name, "acceleration_adjoint_load")),
205  implicit_sensitivity_displacement_start_of_step_(displacement_.space(), "total_deriv_wrt_displacement."),
206  implicit_sensitivity_velocity_start_of_step_(displacement_.space(), "total_deriv_wrt_velocity."),
207  reactions_(StateManager::newDual(H1<order, dim>{}, detail::addPrefix(physics_name, "reactions"), mesh_->tag())),
208  reactions_adjoint_bcs_(reactions_.space(), "reactions_shape_sensitivity"),
209  nonlin_solver_(std::move(solver)),
210  ode2_(displacement_.space().TrueVSize(),
211  {.time = time_, .c0 = c0_, .c1 = c1_, .u = u_, .du_dt = v_, .d2u_dt2 = acceleration_}, *nonlin_solver_,
212  bcs_),
213  use_warm_start_(use_warm_start)
214  {
216  SLIC_ERROR_ROOT_IF(mfemParMesh().Dimension() != dim,
217  axom::fmt::format("Compile time dimension, {0}, and runtime mesh dimension, {1}, mismatch", dim,
218  mfemParMesh().Dimension()));
219 
220  SLIC_ERROR_ROOT_IF(!nonlin_solver_,
221  "EquationSolver argument is nullptr in SolidMechanics constructor. It is possible that it was "
222  "previously moved.");
223 
224  // Check for dynamic mode
225  if (timestepping_opts.timestepper != TimestepMethod::QuasiStatic) {
226  ode2_.SetTimestepper(timestepping_opts.timestepper);
227  ode2_.SetEnforcementMethod(timestepping_opts.enforcement_method);
228  is_quasistatic_ = false;
229  } else {
230  is_quasistatic_ = true;
231  }
232 
233  states_.push_back(&displacement_);
234  if (!is_quasistatic_) {
235  states_.push_back(&velocity_);
236  states_.push_back(&acceleration_);
237  }
238 
239  adjoints_.push_back(&adjoint_displacement_);
240  duals_.push_back(&reactions_);
241  dual_adjoints_.push_back(&reactions_adjoint_bcs_);
242 
243  // Create a pack of the primal field and parameter finite element spaces
244  const mfem::ParFiniteElementSpace* test_space = &displacement_.space();
245  const mfem::ParFiniteElementSpace* shape_space = &mesh_->shapeDisplacementSpace();
246 
247  std::array<const mfem::ParFiniteElementSpace*, NUM_STATE_VARS + sizeof...(parameter_space)> trial_spaces;
248  trial_spaces[0] = &displacement_.space();
249  trial_spaces[1] = &displacement_.space();
250 
251  SLIC_ERROR_ROOT_IF(
252  sizeof...(parameter_space) != parameter_names.size(),
253  axom::fmt::format("{} parameter spaces given in the template argument but {} parameter names were supplied.",
254  sizeof...(parameter_space), parameter_names.size()));
255 
256  if constexpr (sizeof...(parameter_space) > 0) {
257  tuple<parameter_space...> types{};
258  for_constexpr<sizeof...(parameter_space)>([&](auto i) {
259  parameters_.emplace_back(mfemParMesh(), get<i>(types), detail::addPrefix(name_, parameter_names[i]));
260  trial_spaces[i + NUM_STATE_VARS] = &(parameters_[i].state->space());
261  });
262  }
263 
264  residual_ = std::make_unique<ShapeAwareFunctional<shape_trial, test(trial, trial, parameter_space...)>>(
265  shape_space, test_space, trial_spaces);
266 
267  // If the user wants the AMG preconditioner with a linear solver, set the pfes
268  // to be the displacement
269  auto* amg_prec = dynamic_cast<mfem::HypreBoomerAMG*>(&nonlin_solver_->preconditioner());
270  if (amg_prec) {
271  // ZRA - Iterative refinement tends to be more expensive than it is worth
272  // We should add a flag allowing users to enable it
273 
274  // bool iterative_refinement = false;
275  // amg_prec->SetElasticityOptions(&displacement_.space(), iterative_refinement);
276 
277  // SetElasticityOptions only works with byVDIM ordering, some evidence that it is not often optimal
278  amg_prec->SetSystemsOptions(displacement_.space().GetVDim(), serac::ordering == mfem::Ordering::byNODES);
279  }
280 
281  int true_size = velocity_.space().TrueVSize();
282 
283  u_.SetSize(true_size);
284  v_.SetSize(true_size);
285  du_.SetSize(true_size);
286  predicted_displacement_.SetSize(true_size);
287 
288  initializeSolidMechanicsStates();
289  }
290 
292  virtual ~SolidMechanics() {}
293 
299  {
300  c0_ = 0.0;
301  c1_ = 0.0;
302 
303  time_end_step_ = 0.0;
304 
305  displacement_ = 0.0;
306  velocity_ = 0.0;
307  acceleration_ = 0.0;
308 
309  adjoint_displacement_ = 0.0;
310  displacement_adjoint_load_ = 0.0;
311  velocity_adjoint_load_ = 0.0;
312  acceleration_adjoint_load_ = 0.0;
313 
314  implicit_sensitivity_displacement_start_of_step_ = 0.0;
315  implicit_sensitivity_velocity_start_of_step_ = 0.0;
316 
317  reactions_ = 0.0;
318  reactions_adjoint_bcs_ = 0.0;
319 
320  u_ = 0.0;
321  v_ = 0.0;
322  du_ = 0.0;
323  predicted_displacement_ = 0.0;
324  }
325 
327  void resetStates(int cycle = 0, double time = 0.0) override
328  {
330  initializeSolidMechanicsStates();
331 
332  if (checkpoint_to_disk_) {
333  outputStateToDisk();
334  } else {
335  checkpoint_states_.clear();
336  auto state_names = stateNames();
337  for (const auto& state_name : state_names) {
338  checkpoint_states_[state_name].push_back(state(state_name));
339  }
340  }
341  }
342 
351  template <typename T>
352  qdata_type<T> createQuadratureDataBuffer(T initial_state, const std::optional<Domain>& optional_domain = std::nullopt)
353  {
354  Domain domain = (optional_domain) ? *optional_domain : EntireDomain(mfemParMesh());
355  return StateManager::newQuadratureDataBuffer(mesh_->tag(), domain, order, dim, initial_state);
356  }
357 
387  template <typename AppliedDisplacementFunction>
388  void setDisplacementBCs(AppliedDisplacementFunction applied_displacement, Domain& domain,
389  Components components = Component::ALL)
390  {
391  for (int i = 0; i < dim; ++i) {
392  if (components[size_t(i)]) {
393  auto mfem_coefficient_function = [applied_displacement, i](const mfem::Vector& X_mfem, double t) {
394  auto X = make_tensor<dim>([&X_mfem](int k) { return X_mfem[k]; });
395  return applied_displacement(X, t)[i];
396  };
397 
398  component_disp_bdr_coef_ = std::make_shared<mfem::FunctionCoefficient>(mfem_coefficient_function);
399 
400  auto dof_list = domain.dof_list(&displacement_.space());
401 
402  // scalar ldofs -> vector ldofs
403  displacement_.space().DofsToVDofs(i, dof_list);
404 
405  bcs_.addEssential(dof_list, component_disp_bdr_coef_, displacement_.space(), i);
406  }
407  }
408  }
409 
419  void setFixedBCs(Domain& domain, Components components = Component::ALL)
420  {
421  auto zero_vector_function = [](tensor<double, dim>, double) { return tensor<double, dim>{}; };
422  setDisplacementBCs(zero_vector_function, domain, components);
423  }
424 
448  template <typename AppliedDisplacementFunction>
449  void setDisplacementBCsByDofList(AppliedDisplacementFunction applied_displacement, const mfem::Array<int> true_dofs)
450  {
451  // std::function<void(const mfem::Vector&, double, mfem::Vector&)> disp
452  auto mfem_vector_coefficient_function = [applied_displacement](const mfem::Vector& X_mfem, double t,
453  mfem::Vector& u_mfem) {
454  auto X = make_tensor<dim>([&X_mfem](int k) { return X_mfem[k]; });
455  auto u = applied_displacement(X, t);
456  for (int i = 0; i < dim; i++) {
457  u_mfem(i) = u[i];
458  }
459  };
460  disp_bdr_coef_ = std::make_shared<mfem::VectorFunctionCoefficient>(dim, mfem_vector_coefficient_function);
461  bcs_.addEssentialByTrueDofs(true_dofs, disp_bdr_coef_, displacement_.space());
462  }
463 
465  const FiniteElementState& state(const std::string& state_name) const override
466  {
467  if (state_name == "displacement") {
468  return displacement_;
469  } else if (state_name == "velocity") {
470  return velocity_;
471  } else if (state_name == "acceleration") {
472  return acceleration_;
473  }
474 
475  SLIC_ERROR_ROOT(axom::fmt::format("State '{}' requested from solid mechanics module '{}', but it doesn't exist",
476  state_name, name_));
477  return displacement_;
478  }
479 
489  void setState(const std::string& state_name, const FiniteElementState& state) override
490  {
491  if (state_name == "displacement") {
492  displacement_ = state;
493  if (!checkpoint_to_disk_) {
494  checkpoint_states_["displacement"][static_cast<size_t>(cycle_)] = displacement_;
495  }
496  return;
497  } else if (state_name == "velocity") {
498  velocity_ = state;
499  if (!checkpoint_to_disk_) {
500  checkpoint_states_["velocity"][static_cast<size_t>(cycle_)] = velocity_;
501  }
502  return;
503  }
504 
505  SLIC_ERROR_ROOT(axom::fmt::format(
506  "setState for state named '{}' requested from solid mechanics module '{}', but it doesn't exist", state_name,
507  name_));
508  }
509 
511  std::vector<std::string> stateNames() const override
512  {
513  if (is_quasistatic_) {
514  return std::vector<std::string>{"displacement"};
515  } else {
516  return std::vector<std::string>{"displacement", "velocity", "acceleration"};
517  }
518  }
519 
541  template <int... active_parameters, typename callable>
543  const std::optional<Domain>& optional_domain = std::nullopt)
544  {
545  Domain domain = (optional_domain) ? *optional_domain : EntireBoundary(mfemParMesh());
546 
547  residual_->AddBoundaryIntegral(Dimension<dim - 1>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
548  qfunction, domain);
549  }
550 
552  std::vector<std::string> adjointNames() const override { return std::vector<std::string>{{"displacement"}}; }
553 
555  const FiniteElementState& adjoint(const std::string& state_name) const override
556  {
557  if (state_name == "displacement") {
558  return adjoint_displacement_;
559  }
560 
561  SLIC_ERROR_ROOT(axom::fmt::format("adjoint '{}' requested from solid mechanics module '{}', but it doesn't exist",
562  state_name, name_));
563  return adjoint_displacement_;
564  }
565 
567  std::vector<std::string> dualNames() const override { return std::vector<std::string>{{"reactions"}}; }
568 
570  const FiniteElementDual& dual(const std::string& dual_name) const override
571  {
572  if (dual_name == "reactions") {
573  return reactions_;
574  }
575 
576  SLIC_ERROR_ROOT(axom::fmt::format("dual '{}' requested from solid mechanics module '{}', but it doesn't exist",
577  dual_name, name_));
578  return reactions_;
579  }
580 
582  const FiniteElementState& dualAdjoint(const std::string& dual_name) const override
583  {
584  if (dual_name == "reactions") {
585  return reactions_adjoint_bcs_;
586  }
587 
588  SLIC_ERROR_ROOT(axom::fmt::format(
589  "dualAdjoint '{}' requested from solid mechanics module '{}', but it doesn't exist", dual_name, name_));
590  return reactions_adjoint_bcs_;
591  }
592 
594  FiniteElementDual loadCheckpointedDual(const std::string& dual_name, int cycle) override
595  {
596  if (dual_name == "reactions") {
597  auto checkpointed_sol = getCheckpointedStates(cycle);
598 
599  FiniteElementDual reactions(reactions_.space(), "reactions_tmp");
600 
601  if (is_quasistatic_) {
602  reactions = (*residual_)(time_, shapeDisplacement(), checkpointed_sol.at("displacement"), acceleration_,
603  *parameters_[parameter_indices].state...);
604  } else {
605  reactions = (*residual_)(time_, shapeDisplacement(), checkpointed_sol.at("displacement"),
606  checkpointed_sol.at("acceleration"), *parameters_[parameter_indices].state...);
607  }
608 
609  return reactions;
610  }
611 
612  SLIC_ERROR_ROOT(
613  axom::fmt::format("loadCheckpointedDual '{}' requested from solid mechanics module '{}', but it doesn't exist",
614  dual_name, name_));
615  return reactions_;
616  }
617 
648  template <int... active_parameters, typename callable, typename StateType = Nothing>
651  {
652  residual_->AddDomainIntegral(Dimension<dim>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{}, qfunction,
653  domain, qdata);
654  }
655 
660  template <typename Material>
661  struct MaterialStressFunctor {
663  MaterialStressFunctor(Material material) : material_(material) {}
664 
666  Material material_;
667 
683  template <typename X, typename State, typename Displacement, typename Acceleration, typename... Params>
684  auto SERAC_HOST_DEVICE operator()(double, X, State& state, Displacement displacement, Acceleration acceleration,
685  Params... params) const
686  {
687  auto du_dX = get<DERIVATIVE>(displacement);
688  auto d2u_dt2 = get<VALUE>(acceleration);
689 
690  auto stress = material_(state, du_dX, params...);
691 
692  return serac::tuple{material_.density * d2u_dt2, stress};
693  }
694  };
695 
721  template <int... active_parameters, typename MaterialType, typename StateType = Empty>
722  void setMaterial(DependsOn<active_parameters...>, const MaterialType& material, Domain& domain,
724  {
725  static_assert(std::is_same_v<StateType, Empty> || std::is_same_v<StateType, typename MaterialType::State>,
726  "invalid quadrature data provided in setMaterial()");
727  MaterialStressFunctor<MaterialType> material_functor(material);
728  residual_->AddDomainIntegral(
729  Dimension<dim>{},
730  DependsOn<0, 1,
731  active_parameters + NUM_STATE_VARS...>{}, // the magic number "+ NUM_STATE_VARS" accounts for the
732  // fact that the displacement, acceleration, and shape
733  // fields are always-on and come first, so the `n`th
734  // parameter will actually be argument `n + NUM_STATE_VARS`
735  std::move(material_functor), domain, qdata);
736  }
737 
739  template <typename MaterialType, typename StateType = Empty>
740  void setMaterial(const MaterialType& material, Domain& domain,
741  std::shared_ptr<QuadratureData<StateType>> qdata = EmptyQData)
742  {
743  setMaterial(DependsOn<>{}, material, domain, qdata);
744  }
745 
752  template <typename Material>
753  struct RateDependentMaterialStressFunctor {
755  RateDependentMaterialStressFunctor(Material material, const double* dt) : material_(material), time_increment_(dt)
756  {
757  }
758 
760  Material material_;
761 
763  const double* time_increment_;
764 
780  template <typename X, typename State, typename Displacement, typename Acceleration, typename... Params>
781  auto SERAC_HOST_DEVICE operator()(double, X, State& state, Displacement displacement, Acceleration acceleration,
782  Params... params) const
783  {
784  auto du_dX = get<DERIVATIVE>(displacement);
785  auto d2u_dt2 = get<VALUE>(acceleration);
786 
787  auto stress = material_(state, *time_increment_, du_dX, params...);
788 
789  return serac::tuple{material_.density * d2u_dt2, stress};
790  }
791  };
792 
796  template <int... active_parameters, typename RateDependentMaterialType, typename StateType = Empty>
797  void setRateDependentMaterial(DependsOn<active_parameters...>, const RateDependentMaterialType& material,
798  Domain& domain, qdata_type<StateType> qdata = EmptyQData)
799  {
800  static_assert(
801  std::is_same_v<StateType, Empty> || std::is_same_v<StateType, typename RateDependentMaterialType::State>,
802  "invalid quadrature data provided in setMaterial()");
803  RateDependentMaterialStressFunctor<RateDependentMaterialType> material_functor(material, &dt_);
804  residual_->AddDomainIntegral(
805  Dimension<dim>{},
806  DependsOn<0, 1,
807  active_parameters + NUM_STATE_VARS...>{}, // the magic number "+ NUM_STATE_VARS" accounts for the
808  // fact that the displacement, acceleration, and shape
809  // fields are always-on and come first, so the `n`th
810  // parameter will actually be argument `n + NUM_STATE_VARS`
811  std::move(material_functor), domain, qdata);
812  }
813 
815  template <typename RateDependentMaterialType, typename StateType = Empty>
816  void setRateDependentMaterial(const RateDependentMaterialType& material, Domain& domain,
817  std::shared_ptr<QuadratureData<StateType>> qdata = EmptyQData)
818  {
819  setRateDependentMaterial(DependsOn<>{}, material, domain, qdata);
820  }
821 
834  template <typename AppliedDisplacementFunction>
835  void setDisplacement(AppliedDisplacementFunction applied_displacement)
836  {
837  displacement_.setFromFieldFunction(applied_displacement);
838  }
839 
841  void setDisplacement(const FiniteElementState& temp) { displacement_ = temp; }
842 
855  template <typename AppliedVelocityFunction>
856  void setVelocity(AppliedVelocityFunction applied_velocity)
857  {
858  velocity_.setFromFieldFunction(applied_velocity);
859  }
860 
862  void setVelocity(const FiniteElementState& temp) { velocity_ = temp; }
863 
868  template <typename BodyForceType>
869  struct BodyForceIntegrand {
871  BodyForceType body_force_;
873  BodyForceIntegrand(BodyForceType body_force) : body_force_(body_force) {}
874 
889  template <typename T, typename Position, typename Displacement, typename Acceleration, typename... Params>
890  auto SERAC_HOST_DEVICE operator()(T t, Position position, Displacement, Acceleration, Params... params) const
891  {
892  return serac::tuple{-1.0 * body_force_(get<VALUE>(position), t, params...), zero{}};
893  }
894  };
895 
915  template <int... active_parameters, typename BodyForceType>
916  void addBodyForce(DependsOn<active_parameters...>, BodyForceType body_force, Domain& domain)
917  {
918  residual_->AddDomainIntegral(Dimension<dim>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
919  BodyForceIntegrand<BodyForceType>(body_force), domain);
920  }
921 
923  template <typename BodyForceType>
924  void addBodyForce(BodyForceType body_force, Domain& domain)
925  {
926  addBodyForce(DependsOn<>{}, body_force, domain);
927  }
928 
952  template <int... active_parameters, typename TractionType>
953  void setTraction(DependsOn<active_parameters...>, TractionType traction_function, Domain& domain)
954  {
955  residual_->AddBoundaryIntegral(
956  Dimension<dim - 1>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
957  [traction_function](double t, auto X, auto /* displacement */, auto /* acceleration */, auto... params) {
958  auto n = cross(get<DERIVATIVE>(X));
959 
960  return -1.0 * traction_function(get<VALUE>(X), normalize(n), t, params...);
961  },
962  domain);
963  }
964 
966  template <typename TractionType>
967  void setTraction(TractionType traction_function, Domain& domain)
968  {
969  setTraction(DependsOn<>{}, traction_function, domain);
970  }
971 
995  template <int... active_parameters, typename PressureType>
996  void setPressure(DependsOn<active_parameters...>, PressureType pressure_function, Domain& domain)
997  {
998  residual_->AddBoundaryIntegral(
999  Dimension<dim - 1>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
1000  [pressure_function](double t, auto X, auto displacement, auto /* acceleration */, auto... params) {
1001  // Calculate the position and normal in the shape perturbed deformed configuration
1002  auto x = X + displacement;
1003 
1004  auto n = cross(get<DERIVATIVE>(x));
1005 
1006  // serac::Functional's boundary integrals multiply the q-function output by
1007  // norm(cross(dX_dxi)) at that quadrature point, but if we impose a shape displacement
1008  // then that weight needs to be corrected. The new weight should be
1009  // norm(cross(dX_dxi + du_dxi + dp_dxi)) where u is displacement and p is shape displacement. This implies:
1010  //
1011  // pressure * normalize(normal_new) * w_new
1012  // = pressure * normalize(normal_new) * (w_new / w_old) * w_old
1013  // = pressure * normalize(normal_new) * (norm(normal_new) / norm(normal_old)) * w_old
1014  // = pressure * (normal_new / norm(normal_new)) * (norm(normal_new) / norm(normal_old)) * w_old
1015  // = pressure * (normal_new / norm(normal_old)) * w_old
1016 
1017  // We always query the pressure function in the undeformed configuration
1018  return pressure_function(get<VALUE>(X), t, params...) * (n / norm(cross(get<DERIVATIVE>(X))));
1019  },
1020  domain);
1021  }
1022 
1024  template <typename PressureType>
1025  void setPressure(PressureType pressure_function, Domain& domain)
1026  {
1027  setPressure(DependsOn<>{}, pressure_function, domain);
1028  }
1029 
1031  virtual std::unique_ptr<mfem_ext::StdFunctionOperator> buildQuasistaticOperator()
1032  {
1033  // the quasistatic case is entirely described by the residual,
1034  // there is no ordinary differential equation
1035  return std::make_unique<mfem_ext::StdFunctionOperator>(
1036  displacement_.space().TrueVSize(),
1037 
1038  // residual function
1039  [this](const mfem::Vector& u, mfem::Vector& r) {
1040  SERAC_MARK_FUNCTION;
1041  const mfem::Vector res =
1042  (*residual_)(time_, shapeDisplacement(), u, acceleration_, *parameters_[parameter_indices].state...);
1043 
1044  // TODO this copy is required as the sundials solvers do not allow move assignments because of their memory
1045  // tracking strategy
1046  // See https://github.com/mfem/mfem/issues/3531
1047  r = res;
1048  r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
1049  },
1050 
1051  // gradient of residual function
1052  [this](const mfem::Vector& u) -> mfem::Operator& {
1053  SERAC_MARK_FUNCTION;
1054  auto [r, drdu] = (*residual_)(time_, shapeDisplacement(), differentiate_wrt(u), acceleration_,
1055  *parameters_[parameter_indices].state...);
1056  J_.reset();
1057  J_ = assemble(drdu);
1058  J_e_.reset();
1059  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
1060  return *J_;
1061  });
1062  }
1063 
1075  std::pair<const mfem::HypreParMatrix&, const mfem::HypreParMatrix&> stiffnessMatrix() const
1076  {
1077  SLIC_ERROR_ROOT_IF(!J_ || !J_e_, "Stiffness matrix has not yet been assembled.");
1078 
1079  return {*J_, *J_e_};
1080  }
1081 
1083  void completeSetup() override
1084  {
1085  if (is_quasistatic_) {
1086  residual_with_bcs_ = buildQuasistaticOperator();
1087  } else {
1088  // the dynamic case is described by a residual function and a second order
1089  // ordinary differential equation. Here, we define the residual function in
1090  // terms of an acceleration.
1091  residual_with_bcs_ = std::make_unique<mfem_ext::StdFunctionOperator>(
1092  displacement_.space().TrueVSize(),
1093 
1094  [this](const mfem::Vector& d2u_dt2, mfem::Vector& r) {
1095  add(1.0, u_, c0_, d2u_dt2, predicted_displacement_);
1096  const mfem::Vector res = (*residual_)(time_, shapeDisplacement(), predicted_displacement_, d2u_dt2,
1097  *parameters_[parameter_indices].state...);
1098 
1099  // TODO this copy is required as the sundials solvers do not allow move assignments because of their memory
1100  // tracking strategy
1101  // See https://github.com/mfem/mfem/issues/3531
1102  r = res;
1103  r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
1104  },
1105 
1106  [this](const mfem::Vector& d2u_dt2) -> mfem::Operator& {
1107  add(1.0, u_, c0_, d2u_dt2, predicted_displacement_);
1108 
1109  // K := dR/du
1110  auto K = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(),
1111  differentiate_wrt(predicted_displacement_), d2u_dt2,
1112  *parameters_[parameter_indices].state...));
1113  std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
1114 
1115  // M := dR/da
1116  auto M = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), predicted_displacement_,
1117  differentiate_wrt(d2u_dt2),
1118  *parameters_[parameter_indices].state...));
1119  std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
1120 
1121  // J = M + c0 * K
1122  J_.reset();
1123  J_.reset(mfem::Add(1.0, *m_mat, c0_, *k_mat));
1124  J_e_.reset();
1125  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
1126 
1127  return *J_;
1128  });
1129  }
1130 
1131 #ifdef SERAC_USE_PETSC
1132  auto* space_dep_pc =
1133  dynamic_cast<serac::mfem_ext::PetscPreconditionerSpaceDependent*>(&nonlin_solver_->preconditioner());
1134  if (space_dep_pc) {
1135  // This call sets the displacement ParFiniteElementSpace used to get the spatial coordinates and to
1136  // generate the near null space for the PCGAMG preconditioner
1137  space_dep_pc->SetFESpace(&displacement_.space());
1138  }
1139 #endif
1140 
1141  nonlin_solver_->setOperator(*residual_with_bcs_);
1142 
1143  if (checkpoint_to_disk_) {
1144  outputStateToDisk();
1145  } else {
1146  checkpoint_states_.clear();
1147  auto state_names = stateNames();
1148  for (const auto& state_name : state_names) {
1149  checkpoint_states_[state_name].push_back(state(state_name));
1150  }
1151  }
1152  }
1153 
1156  {
1157  for (const auto& essential : bcs_.essentials()) {
1158  field.SetSubVector(essential.getTrueDofList(), 0.0);
1159  }
1160  }
1161 
1163  void advanceTimestep(double dt) override
1164  {
1166  SLIC_ERROR_ROOT_IF(!residual_, "completeSetup() must be called prior to advanceTimestep(dt) in SolidMechanics.");
1167 
1168  // If this is the first call, initialize the previous parameter values as the initial values
1169  if (cycle_ == 0) {
1170  for (auto& parameter : parameters_) {
1171  *parameter.previous_state = *parameter.state;
1172  }
1173  }
1174 
1175  if (is_quasistatic_) {
1176  quasiStaticSolve(dt);
1177  } else {
1178  // The current ode interface tracks 2 times, one internally which we have a handle to via time_,
1179  // and one here via the step interface.
1180  // We are ignoring this one, and just using the internal version for now.
1181  // This may need to be revisited when more complex time integrators are required,
1182  // but at the moment, the double times creates a lot of confusion, so
1183  // we short circuit the extra time here by passing a dummy time and ignoring it.
1184  double time_tmp = time_;
1185  dt_ = dt;
1186  ode2_.Step(displacement_, velocity_, time_tmp, dt);
1187  }
1188 
1189  cycle_ += 1;
1190 
1191  if (checkpoint_to_disk_) {
1192  outputStateToDisk();
1193  } else {
1194  for (const auto& state_name : stateNames()) {
1195  checkpoint_states_[state_name].push_back(state(state_name));
1196  }
1197  }
1198 
1199  {
1200  // after finding displacements that satisfy equilibrium,
1201  // compute the residual one more time, this time enabling
1202  // the material state buffers to be updated
1203  residual_->updateQdata(true);
1204 
1205  reactions_ = (*residual_)(time_, shapeDisplacement(), displacement_, acceleration_,
1206  *parameters_[parameter_indices].state...);
1207 
1208  residual_->updateQdata(false);
1209  }
1210 
1211  if (cycle_ > max_cycle_) {
1212  timesteps_.push_back(dt);
1213  max_cycle_ = cycle_;
1214  max_time_ = time_;
1215  }
1216  }
1217 
1232  virtual void setAdjointLoad(std::unordered_map<std::string, const serac::FiniteElementDual&> loads) override
1233  {
1234  SLIC_ERROR_ROOT_IF(loads.size() == 0, "Adjoint load container size must be greater than 0 in the solid mechanics.");
1235 
1236  auto disp_adjoint_load = loads.find("displacement");
1237 
1238  SLIC_ERROR_ROOT_IF(disp_adjoint_load == loads.end(), "Adjoint load for \"displacement\" not found.");
1239 
1240  if (disp_adjoint_load != loads.end()) {
1241  displacement_adjoint_load_ = disp_adjoint_load->second;
1242  // Add the sign correction to move the term to the RHS
1243  displacement_adjoint_load_ *= -1.0;
1244  }
1245 
1246  auto velo_adjoint_load = loads.find("velocity");
1247 
1248  if (velo_adjoint_load != loads.end()) {
1249  velocity_adjoint_load_ = velo_adjoint_load->second;
1250  // Add the sign correction to move the term to the RHS
1251  velocity_adjoint_load_ *= -1.0;
1252  }
1253 
1254  auto accel_adjoint_load = loads.find("acceleration");
1255 
1256  if (accel_adjoint_load != loads.end()) {
1257  acceleration_adjoint_load_ = accel_adjoint_load->second;
1258  // Add the sign correction to move the term to the RHS
1259  acceleration_adjoint_load_ *= -1.0;
1260  }
1261  }
1262 
1263  virtual void setDualAdjointBcs(std::unordered_map<std::string, const serac::FiniteElementState&> bcs) override
1264  {
1265  SLIC_ERROR_ROOT_IF(bcs.size() == 0, "Adjoint load container size must be greater than 0 in the solid mechanics.");
1266 
1267  auto reaction_adjoint_load = bcs.find("reactions");
1268 
1269  SLIC_ERROR_ROOT_IF(reaction_adjoint_load == bcs.end(), "Adjoint load for \"reaction\" not found.");
1270 
1271  if (reaction_adjoint_load != bcs.end()) {
1272  reactions_adjoint_bcs_ = reaction_adjoint_load->second;
1273  }
1274  }
1275 
1277  void reverseAdjointTimestep() override
1278  {
1280  SLIC_ERROR_ROOT_IF(cycle_ <= min_cycle_,
1281  "Maximum number of adjoint timesteps exceeded! The number of adjoint timesteps must equal the "
1282  "number of forward timesteps");
1283 
1284  cycle_--; // cycle is now at n \in [0,N-1]
1285 
1286  double dt_np1_to_np2 = getCheckpointedTimestep(cycle_ + 1);
1287  const double dt_n_to_np1 = getCheckpointedTimestep(cycle_);
1288 
1289  auto end_step_solution = getCheckpointedStates(cycle_ + 1);
1290 
1291  displacement_ = end_step_solution.at("displacement");
1292 
1293  if (is_quasistatic_) {
1294  quasiStaticAdjointSolve(dt_n_to_np1);
1295  } else {
1296  SLIC_ERROR_ROOT_IF(ode2_.GetTimestepper() != TimestepMethod::Newmark,
1297  "Only Newmark implemented for transient adjoint solid mechanics.");
1298 
1299  // Load the end of step velo, accel from the previous cycle
1300 
1301  velocity_ = end_step_solution.at("velocity");
1302  acceleration_ = end_step_solution.at("acceleration");
1303 
1304  // K := dR/du
1305  auto K = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), differentiate_wrt(displacement_),
1306  acceleration_, *parameters_[parameter_indices].state...));
1307  std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
1308 
1309  // M := dR/da
1310  auto M = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), displacement_,
1311  differentiate_wrt(acceleration_),
1312  *parameters_[parameter_indices].state...));
1313  std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
1314 
1315  auto& lin_solver = nonlin_solver_->linearSolver();
1316  solid_mechanics::detail::adjoint_integrate(
1317  dt_n_to_np1, dt_np1_to_np2, m_mat.get(), k_mat.get(), displacement_adjoint_load_, velocity_adjoint_load_,
1318  acceleration_adjoint_load_, adjoint_displacement_, implicit_sensitivity_displacement_start_of_step_,
1319  implicit_sensitivity_velocity_start_of_step_, reactions_adjoint_bcs_, bcs_, lin_solver);
1320  }
1321 
1322  time_end_step_ = time_;
1323  time_ -= dt_n_to_np1;
1324  }
1325 
1327  FiniteElementDual computeTimestepSensitivity(size_t parameter_field) override
1328  {
1329  SLIC_ASSERT_MSG(parameter_field < sizeof...(parameter_indices),
1330  axom::fmt::format("Invalid parameter index '{}' requested for sensitivity."));
1331 
1332  auto drdparam = serac::get<DERIVATIVE>(d_residual_d_[parameter_field](time_end_step_));
1333  auto drdparam_mat = assemble(drdparam);
1334 
1335  drdparam_mat->MultTranspose(adjoint_displacement_, *parameters_[parameter_field].sensitivity);
1336 
1337  return *parameters_[parameter_field].sensitivity;
1338  }
1339 
1342  {
1343  auto drdshape =
1344  serac::get<DERIVATIVE>((*residual_)(time_end_step_, differentiate_wrt(shapeDisplacement()), displacement_,
1345  acceleration_, *parameters_[parameter_indices].state...));
1346 
1347  auto drdshape_mat = assemble(drdshape);
1348 
1349  drdshape_mat->MultTranspose(adjoint_displacement_, shape_displacement_dual_);
1350 
1351  return shapeDisplacementSensitivity();
1352  }
1353 
1355  const std::unordered_map<std::string, const serac::FiniteElementDual&> computeInitialConditionSensitivity()
1356  const override
1357  {
1358  return {{"displacement", implicit_sensitivity_displacement_start_of_step_},
1359  {"velocity", implicit_sensitivity_velocity_start_of_step_}};
1360  }
1361 
1367  const serac::FiniteElementState& displacement() const { return displacement_; };
1368 
1374  const serac::FiniteElementState& velocity() const { return velocity_; };
1375 
1381  const serac::FiniteElementState& acceleration() const { return acceleration_; };
1382 
1384  const serac::FiniteElementDual& reactions() const { return reactions_; };
1385 
1386  protected:
1389 
1392 
1396 
1399 
1402 
1405 
1406  // In the case of transient dynamics, this is more like an adjoint_acceleration
1409 
1412 
1415 
1418 
1421 
1424 
1427 
1430 
1432  std::unique_ptr<ShapeAwareFunctional<shape_trial, test(trial, trial, parameter_space...)>> residual_;
1433 
1435  std::unique_ptr<mfem_ext::StdFunctionOperator> residual_with_bcs_;
1436 
1438  std::unique_ptr<EquationSolver> nonlin_solver_;
1439 
1446 
1448  std::unique_ptr<mfem::HypreParMatrix> J_;
1449 
1452  std::unique_ptr<mfem::HypreParMatrix> J_e_;
1453 
1456 
1458  mfem::Vector du_;
1459 
1461  mfem::Vector u_;
1462 
1464  mfem::Vector v_;
1465 
1467  double c0_;
1468 
1470  double c1_;
1471 
1475 
1478 
1480  std::shared_ptr<mfem::VectorCoefficient> disp_bdr_coef_;
1481 
1483  std::shared_ptr<mfem::Coefficient> component_disp_bdr_coef_;
1484 
1488  std::array<std::function<decltype((*residual_)(DifferentiateWRT<1>{}, 0.0, shape_displacement_, displacement_,
1489  acceleration_, *parameters_[parameter_indices].state...))(double)>,
1490  sizeof...(parameter_indices)>
1491  d_residual_d_ = {[&](double _t) {
1492  return (*residual_)(DifferentiateWRT<NUM_STATE_VARS + 1 + parameter_indices>{}, _t, shape_displacement_,
1493  displacement_, acceleration_, *parameters_[parameter_indices].state...);
1494  }...};
1495 
1497  virtual void quasiStaticSolve(double dt)
1498  {
1499  // warm start must be called prior to the time update so that the previous Jacobians can be used consistently
1500  // throughout.
1501  warmStartDisplacement(dt);
1502  time_ += dt;
1503 
1504  // this method is essentially equivalent to the 1-liner
1505  // u += dot(inv(J), dot(J_elim[:, dofs], (U(t + dt) - u)[dofs]));
1506  nonlin_solver_->solve(displacement_);
1507  }
1508 
1510  virtual void quasiStaticAdjointSolve(double /*dt*/)
1511  {
1512  auto [_, drdu] = (*residual_)(time_, shapeDisplacement(), differentiate_wrt(displacement_), acceleration_,
1513  *parameters_[parameter_indices].state...);
1514  J_.reset();
1515  J_ = assemble(drdu);
1516 
1517  auto J_T = std::unique_ptr<mfem::HypreParMatrix>(J_->Transpose());
1518 
1519  J_e_.reset();
1520  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_T);
1521 
1522  auto& constrained_dofs = bcs_.allEssentialTrueDofs();
1523 
1524  mfem::EliminateBC(*J_T, *J_e_, constrained_dofs, reactions_adjoint_bcs_, displacement_adjoint_load_);
1525  for (int i = 0; i < constrained_dofs.Size(); i++) {
1526  int j = constrained_dofs[i];
1527  displacement_adjoint_load_[j] = reactions_adjoint_bcs_[j];
1528  }
1529 
1530  auto& lin_solver = nonlin_solver_->linearSolver();
1531  lin_solver.SetOperator(*J_T);
1532  lin_solver.Mult(displacement_adjoint_load_, adjoint_displacement_);
1533 
1534  // Reset the equation solver to use the full nonlinear residual operator. MRT, is this needed?
1535  nonlin_solver_->setOperator(*residual_with_bcs_);
1536  }
1537 
1547  mfem::Array<int> calculateConstrainedDofs(std::function<bool(const mfem::Vector&)> is_node_constrained,
1548  std::optional<int> component = {}) const
1549  {
1550  // Get the nodal positions for the displacement vector in grid function form
1551  mfem::ParGridFunction nodal_positions(
1552  const_cast<mfem::ParFiniteElementSpace*>(&displacement_.space())); // MRT mfem const correctness issue
1553  mfemParMesh().GetNodes(nodal_positions);
1554 
1555  const int num_nodes = nodal_positions.Size() / dim;
1556  mfem::Array<int> constrained_dofs;
1557 
1558  for (int i = 0; i < num_nodes; i++) {
1559  // Determine if this "local" node (L-vector node) is in the local true vector. I.e. ensure this node is not a
1560  // shared node owned by another processor
1561  int idof = mfem::Ordering::Map<serac::ordering>(nodal_positions.FESpace()->GetNDofs(),
1562  nodal_positions.FESpace()->GetVDim(), i, 0);
1563  if (nodal_positions.ParFESpace()->GetLocalTDofNumber(idof) >= 0) {
1564  mfem::Vector node_coords(dim);
1565  mfem::Array<int> node_dofs;
1566  for (int d = 0; d < dim; d++) {
1567  // Get the local dof number for the prescribed component
1568  int local_vector_dof = mfem::Ordering::Map<serac::ordering>(nodal_positions.FESpace()->GetNDofs(),
1569  nodal_positions.FESpace()->GetVDim(), i, d);
1570 
1571  // Save the spatial position for this coordinate dof
1572  node_coords(d) = nodal_positions(local_vector_dof);
1573 
1574  // Check if this component of the displacement vector is constrained
1575  bool is_active_component = true;
1576  if (component) {
1577  if (*component != d) {
1578  is_active_component = false;
1579  }
1580  }
1581 
1582  if (is_active_component) {
1583  // Add the true dof for this component to the related dof list
1584  node_dofs.Append(nodal_positions.ParFESpace()->GetLocalTDofNumber(local_vector_dof));
1585  }
1586  }
1587 
1588  // Do the user-defined spatial query to determine if these dofs should be constrained
1589  if (is_node_constrained(node_coords)) {
1590  constrained_dofs.Append(node_dofs);
1591  }
1592 
1593  // Reset the temporary container for the dofs associated with a particular node
1594  node_dofs.DeleteAll();
1595  }
1596  }
1597  return constrained_dofs;
1598  }
1599 
1633  void warmStartDisplacement(double dt)
1634  {
1636 
1637  du_ = 0.0;
1638  for (auto& bc : bcs_.essentials()) {
1639  // apply the future boundary conditions, but use the most recent Jacobians stiffness.
1640  bc.setDofs(du_, time_ + dt);
1641  }
1642 
1643  auto& constrained_dofs = bcs_.allEssentialTrueDofs();
1644  for (int i = 0; i < constrained_dofs.Size(); i++) {
1645  int j = constrained_dofs[i];
1646  du_[j] -= displacement_(j);
1647  }
1648 
1649  if (use_warm_start_) {
1650  // Update external forcing
1651  auto r = (*residual_)(time_ + dt, shapeDisplacement(), displacement_, acceleration_,
1652  *parameters_[parameter_indices].state...);
1653 
1654  // use the most recently evaluated Jacobian
1655  auto [_, drdu] = (*residual_)(time_, shapeDisplacement(), differentiate_wrt(displacement_), acceleration_,
1656  *parameters_[parameter_indices].previous_state...);
1657  J_.reset();
1658  J_ = assemble(drdu);
1659  J_e_.reset();
1660  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
1661 
1662  r *= -1.0;
1663 
1664  mfem::EliminateBC(*J_, *J_e_, constrained_dofs, du_, r);
1665  for (int i = 0; i < constrained_dofs.Size(); i++) {
1666  int j = constrained_dofs[i];
1667  r[j] = du_[j];
1668  }
1669 
1670  auto& lin_solver = nonlin_solver_->linearSolver();
1671 
1672  lin_solver.SetOperator(*J_);
1673 
1674  lin_solver.Mult(r, du_);
1675  }
1676 
1677  displacement_ += du_;
1678  dt_ = dt;
1679  }
1680 };
1681 
1682 } // namespace serac
This file contains the interface used for initializing/terminating any hardware accelerator-related f...
#define SERAC_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
Definition: accelerator.hpp:38
The base interface class for a generic PDE solver.
This file contains the declaration of the boundary condition manager class.
This is the abstract base class for a generic forward solver.
void initializeBasePhysicsStates(int cycle, double time)
Protected, non-virtual method to reset physics states to zero. This does not reset design parameters ...
A set to flag components of a vector field.
Definition: components.hpp:29
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
Class for encapsulating the critical MFEM components of a primal finite element field.
Class for encapsulating the data associated with a vector derived from a MFEM finite element space....
mfem::ParFiniteElementSpace & space()
Returns a non-owning reference to the internal FESpace.
const serac::FiniteElementDual & reactions() const
getter for nodal forces (before zeroing-out essential dofs)
void addCustomBoundaryIntegral(DependsOn< active_parameters... >, callable qfunction, const std::optional< Domain > &optional_domain=std::nullopt)
register a custom boundary integral calculation as part of the residual
void setMaterial(DependsOn< active_parameters... >, const MaterialType &material, Domain &domain, qdata_type< StateType > qdata=EmptyQData)
Set the material stress response and mass properties for the physics module.
void setMaterial(const MaterialType &material, Domain &domain, std::shared_ptr< QuadratureData< StateType >> qdata=EmptyQData)
This is an overloaded member function, provided for convenience. It differs from the above function o...
mfem_ext::SecondOrderODE ode2_
the ordinary differential equation that describes how to solve for the second time derivative of disp...
const FiniteElementState & dualAdjoint(const std::string &dual_name) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< mfem::HypreParMatrix > J_
Assembled sparse matrix for the Jacobian df/du (11 block if using Lagrange multiplier contact)
mfem::Vector u_
used to communicate the ODE solver's predicted displacement to the residual operator
void addBodyForce(BodyForceType body_force, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void zeroEssentials(FiniteElementVector &field) const
Set field to zero wherever their are essential boundary conditions applies.
FiniteElementDual implicit_sensitivity_velocity_start_of_step_
The total/implicit sensitivity of the qoi with respect to the start of the previous timestep's veloci...
std::unique_ptr< mfem_ext::StdFunctionOperator > residual_with_bcs_
mfem::Operator that calculates the residual after applying essential boundary conditions
void setTraction(DependsOn< active_parameters... >, TractionType traction_function, Domain &domain)
Set the traction boundary condition.
const FiniteElementDual & computeTimestepShapeSensitivity() override
This is an overloaded member function, provided for convenience. It differs from the above function o...
mfem::Array< int > calculateConstrainedDofs(std::function< bool(const mfem::Vector &)> is_node_constrained, std::optional< int > component={}) const
Calculate a list of constrained dofs in the true displacement vector from a function that returns tru...
void resetStates(int cycle=0, double time=0.0) override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void setDisplacementBCsByDofList(AppliedDisplacementFunction applied_displacement, const mfem::Array< int > true_dofs)
Set the displacement essential boundary conditions on a set of true degrees of freedom.
void setPressure(DependsOn< active_parameters... >, PressureType pressure_function, Domain &domain)
Apply a pressure-type follower load.
virtual void setAdjointLoad(std::unordered_map< std::string, const serac::FiniteElementDual & > loads) override
Set the loads for the adjoint reverse timestep solve.
virtual void setDualAdjointBcs(std::unordered_map< std::string, const serac::FiniteElementState & > bcs) override
Set the dual loads (dirichlet values) for the adjoint reverse timestep solve This must be called afte...
void setTraction(TractionType traction_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void completeSetup() override
This is an overloaded member function, provided for convenience. It differs from the above function o...
FiniteElementDual velocity_adjoint_load_
The adjoint load (RHS) for the velocity adjoint system solve (downstream -dQOI/d velocity)
std::shared_ptr< mfem::Coefficient > component_disp_bdr_coef_
Coefficient containing the essential boundary values.
mfem::Vector predicted_displacement_
an intermediate variable used to store the predicted end-step displacement
double time_end_step_
End of step time used in reverse mode so that the time can be decremented on reverse steps.
void setRateDependentMaterial(const RateDependentMaterialType &material, Domain &domain, std::shared_ptr< QuadratureData< StateType >> qdata=EmptyQData)
This is an overloaded member function, provided for convenience. It differs from the above function o...
mfem::Vector v_
used to communicate the ODE solver's predicted velocity to the residual operator
void setVelocity(AppliedVelocityFunction applied_velocity)
Set the underlying finite element state to a prescribed velocity field.
void initializeSolidMechanicsStates()
Non virtual method to reset thermal states to zero. This does not reset design parameters or shape.
virtual std::unique_ptr< mfem_ext::StdFunctionOperator > buildQuasistaticOperator()
Build the quasi-static operator corresponding to the total Lagrangian formulation.
void setRateDependentMaterial(DependsOn< active_parameters... >, const RateDependentMaterialType &material, Domain &domain, qdata_type< StateType > qdata=EmptyQData)
void setDisplacement(const FiniteElementState &temp)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void setState(const std::string &state_name, const FiniteElementState &state) override
Set the primal solution field (displacement, velocity) for the underlying solid mechanics solver.
void setDisplacement(AppliedDisplacementFunction applied_displacement)
Set the underlying finite element state to a prescribed displacement field.
void setPressure(PressureType pressure_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
double c1_
coefficient used to calculate predicted velocity: dudt_p := dudt + c1 * d2u_dt2
std::shared_ptr< mfem::VectorCoefficient > disp_bdr_coef_
Coefficient containing the essential boundary values.
const FiniteElementState & state(const std::string &state_name) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::vector< std::string > adjointNames() const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void advanceTimestep(double dt) override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::pair< const mfem::HypreParMatrix &, const mfem::HypreParMatrix & > stiffnessMatrix() const
Return the assembled stiffness matrix.
FiniteElementDual implicit_sensitivity_displacement_start_of_step_
The total/implicit sensitivity of the qoi with respect to the start of the previous timestep's displa...
SolidMechanics(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< serac::Mesh > serac_mesh, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false, bool use_warm_start=true)
Construct a new SolidMechanics object.
void reverseAdjointTimestep() override
This is an overloaded member function, provided for convenience. It differs from the above function o...
const FiniteElementDual & dual(const std::string &dual_name) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::vector< std::string > dualNames() const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
FiniteElementDual loadCheckpointedDual(const std::string &dual_name, int cycle) override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addBodyForce(DependsOn< active_parameters... >, BodyForceType body_force, Domain &domain)
Set the body forcefunction.
FiniteElementDual displacement_adjoint_load_
The adjoint load (RHS) for the displacement adjoint system solve (downstream -dQOI/d displacement)
std::unique_ptr< ShapeAwareFunctional< shape_trial, test(trial, trial, parameter_space...)> > residual_
serac::Functional that is used to calculate the residual and its derivatives
void setVelocity(const FiniteElementState &temp)
This is an overloaded member function, provided for convenience. It differs from the above function o...
const std::unordered_map< std::string, const serac::FiniteElementDual & > computeInitialConditionSensitivity() const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void addCustomDomainIntegral(DependsOn< active_parameters... >, callable qfunction, Domain &domain, qdata_type< StateType > qdata=NoQData)
register a custom domain integral calculation as part of the residual
void setFixedBCs(Domain &domain, Components components=Component::ALL)
Shortcut to set selected components of displacements to zero for all time.
std::vector< std::string > stateNames() const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
const FiniteElementState & adjoint(const std::string &state_name) const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
void warmStartDisplacement(double dt)
Sets the Dirichlet BCs for the current time and computes an initial guess for parameters and displace...
FiniteElementDual computeTimestepSensitivity(size_t parameter_field) override
This is an overloaded member function, provided for convenience. It differs from the above function o...
FiniteElementDual acceleration_adjoint_load_
The adjoint load (RHS) for the adjoint system solve (downstream -dQOI/d acceleration)
std::unique_ptr< EquationSolver > nonlin_solver_
the specific methods and tolerances specified to solve the nonlinear residual equations
std::shared_ptr< QuadratureData< T > > qdata_type
a container holding quadrature point data of the specified type
qdata_type< T > createQuadratureDataBuffer(T initial_state, const std::optional< Domain > &optional_domain=std::nullopt)
Create a shared ptr to a quadrature data buffer for the given material type.
SolidMechanics(std::unique_ptr< serac::EquationSolver > solver, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< serac::Mesh > serac_mesh, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false, bool use_warm_start=true)
Construct a new SolidMechanics object.
void setDisplacementBCs(AppliedDisplacementFunction applied_displacement, Domain &domain, Components components=Component::ALL)
Set essential displacement boundary conditions on selected components.
static std::shared_ptr< QuadratureData< T > > newQuadratureDataBuffer(const std::string &mesh_tag, const Domain &domain, int order, int dim, T initial_state)
Create a shared ptr to a quadrature data buffer for the given material type.
static FiniteElementState newState(FunctionSpace space, const std::string &state_name, const std::string &mesh_tag)
Factory method for creating a new FEState object.
static FiniteElementDual newDual(FunctionSpace space, const std::string &dual_name, const std::string &mesh_tag)
Factory method for creating a new FEDual object.
SecondOrderODE is a class wrapping mfem::SecondOrderTimeDependentOperator so that the user can use st...
Definition: odes.hpp:36
A file defining some enums and structs that are used by the different physics modules.
many of the functions in this file amount to extracting element indices from an mesh_t like
This file contains the declaration of an equation solver wrapper.
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 ...
This file contains the declaration of structure that manages vectors derived from an MFEM finite elem...
Implementation of the quadrature-function-based functional enabling rapid development of FEM formulat...
Serac mesh class which assists in constructing the appropriate parallel mfem meshes and registering a...
constexpr SERAC_HOST_DEVICE void for_constexpr(const lambda &f)
multidimensional loop tool that evaluates the lambda body inside the innermost loop.
const TimesteppingOptions default_timestepping_options
default implicit dynamic timestepping options for solid mechanics
const LinearSolverOptions direct_linear_options
the default direct solver option for solving the linear stiffness equations
const NonlinearSolverOptions default_nonlinear_options
default iteration limits, tolerances and verbosity for solving the systems of nonlinear equations tha...
const TimesteppingOptions default_quasistatic_options
default quasistatic timestepping options for solid mechanics
const LinearSolverOptions default_linear_options
default method and tolerances for solving the systems of linear equations that show up in implicit so...
Accelerator functionality.
Definition: serac.cpp:36
std::shared_ptr< QuadratureData< Empty > > EmptyQData
a single instance of a QuadratureData container of Emptys, since they are all interchangeable
std::shared_ptr< QuadratureData< Nothing > > NoQData
a single instance of a QuadratureData container of Nothings, since they are all interchangeable
auto differentiate_wrt(const mfem::Vector &v)
this function is intended to only be used in combination with serac::Functional::operator(),...
Domain EntireBoundary(const mesh_t &mesh)
constructs a domain from all the boundary elements in a mesh
Definition: domain.cpp:609
tuple(T...) -> tuple< T... >
Class template argument deduction rule for tuples.
constexpr SERAC_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
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:961
SERAC_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:1117
Domain EntireDomain(const mesh_t &mesh)
constructs a domain from all the elements in a mesh
Definition: domain.cpp:594
Wrapper classes for using MFEM's ODE solver objects.
Various helper functions and macros for profiling using Caliper.
#define SERAC_MARK_FUNCTION
Definition: profiling.hpp:90
This file contains the declaration of the structures that manage quadrature point data.
Wrapper of serac::Functional for evaluating integrals and derivatives of quantities with shape displa...
The material and load types for the solid functional physics module.
An object containing all input file options for the solver for total Lagrangian finite deformation so...
This file contains enumerations and record types for physics solver configuration.
This file contains the declaration of the StateManager class.
Classes for defining an mfem::Operator with std::functions.
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
mfem::Array< int > dof_list(const fes_t *fes) const
get mfem degree of freedom list for a given FiniteElementSpace
Definition: domain.cpp:466
see Nothing for a complete description of this class and when to use it
H1 elements of order p.
Parameters for an iterative linear solution scheme.
LinearSolver linear_solver
Linear solver selection.
Nonlinear solution scheme parameters.
NonlinearSolver nonlin_solver
Nonlinear solver selection.
these classes are a little confusing. These two special types represent the similar (but different) c...
a struct that is used in the physics modules to clarify which template arguments are user-controlled ...
Definition: common.hpp:21
A class for storing and access user-defined types at quadrature points.
auto SERAC_HOST_DEVICE operator()(T t, Position position, Displacement, Acceleration, Params... params) const
Body force call.
auto SERAC_HOST_DEVICE operator()(double, X, State &state, Displacement displacement, Acceleration acceleration, Params... params) const
Material stress response call.
auto SERAC_HOST_DEVICE operator()(double, X, State &state, Displacement displacement, Acceleration acceleration, Params... params) const
Material stress response call.
A timestep and boundary condition enforcement method for a dynamic solver.
TimestepMethod timestepper
The timestepping method to be applied.
DirichletEnforcementMethod enforcement_method
The essential boundary enforcement method to use.
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
Implementation of the tensor class used by Functional.
Implements a std::tuple-like object that works in CUDA kernels.