18 #include <string_view>
19 #include <system_error>
20 #include <unordered_map>
28 template <
int order,
int dim,
typename parameters = Parameters<>,
29 typename parameter_indices = std::make_
integer_sequence<
int, parameters::n>>
41 template <
int order,
int dim,
typename... parameter_space,
int... parameter_indices>
43 std::integer_sequence<int, parameter_indices...>>
44 :
public SolidMechanics<order, dim, Parameters<parameter_space...>,
45 std::integer_sequence<int, parameter_indices...>> {
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)
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),
109 std::shared_ptr<smith::Mesh> smith_mesh,
int cycle = 0,
double time = 0.0)
117 SolidMechanicsBase::resetStates(cycle, time);
118 for (
auto& [_, force] : contact_interaction_forces_) {
123 for (
auto& [_, seed] : contact_interaction_force_adjoint_bcs_) {
130 mfem::Vector p(contact_.numPressureDofs());
133 updateContactForceOutputs();
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());
142 *parameters_[parameter_indices].state...);
147 mfem::Vector r_blk(r, 0, displacement_.Size());
151 r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
155 if (contact_.haveLagrangeMultipliers()) {
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,
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...);
168 J_constraint_ = contact_.jacobianFunction(assemble(drdu));
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)));
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());
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());
191 J_operator_ = J_constraint_.get();
192 return *J_constraint_;
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...);
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)));
207 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
209 J_operator_ = J_.get();
224 const std::set<int>& bdry_attr_surf2,
ContactOptions contact_opts)
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.");
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_));
236 contact_.addContactInteraction(interaction_id, bdry_attr_surf1, bdry_attr_surf2, contact_opts);
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));
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));
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()));
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()) {
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()) {
290 return SolidMechanicsBase::dualAdjoint(dual_name);
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.");
300 const FiniteElementState checkpointed_displacement = this->loadCheckpointedState(
"displacement", cycle);
301 double dt = this->getCheckpointedTimestep(cycle);
303 updateContactForceOutputs();
305 const auto interaction_id = parseContactInteractionForceId(dual_name);
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);
311 it == contact_interaction_forces_.end(),
312 std::format(
"Requested checkpointed dual '{}' does not exist in physics module '{}'.", dual_name, name_));
316 return SolidMechanicsBase::loadCheckpointedDual(dual_name, cycle);
326 if (!contact_dof_prolongation_) {
327 contact_dof_prolongation_ = contact_.contactSubspaceTransferOperator();
339 mfem::Vector p = pressure();
341 updateContactForceOutputs();
343 SolidMechanicsBase::completeSetup();
351 mfem::HypreParVector
pressure()
const {
return contact_.mergedPressures(); }
353 #ifdef SMITH_USE_TRIBOL
360 const ContactInteraction& contactInteraction(
int interaction_id)
const
362 for (
const auto& interaction : contact_.getContactInteractions()) {
363 if (interaction.getInteractionId() == interaction_id) {
367 SLIC_ERROR_ROOT(std::format(
"No contact interaction found with interaction_id={}", interaction_id));
368 return contact_.getContactInteractions().front();
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));
386 return SolidMechanicsBase::trySetDualAdjointBc(dual_name, bc);
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};
396 if (normalized_name_view.substr(0, prefix.size()) != prefix) {
400 return parseInteractionId(normalized_name_view.substr(prefix.size()));
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) {
413 return interaction_id;
422 warmStartDisplacementContact(dt);
429 mfem::Vector augmented_solution(displacement_.Size() + contact_.numPressureDofs());
430 augmented_solution.SetVector(displacement_, 0);
431 augmented_solution.SetVector(contact_.mergedPressures(), displacement_.Size());
434 nonlin_solver_->solve(augmented_solution);
435 displacement_.Set(1.0, mfem::Vector(augmented_solution, 0, displacement_.Size()));
436 updateContactForceOutputs();
477 for (
auto& bc : bcs_.essentials()) {
479 bc.setDofs(du_, time_ + dt);
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);
488 auto amgf_prec =
dynamic_cast<mfem::AMGFSolver*
>(&nonlin_solver_->preconditioner());
491 computeContactSubspaceTransferOperator();
493 amgf_prec->SetFilteredSubspaceTransferOperator(*(contact_dof_prolongation_.get()));
496 int filter_solver_print_level = 0;
498 std::make_unique<StrumpackSolver>(filter_solver_print_level, contact_dof_prolongation_->GetComm());
499 amgf_prec->SetFilteredSubspaceSolver(*filter_solver_.get());
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");
507 if (use_warm_start_) {
509 mfem::Vector augmented_residual(displacement_.Size() + contact_.numPressureDofs());
510 augmented_residual = 0.0;
512 *parameters_[parameter_indices].state...);
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());
518 mfem::Vector p_blk(augmented_solution, displacement_.Size(), contact_.numPressureDofs());
524 mfem::Vector r_blk(augmented_residual, 0, displacement_.space().TrueVSize());
526 r_blk += contact_.forces();
528 mfem::Vector g_blk(augmented_residual, displacement_.Size(), contact_.numPressureDofs());
529 g_blk.Set(1.0, contact_.mergedGaps(
true));
531 r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
535 acceleration_, *parameters_[parameter_indices].previous_state...);
537 if (contact_.haveLagrangeMultipliers()) {
538 J_offsets_ = mfem::Array<int>({0, displacement_.Size(), displacement_.Size() + contact_.numPressureDofs()});
539 J_constraint_ = contact_.jacobianFunction(assemble(drdu));
542 J_constraint_->owns_blocks =
false;
543 J_ = std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(0, 0)));
545 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(0, 1)));
547 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(1, 0)));
549 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(1, 1)));
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());
555 J_operator_ = J_constraint_.get();
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)));
562 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
564 J_operator_ = J_.get();
567 augmented_residual *= -1.0;
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];
576 auto& lin_solver = nonlin_solver_->linearSolver();
577 lin_solver.SetOperator(*J_operator_);
578 lin_solver.Mult(augmented_residual, augmented_solution);
583 displacement_ += du_;
589 SLIC_ERROR_ROOT_IF(contact_.haveLagrangeMultipliers(),
590 "Lagrange multiplier contact does not currently support sensitivities/adjoints.");
595 acceleration_, *parameters_[parameter_indices].state...);
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());
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;
612 for (
const auto& [interaction_id, force_seed] : contact_interaction_force_adjoint_bcs_) {
618 if (force_seed->Norml2() == 0.0) {
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.");
626 FiniteElementDual tmp(displacement_.space(),
"contact_force_dual_adjoint_load_tmp");
628 J00->MultTranspose(*force_seed, tmp);
629 contact_force_load.Add(1.0, tmp);
632 displacement_adjoint_load_.Add(-1.0, contact_force_load);
636 auto J_e_T = bcs_.eliminateAllEssentialDofsFromMatrix(*J_T);
637 auto& constrained_dofs = bcs_.allEssentialTrueDofs();
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];
645 auto& lin_solver = nonlin_solver_->linearSolver();
646 lin_solver.SetOperator(*J_T);
647 lin_solver.Mult(displacement_adjoint_load_, adjoint_displacement_);
653 double dt = this->getCheckpointedTimestep(cycle_);
658 displacement_, acceleration_, *parameters_[parameter_indices].state...));
660 auto block_J = contact_.jacobianFunction(assemble(drdshape));
661 block_J->owns_blocks =
false;
663 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&block_J->GetBlock(0, 0)));
665 drdshape_mat->MultTranspose(adjoint_displacement_, shape_displacement_dual_);
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;
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;
715 std::unique_ptr<mfem::HypreParMatrix>
J_21_;
718 std::unique_ptr<mfem::HypreParMatrix>
J_12_;
721 std::unique_ptr<mfem::HypreParMatrix>
J_22_;
730 std::unique_ptr<mfem::HypreParMatrix>
J_e_21_;
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.
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.
The nonlinear solid solver class.
static FiniteElementDual newDual(FunctionSpace space, const std::string &dual_name, const std::string &mesh_tag)
Factory method for creating a new FEDual object.
Accelerator functionality.
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
Tools for tagging a set of components of a vector field for boundary condition enforcement.
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 ...
A timestep and boundary condition enforcement method for a dynamic solver.