Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
solid_mechanics_contact.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 
13 #pragma once
14 
15 #include <charconv>
16 #include <memory>
17 #include <optional>
18 #include <string_view>
19 #include <system_error>
20 #include <unordered_map>
21 #include <vector>
22 
25 
26 namespace smith {
27 
28 template <int order, int dim, typename parameters = Parameters<>,
29  typename parameter_indices = std::make_integer_sequence<int, parameters::n>>
31 
41 template <int order, int dim, typename... parameter_space, int... parameter_indices>
42 class SolidMechanicsContact<order, dim, Parameters<parameter_space...>,
43  std::integer_sequence<int, parameter_indices...>>
44  : public SolidMechanics<order, dim, Parameters<parameter_space...>,
45  std::integer_sequence<int, parameter_indices...>> {
46  using SolidMechanicsBase =
47  SolidMechanics<order, dim, Parameters<parameter_space...>, std::integer_sequence<int, parameter_indices...>>;
48 
49  public:
66  const smith::TimesteppingOptions timestepping_opts, const std::string& physics_name,
67  std::shared_ptr<smith::Mesh> smith_mesh, std::vector<std::string> parameter_names = {},
68  int cycle = 0, double time = 0.0, bool checkpoint_to_disk = false, bool use_warm_start = true)
69  : SolidMechanicsContact(std::make_unique<EquationSolver>(nonlinear_opts, lin_opts, smith_mesh->getComm()),
70  timestepping_opts, physics_name, smith_mesh, parameter_names, cycle, time,
71  checkpoint_to_disk, use_warm_start)
72  {
73  }
74 
89  SolidMechanicsContact(std::unique_ptr<smith::EquationSolver> solver,
90  const smith::TimesteppingOptions timestepping_opts, const std::string& physics_name,
91  std::shared_ptr<smith::Mesh> smith_mesh, std::vector<std::string> parameter_names = {},
92  int cycle = 0, double time = 0.0, bool checkpoint_to_disk = false, bool use_warm_start = true)
93  : SolidMechanicsBase(std::move(solver), timestepping_opts, physics_name, smith_mesh, parameter_names, cycle, time,
94  checkpoint_to_disk, use_warm_start),
95  contact_(BasePhysics::mfemParMesh())
96  {
97  }
98 
108  SolidMechanicsContact(const SolidMechanicsInputOptions& input_options, const std::string& physics_name,
109  std::shared_ptr<smith::Mesh> smith_mesh, int cycle = 0, double time = 0.0)
110  : SolidMechanicsBase(input_options, physics_name, smith_mesh, cycle, time), contact_(BasePhysics::mfemParMesh())
111  {
112  }
113 
115  void resetStates(int cycle = 0, double time = 0.0) override
116  {
117  SolidMechanicsBase::resetStates(cycle, time);
118  for (auto& [_, force] : contact_interaction_forces_) {
119  if (force) {
120  *force = 0.0;
121  }
122  }
123  for (auto& [_, seed] : contact_interaction_force_adjoint_bcs_) {
124  if (seed) {
125  *seed = 0.0;
126  }
127  }
128  contact_.reset();
129  double dt = 0.0;
130  mfem::Vector p(contact_.numPressureDofs());
131  p = 0.0;
132  contact_.updateForcesAndJacobian(cycle, time, dt, BasePhysics::shapeDisplacement(), displacement_, p);
133  updateContactForceOutputs();
134  }
135 
137  std::unique_ptr<mfem_ext::StdFunctionOperator> buildQuasistaticOperator() override
138  {
139  auto residual_fn = [this](const mfem::Vector& u, mfem::Vector& r) {
140  const mfem::Vector u_blk(const_cast<mfem::Vector&>(u), 0, displacement_.Size());
141  const mfem::Vector res = (*residual_)(time_, BasePhysics::shapeDisplacement(), u_blk, acceleration_,
142  *parameters_[parameter_indices].state...);
143 
144  // NOTE this copy is required as the sundials solvers do not allow move assignments because of their memory
145  // tracking strategy
146  // See https://github.com/mfem/mfem/issues/3531
147  mfem::Vector r_blk(r, 0, displacement_.Size());
148  r_blk = res;
149 
150  contact_.residualFunction(BasePhysics::shapeDisplacement(), u, r);
151  r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
152  };
153  // This if-block below breaks up building the Jacobian operator depending if there is Lagrange multiplier
154  // enforcement or not
155  if (contact_.haveLagrangeMultipliers()) {
156  // The quasistatic operator has blocks if any of the contact interactions are enforced using Lagrange multipliers.
157  // Jacobian operator is an mfem::BlockOperator
158  J_offsets_ = mfem::Array<int>({0, displacement_.Size(), displacement_.Size() + contact_.numPressureDofs()});
159  return std::make_unique<mfem_ext::StdFunctionOperator>(
160  displacement_.space().TrueVSize() + contact_.numPressureDofs(), residual_fn,
161  // gradient of residual function
162  [this](const mfem::Vector& u) -> mfem::Operator& {
163  const mfem::Vector u_blk(const_cast<mfem::Vector&>(u), 0, displacement_.Size());
164  auto [r, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(u_blk),
165  acceleration_, *parameters_[parameter_indices].state...);
166 
167  // create block operator holding jacobian contributions
168  J_constraint_ = contact_.jacobianFunction(assemble(drdu));
169 
170  // take ownership of blocks
171  J_constraint_->owns_blocks = false;
172  J_ = std::unique_ptr<mfem::HypreParMatrix>(
173  static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 0)));
174  J_12_ = std::unique_ptr<mfem::HypreParMatrix>(
175  static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 1)));
176  J_21_ = std::unique_ptr<mfem::HypreParMatrix>(
177  static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 0)));
178  J_22_ = std::unique_ptr<mfem::HypreParMatrix>(
179  static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 1)));
180 
181  // eliminate bcs and compute eliminated blocks
182  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
183  J_e_21_ = std::unique_ptr<mfem::HypreParMatrix>(J_21_->EliminateCols(bcs_.allEssentialTrueDofs()));
184  J_12_->EliminateRows(bcs_.allEssentialTrueDofs());
185 
186  // create block operator for constraints
187  J_constraint_e_ = std::make_unique<mfem::BlockOperator>(J_offsets_);
188  J_constraint_e_->SetBlock(0, 0, J_e_.get());
189  J_constraint_e_->SetBlock(1, 0, J_e_21_.get());
190 
191  J_operator_ = J_constraint_.get();
192  return *J_constraint_;
193  });
194  } else {
195  // If all of the contact interactions are penalty, then there will be no blocks. Jacobian operator is a single
196  // mfem::HypreParMatrix
197  return std::make_unique<mfem_ext::StdFunctionOperator>(
198  displacement_.space().TrueVSize(), residual_fn, [this](const mfem::Vector& u) -> mfem::Operator& {
199  auto [r, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(u), acceleration_,
200  *parameters_[parameter_indices].state...);
201 
202  // get 11-block holding Jacobian contributions
203  auto block_J = contact_.jacobianFunction(assemble(drdu));
204  block_J->owns_blocks = false;
205  J_ = std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&block_J->GetBlock(0, 0)));
206 
207  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
208 
209  J_operator_ = J_.get();
210  return *J_;
211  });
212  }
213  }
214 
223  void addContactInteraction(int interaction_id, const std::set<int>& bdry_attr_surf1,
224  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts)
225  {
226  SLIC_ERROR_ROOT_IF(!is_quasistatic_, "Contact can only be applied to quasistatic problems.");
227  SLIC_ERROR_ROOT_IF(order > 1, "Contact can only be applied to linear (order = 1) meshes.");
228 
229  // Disallow duplicates for this physics module. Tribol ids are global, so reusing the same interaction_id is an
230  // error even if ContactData rejects it later.
231  SLIC_ERROR_ROOT_IF(
232  contact_interaction_forces_.find(interaction_id) != contact_interaction_forces_.end(),
233  std::format("Contact interaction id {} already exists for physics module '{}'.", interaction_id, name_));
234 
235  // Register with Tribol via ContactData (includes global uniqueness check).
236  contact_.addContactInteraction(interaction_id, bdry_attr_surf1, bdry_attr_surf2, contact_opts);
237 
238  const auto interaction_force_name = detail::addPrefix(name_, std::format("contact_force_{}", interaction_id));
239  auto interaction_force_dual =
240  std::make_unique<FiniteElementDual>(StateManager::newDual(displacement_.space(), interaction_force_name));
241  *interaction_force_dual = 0.0;
242  duals_.push_back(interaction_force_dual.get());
243  contact_interaction_forces_.emplace(interaction_id, std::move(interaction_force_dual));
244 
245  const auto force_adjoint_bcs_name =
246  detail::addPrefix(name_, std::format("contact_force_adjoint_bcs_{}", interaction_id));
247  auto force_adjoint_bcs = std::make_unique<FiniteElementState>(displacement_.space(), force_adjoint_bcs_name);
248  *force_adjoint_bcs = 0.0;
249  this->dual_adjoints_.push_back(force_adjoint_bcs.get());
250  contact_interaction_force_adjoint_bcs_.emplace(interaction_id, std::move(force_adjoint_bcs));
251  }
252 
254  std::vector<std::string> dualNames() const override
255  {
256  auto dual_names = SolidMechanicsBase::dualNames();
257 #ifdef SMITH_USE_TRIBOL
258  for (const auto& interaction : contact_.getContactInteractions()) {
259  dual_names.push_back(std::format("contact_force_{}", interaction.getInteractionId()));
260  }
261 #endif
262  return dual_names;
263  }
264 
266  const FiniteElementDual& dual(const std::string& dual_name) const override
267  {
268  const auto interaction_id = parseContactInteractionForceId(dual_name);
269  if (interaction_id.has_value()) {
270  auto it = contact_interaction_forces_.find(*interaction_id);
271  if (it != contact_interaction_forces_.end()) {
272  return *it->second;
273  }
274  }
275 
276  return SolidMechanicsBase::dual(dual_name);
277  }
278 
280  const FiniteElementState& dualAdjoint(const std::string& dual_name) const override
281  {
282  const auto interaction_id = parseContactInteractionForceId(dual_name);
283  if (interaction_id.has_value()) {
284  auto it = contact_interaction_force_adjoint_bcs_.find(*interaction_id);
285  if (it != contact_interaction_force_adjoint_bcs_.end()) {
286  return *it->second;
287  }
288  }
289 
290  return SolidMechanicsBase::dualAdjoint(dual_name);
291  }
292 
294  FiniteElementDual loadCheckpointedDual(const std::string& dual_name, int cycle) override
295  {
296  if (parseContactInteractionForceId(dual_name).has_value()) {
297  SLIC_ERROR_ROOT_IF(contact_.haveLagrangeMultipliers(),
298  "Checkpointed retrieval of contact forces is not supported for Lagrange multiplier contact.");
299 
300  const FiniteElementState checkpointed_displacement = this->loadCheckpointedState("displacement", cycle);
301  double dt = this->getCheckpointedTimestep(cycle);
302  contact_.updateForcesAndJacobian(cycle, time_, dt, BasePhysics::shapeDisplacement(), checkpointed_displacement);
303  updateContactForceOutputs();
304 
305  const auto interaction_id = parseContactInteractionForceId(dual_name);
306  SLIC_ERROR_ROOT_IF(
307  !interaction_id.has_value(),
308  std::format("Requested checkpointed dual '{}' does not exist in physics module '{}'.", dual_name, name_));
309  auto it = contact_interaction_forces_.find(*interaction_id);
310  SLIC_ERROR_ROOT_IF(
311  it == contact_interaction_forces_.end(),
312  std::format("Requested checkpointed dual '{}' does not exist in physics module '{}'.", dual_name, name_));
313  return *it->second;
314  }
315 
316  return SolidMechanicsBase::loadCheckpointedDual(dual_name, cycle);
317  }
318 
323  {
324  // compute contact dof --> displacement dof prolongation operator
325  // if not previously computed
326  if (!contact_dof_prolongation_) {
327  contact_dof_prolongation_ = contact_.contactSubspaceTransferOperator();
328  }
329  }
330 
336  void completeSetup() override
337  {
338  double dt = 0.0;
339  mfem::Vector p = pressure();
340  contact_.updateForcesAndJacobian(cycle_, time_, dt, BasePhysics::shapeDisplacement(), displacement_, p);
341  updateContactForceOutputs();
342 
343  SolidMechanicsBase::completeSetup();
344  }
345 
351  mfem::HypreParVector pressure() const { return contact_.mergedPressures(); }
352 
353 #ifdef SMITH_USE_TRIBOL
360  const ContactInteraction& contactInteraction(int interaction_id) const
361  {
362  for (const auto& interaction : contact_.getContactInteractions()) {
363  if (interaction.getInteractionId() == interaction_id) {
364  return interaction;
365  }
366  }
367  SLIC_ERROR_ROOT(std::format("No contact interaction found with interaction_id={}", interaction_id));
368  return contact_.getContactInteractions().front();
369  }
370 #endif
371 
372  protected:
374  bool trySetDualAdjointBc(const std::string& dual_name, const smith::FiniteElementState& bc) override
375  {
376  const auto interaction_id = parseContactInteractionForceId(dual_name);
377  if (interaction_id.has_value()) {
378  auto it = contact_interaction_force_adjoint_bcs_.find(*interaction_id);
379  SLIC_ERROR_ROOT_IF(it == contact_interaction_force_adjoint_bcs_.end(),
380  std::format("No contact force adjoint BC registered for interaction_id={}", *interaction_id));
381  *it->second = bc;
382  return true;
383  }
384 
385  // No match found for dual_name in SolidMechanicsContact. Try the base class.
386  return SolidMechanicsBase::trySetDualAdjointBc(dual_name, bc);
387  }
388 
390  std::optional<int> parseContactInteractionForceId(std::string_view dual_name) const
391  {
392  constexpr std::string_view prefix = "contact_force_";
393  const std::string normalized_name = detail::removePrefix(this->name_, std::string(dual_name));
394  const std::string_view normalized_name_view{normalized_name};
395 
396  if (normalized_name_view.substr(0, prefix.size()) != prefix) {
397  return std::nullopt;
398  }
399 
400  return parseInteractionId(normalized_name_view.substr(prefix.size()));
401  }
402 
404  static std::optional<int> parseInteractionId(std::string_view id_str)
405  {
406  int interaction_id = -1;
407  const auto* begin = id_str.data();
408  const auto* end = id_str.data() + id_str.size();
409  auto [ptr, ec] = std::from_chars(begin, end, interaction_id);
410  if (ec != std::errc{} || ptr != end) {
411  return std::nullopt;
412  }
413  return interaction_id;
414  }
415 
417  void quasiStaticSolve(double dt) override
418  {
419  // warm start must be called prior to the time update so that the previous Jacobians can be used consistently
420  // throughout. warm start for contact needs to include the previous stiffness terms associated with contact
421  // otherwise the system will interpenetrate instantly on warm-starting.
422  warmStartDisplacementContact(dt);
423  time_ += dt;
424 
425  // In general, the solution vector is a stacked (block) vector:
426  // | displacement |
427  // | contact pressure |
428  // Contact pressure is only active when solving a contact problem with Lagrange multipliers.
429  mfem::Vector augmented_solution(displacement_.Size() + contact_.numPressureDofs());
430  augmented_solution.SetVector(displacement_, 0);
431  augmented_solution.SetVector(contact_.mergedPressures(), displacement_.Size());
432 
433  // solve the non-linear system resid = 0 and pressure * gap = 0
434  nonlin_solver_->solve(augmented_solution);
435  displacement_.Set(1.0, mfem::Vector(augmented_solution, 0, displacement_.Size()));
436  updateContactForceOutputs();
437  }
438 
473  {
475 
476  du_ = 0.0;
477  for (auto& bc : bcs_.essentials()) {
478  // apply the future boundary conditions, but use the most recent Jacobians stiffness.
479  bc.setDofs(du_, time_ + dt);
480  }
481 
482  auto& constrained_dofs = bcs_.allEssentialTrueDofs();
483  for (int i = 0; i < constrained_dofs.Size(); i++) {
484  int j = constrained_dofs[i];
485  du_[j] -= displacement_(j);
486  }
487 
488  auto amgf_prec = dynamic_cast<mfem::AMGFSolver*>(&nonlin_solver_->preconditioner());
489  if (amgf_prec) {
490  // compute contact_dof_prolongation
491  computeContactSubspaceTransferOperator();
492  // set AMGF subspace transfer operator
493  amgf_prec->SetFilteredSubspaceTransferOperator(*(contact_dof_prolongation_.get()));
494  // set the filteredsubspace solver component of AMGF
495  // better solution: retrieve print level from .preconditioner_print_level from linear_solver_options
496  int filter_solver_print_level = 0;
497  filter_solver_ =
498  std::make_unique<StrumpackSolver>(filter_solver_print_level, contact_dof_prolongation_->GetComm());
499  amgf_prec->SetFilteredSubspaceSolver(*filter_solver_.get());
500 
501  auto& lin_solver = nonlin_solver_->linearSolver();
502  auto iterative_solver = dynamic_cast<mfem::IterativeSolver*>(&lin_solver);
503  SLIC_WARNING_ROOT_IF(!iterative_solver,
504  "AMGFContact should only be used as a preconditioner for an iterative solver");
505  }
506 
507  if (use_warm_start_) {
508  // Update the system residual
509  mfem::Vector augmented_residual(displacement_.Size() + contact_.numPressureDofs());
510  augmented_residual = 0.0;
511  const mfem::Vector res = (*residual_)(time_ + dt, BasePhysics::shapeDisplacement(), displacement_, acceleration_,
512  *parameters_[parameter_indices].state...);
513 
514  mfem::Vector augmented_solution(displacement_.space().TrueVSize() + contact_.numPressureDofs());
515  augmented_solution = 0.0;
516  mfem::Vector du(augmented_solution, 0, displacement_.space().TrueVSize());
517  du = displacement_;
518  mfem::Vector p_blk(augmented_solution, displacement_.Size(), contact_.numPressureDofs());
519 
520  // Perform a single update for the warm start evaluation.
521  // Note: we use time_ to match the previous Jacobian evaluation point.
522  contact_.updateForcesAndJacobian(cycle_, time_, dt, BasePhysics::shapeDisplacement(), displacement_, p_blk);
523 
524  mfem::Vector r_blk(augmented_residual, 0, displacement_.space().TrueVSize());
525  r_blk = res;
526  r_blk += contact_.forces();
527 
528  mfem::Vector g_blk(augmented_residual, displacement_.Size(), contact_.numPressureDofs());
529  g_blk.Set(1.0, contact_.mergedGaps(true));
530 
531  r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
532 
533  // use the most recently evaluated Jacobian
534  auto [_, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(displacement_),
535  acceleration_, *parameters_[parameter_indices].previous_state...);
536 
537  if (contact_.haveLagrangeMultipliers()) {
538  J_offsets_ = mfem::Array<int>({0, displacement_.Size(), displacement_.Size() + contact_.numPressureDofs()});
539  J_constraint_ = contact_.jacobianFunction(assemble(drdu));
540 
541  // take ownership of blocks
542  J_constraint_->owns_blocks = false;
543  J_ = std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 0)));
544  J_12_ =
545  std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 1)));
546  J_21_ =
547  std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 0)));
548  J_22_ =
549  std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 1)));
550 
551  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
552  J_e_21_ = std::unique_ptr<mfem::HypreParMatrix>(J_21_->EliminateCols(bcs_.allEssentialTrueDofs()));
553  J_12_->EliminateRows(bcs_.allEssentialTrueDofs());
554 
555  J_operator_ = J_constraint_.get();
556  } else {
557  // get 11-block holding jacobian contributions
558  auto block_J = contact_.jacobianFunction(assemble(drdu));
559  block_J->owns_blocks = false;
560  J_ = std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&block_J->GetBlock(0, 0)));
561 
562  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
563 
564  J_operator_ = J_.get();
565  }
566 
567  augmented_residual *= -1.0;
568 
569  du = du_;
570  mfem::EliminateBC(*J_, *J_e_, constrained_dofs, du, r_blk);
571  for (int i = 0; i < constrained_dofs.Size(); i++) {
572  int j = constrained_dofs[i];
573  r_blk[j] = du[j];
574  }
575 
576  auto& lin_solver = nonlin_solver_->linearSolver();
577  lin_solver.SetOperator(*J_operator_);
578  lin_solver.Mult(augmented_residual, augmented_solution);
579 
580  du_ = du;
581  }
582 
583  displacement_ += du_;
584  }
585 
587  void quasiStaticAdjointSolve(double dt) override
588  {
589  SLIC_ERROR_ROOT_IF(contact_.haveLagrangeMultipliers(),
590  "Lagrange multiplier contact does not currently support sensitivities/adjoints.");
591 
592  contact_.updateForcesAndJacobian(cycle_, time_, dt, BasePhysics::shapeDisplacement(), displacement_);
593 
594  auto [_, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(displacement_),
595  acceleration_, *parameters_[parameter_indices].state...);
596 
597  auto block_J = contact_.jacobianFunction(assemble(drdu));
598  block_J->owns_blocks = false;
599  auto jacobian = std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&block_J->GetBlock(0, 0)));
600  auto J_T = std::unique_ptr<mfem::HypreParMatrix>(jacobian->Transpose());
601 
602  // If a QoI depends on per-interaction contact force duals, the dual-adjoint seeds define an additional contribution
603  // to the adjoint load:
604  // dJ/du += (df_i/du)^T * dJ/df_i
605  // Following SolidMechanics::setAdjointLoad() sign convention, displacement_adjoint_load_ stores the negative of the
606  // provided dJ/du, so we subtract these contributions here.
607 #ifdef SMITH_USE_TRIBOL
608  if (!contact_interaction_force_adjoint_bcs_.empty()) {
609  FiniteElementDual contact_force_load(displacement_.space(), "contact_force_dual_adjoint_load");
610  contact_force_load = 0.0;
611 
612  for (const auto& [interaction_id, force_seed] : contact_interaction_force_adjoint_bcs_) {
613  if (!force_seed) {
614  continue;
615  }
616 
617  // Only apply if the seed is nonzero.
618  if (force_seed->Norml2() == 0.0) {
619  continue;
620  }
621 
622  const auto interaction_J = contactInteraction(interaction_id).jacobianContribution();
623  auto* J00 = dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(0, 0));
624  SLIC_ERROR_ROOT_IF(!J00, "Expected HypreParMatrix (0,0) block for contact interaction Jacobian.");
625 
626  FiniteElementDual tmp(displacement_.space(), "contact_force_dual_adjoint_load_tmp");
627  tmp = 0.0;
628  J00->MultTranspose(*force_seed, tmp);
629  contact_force_load.Add(1.0, tmp);
630  }
631 
632  displacement_adjoint_load_.Add(-1.0, contact_force_load);
633  }
634 #endif
635 
636  auto J_e_T = bcs_.eliminateAllEssentialDofsFromMatrix(*J_T);
637  auto& constrained_dofs = bcs_.allEssentialTrueDofs();
638 
639  mfem::EliminateBC(*J_T, *J_e_T, constrained_dofs, reactions_adjoint_bcs_, displacement_adjoint_load_);
640  for (int i = 0; i < constrained_dofs.Size(); i++) {
641  int j = constrained_dofs[i];
642  displacement_adjoint_load_[j] = reactions_adjoint_bcs_[j];
643  }
644 
645  auto& lin_solver = nonlin_solver_->linearSolver();
646  lin_solver.SetOperator(*J_T);
647  lin_solver.Mult(displacement_adjoint_load_, adjoint_displacement_);
648  }
649 
652  {
653  double dt = this->getCheckpointedTimestep(cycle_);
654  contact_.updateForcesAndJacobian(cycle_ + 1, time_end_step_, dt, BasePhysics::shapeDisplacement(), displacement_);
655 
656  auto drdshape =
657  smith::get<DERIVATIVE>((*residual_)(time_end_step_, differentiate_wrt(BasePhysics::shapeDisplacement()),
658  displacement_, acceleration_, *parameters_[parameter_indices].state...));
659 
660  auto block_J = contact_.jacobianFunction(assemble(drdshape));
661  block_J->owns_blocks = false;
662  auto drdshape_mat =
663  std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&block_J->GetBlock(0, 0)));
664 
665  drdshape_mat->MultTranspose(adjoint_displacement_, shape_displacement_dual_);
666 
668  }
669 
672  {
673 #ifdef SMITH_USE_TRIBOL
674  for (const auto& interaction : contact_.getContactInteractions()) {
675  auto it = contact_interaction_forces_.find(interaction.getInteractionId());
676  if (it != contact_interaction_forces_.end()) {
677  const auto interaction_forces = interaction.forces();
678  *it->second = interaction_forces;
679  }
680  }
681 #endif
682  }
683 
684  using BasePhysics::bcs_;
685  using BasePhysics::cycle_;
686  using BasePhysics::duals_;
688  using BasePhysics::mesh_;
689  using BasePhysics::name_;
691  using BasePhysics::states_;
692  using BasePhysics::time_;
693  using SolidMechanicsBase::acceleration_;
694  using SolidMechanicsBase::adjoint_displacement_;
695  using SolidMechanicsBase::d_residual_d_;
696  using SolidMechanicsBase::DERIVATIVE;
697  using SolidMechanicsBase::displacement_;
698  using SolidMechanicsBase::displacement_adjoint_load_;
699  using SolidMechanicsBase::du_;
700  using SolidMechanicsBase::J_;
701  using SolidMechanicsBase::J_e_;
702  using SolidMechanicsBase::nonlin_solver_;
703  using SolidMechanicsBase::reactions_adjoint_bcs_;
704  using SolidMechanicsBase::residual_;
705  using SolidMechanicsBase::residual_with_bcs_;
706  using SolidMechanicsBase::shape_displacement_dual_;
707  using SolidMechanicsBase::time_end_step_;
708  using SolidMechanicsBase::use_warm_start_;
709  using SolidMechanicsBase::warmStartDisplacement;
710 
712  mfem::Operator* J_operator_;
713 
715  std::unique_ptr<mfem::HypreParMatrix> J_21_;
716 
718  std::unique_ptr<mfem::HypreParMatrix> J_12_;
719 
721  std::unique_ptr<mfem::HypreParMatrix> J_22_;
722 
724  mfem::Array<int> J_offsets_;
725 
727  std::unique_ptr<mfem::BlockOperator> J_constraint_;
728 
730  std::unique_ptr<mfem::HypreParMatrix> J_e_21_;
731 
734  std::unique_ptr<mfem::BlockOperator> J_constraint_e_;
735 
738 
740  std::unordered_map<int, std::unique_ptr<FiniteElementDual>> contact_interaction_forces_;
741 
743  std::unordered_map<int, std::unique_ptr<FiniteElementState>> contact_interaction_force_adjoint_bcs_;
744 
746  std::unique_ptr<mfem::HypreParMatrix> contact_dof_prolongation_;
747 
749  std::unique_ptr<mfem::Solver> filter_solver_;
750 };
751 
752 } // namespace smith
This is the abstract base class for a generic forward solver.
std::string name_
Name of the physics module.
std::shared_ptr< smith::Mesh > mesh_
The primary mesh.
std::vector< const smith::FiniteElementState * > states_
List of finite element primal states associated with this physics module.
int cycle_
Current cycle (forward pass time iteration count)
std::vector< const smith::FiniteElementDual * > duals_
List of finite element duals associated with this physics module.
std::vector< ParameterInfo > parameters_
A vector of the parameters associated with this physics module.
virtual const FiniteElementState & shapeDisplacement() const
Accessor for getting the shape displacement field from the physics modules.
double time_
Current time for the forward pass.
BoundaryConditionManager bcs_
Boundary condition manager instance.
bool is_quasistatic_
Whether the simulation is time-independent.
const FiniteElementDual & shapeDisplacementSensitivity() const
Internally used accessor for getting the shape displacement sensitivity from the physics modules.
const mfem::ParMesh & mfemParMesh() const
Returns a reference to the mfem ParMesh object.
This class stores all ContactInteractions for a problem, calls Tribol functions that act on all conta...
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.
std::unique_ptr< mfem::HypreParMatrix > J_12_
12 Jacobian block if using Lagrange multiplier contact (df/dp)
bool trySetDualAdjointBc(const std::string &dual_name, const smith::FiniteElementState &bc) override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unordered_map< int, std::unique_ptr< FiniteElementState > > contact_interaction_force_adjoint_bcs_
per-interaction dual-adjoint (BC) fields for contact force duals
static std::optional< int > parseInteractionId(std::string_view id_str)
Parses the interaction id from the suffix of a contact force dual name.
std::unique_ptr< mfem::HypreParMatrix > J_21_
21 Jacobian block if using Lagrange multiplier contact (dg/dx)
SolidMechanicsContact(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts, const smith::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_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 SolidMechanicsContact object.
std::unique_ptr< mfem::HypreParMatrix > J_e_21_
Columns of J_21_ that have been separated out because they are associated with essential boundary con...
std::vector< std::string > dualNames() const override
This is an overloaded member function, provided for convenience. It differs from the above function o...
SolidMechanicsContact(const SolidMechanicsInputOptions &input_options, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_mesh, int cycle=0, double time=0.0)
Construct a new Nonlinear SolidMechanicsContact Solver object.
std::unique_ptr< mfem::HypreParMatrix > J_22_
22 Jacobian block if using Lagrange multiplier contact (ones on diagonal for inactive t-dofs)
mfem::HypreParVector pressure() const
Get the contact pressures from all contact interactions, merged into a single HypreParVector.
mfem::Operator * J_operator_
Pointer to the Jacobian operator (J_ if no Lagrange multiplier contact, J_constraint_ otherwise)
const FiniteElementDual & computeTimestepShapeSensitivity() override
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< mfem::BlockOperator > J_constraint_
Assembled sparse matrix for the Jacobian with constraint blocks.
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...
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...
std::unordered_map< int, std::unique_ptr< FiniteElementDual > > contact_interaction_forces_
per-interaction contact forces for output
mfem::Array< int > J_offsets_
Block offsets for the J_constraint_ BlockOperator (must be owned outside J_constraint_)
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...
void warmStartDisplacementContact(double dt)
Sets the Dirichlet BCs for the current time and computes an initial guess for parameters and displace...
SolidMechanicsContact(std::unique_ptr< smith::EquationSolver > solver, const smith::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_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 SolidMechanicsContact object.
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::unique_ptr< mfem_ext::StdFunctionOperator > buildQuasistaticOperator() override
Build the quasi-static operator corresponding to the total Lagrangian formulation.
void addContactInteraction(int interaction_id, const std::set< int > &bdry_attr_surf1, const std::set< int > &bdry_attr_surf2, ContactOptions contact_opts)
Add a mortar contact boundary condition.
std::optional< int > parseContactInteractionForceId(std::string_view dual_name) const
Converts a dual name into an interaction id (if it exists)
static FiniteElementDual newDual(FunctionSpace space, const std::string &dual_name, const std::string &mesh_tag)
Factory method for creating a new FEDual object.
Class for storing contact data.
Accelerator functionality.
Definition: smith.cpp:36
dual(double, T) -> dual< T >
class template argument deduction guide for type dual.
auto differentiate_wrt(const mfem::Vector &v)
this function is intended to only be used in combination with smith::Functional::operator(),...
#define SMITH_MARK_FUNCTION
Definition: profiling.hpp:90
Tools for tagging a set of components of a vector field for boundary condition enforcement.
Stores the options for a contact pair.
Parameters for an iterative linear solution scheme.
Nonlinear solution scheme parameters.
a struct that is used in the physics modules to clarify which template arguments are user-controlled ...
Definition: common.hpp:45
Stores all information held in the input file that is used to configure the solver.
A timestep and boundary condition enforcement method for a dynamic solver.