20 template <
int order,
int dim,
typename parameters = Parameters<>,
21 typename parameter_indices = std::make_
integer_sequence<
int, parameters::n>>
33 template <
int order,
int dim,
typename... parameter_space,
int... parameter_indices>
35 std::integer_sequence<int, parameter_indices...>>
36 :
public SolidMechanics<order, dim, Parameters<parameter_space...>,
37 std::integer_sequence<int, parameter_indices...>> {
59 std::shared_ptr<smith::Mesh> smith_mesh, std::vector<std::string> parameter_names = {},
60 int cycle = 0,
double time = 0.0,
bool checkpoint_to_disk =
false,
bool use_warm_start =
true)
61 :
SolidMechanicsContact(std::make_unique<EquationSolver>(nonlinear_opts, lin_opts, smith_mesh->getComm()),
62 timestepping_opts, physics_name, smith_mesh, parameter_names, cycle, time,
63 checkpoint_to_disk, use_warm_start)
83 std::shared_ptr<smith::Mesh> smith_mesh, std::vector<std::string> parameter_names = {},
84 int cycle = 0,
double time = 0.0,
bool checkpoint_to_disk =
false,
bool use_warm_start =
true)
85 : SolidMechanicsBase(std::move(solver), timestepping_opts, physics_name, smith_mesh, parameter_names, cycle, time,
86 checkpoint_to_disk, use_warm_start),
88 forces_(
StateManager::newDual(displacement_.space(), detail::addPrefix(physics_name,
"contact_forces")))
91 duals_.push_back(&forces_);
104 std::shared_ptr<smith::Mesh> smith_mesh,
int cycle = 0,
double time = 0.0)
107 forces_(
StateManager::newDual(displacement_.space(), detail::addPrefix(physics_name,
"contact_forces")))
110 duals_.push_back(&forces_);
116 SolidMechanicsBase::resetStates(cycle, time);
121 contact_.update(cycle, time, dt);
127 auto residual_fn = [
this](
const mfem::Vector& u, mfem::Vector& r) {
128 const mfem::Vector u_blk(
const_cast<mfem::Vector&
>(u), 0, displacement_.Size());
130 *parameters_[parameter_indices].state...);
135 mfem::Vector r_blk(r, 0, displacement_.Size());
139 r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
143 if (contact_.haveLagrangeMultipliers()) {
146 J_offsets_ = mfem::Array<int>({0, displacement_.Size(), displacement_.Size() + contact_.numPressureDofs()});
147 return std::make_unique<mfem_ext::StdFunctionOperator>(
148 displacement_.space().TrueVSize() + contact_.numPressureDofs(), residual_fn,
150 [
this](
const mfem::Vector& u) -> mfem::Operator& {
151 const mfem::Vector u_blk(const_cast<mfem::Vector&>(u), 0, displacement_.Size());
152 auto [r, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(u_blk),
153 acceleration_, *parameters_[parameter_indices].state...);
157 J_constraint_ = contact_.jacobianFunction(J_.release());
160 J_constraint_->owns_blocks = false;
161 J_ = std::unique_ptr<mfem::HypreParMatrix>(
162 static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 0)));
163 J_12_ = std::unique_ptr<mfem::HypreParMatrix>(
164 static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(0, 1)));
165 J_21_ = std::unique_ptr<mfem::HypreParMatrix>(
166 static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 0)));
167 J_22_ = std::unique_ptr<mfem::HypreParMatrix>(
168 static_cast<mfem::HypreParMatrix*>(&J_constraint_->GetBlock(1, 1)));
171 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
172 J_e_21_ = std::unique_ptr<mfem::HypreParMatrix>(J_21_->EliminateCols(bcs_.allEssentialTrueDofs()));
173 J_12_->EliminateRows(bcs_.allEssentialTrueDofs());
176 J_constraint_e_ = std::make_unique<mfem::BlockOperator>(J_offsets_);
177 J_constraint_e_->SetBlock(0, 0, J_e_.get());
178 J_constraint_e_->SetBlock(1, 0, J_e_21_.get());
180 J_operator_ = J_constraint_.get();
181 return *J_constraint_;
186 return std::make_unique<mfem_ext::StdFunctionOperator>(
187 displacement_.space().TrueVSize(), residual_fn, [
this](
const mfem::Vector& u) -> mfem::Operator& {
188 auto [r, drdu] = (*residual_)(time_, BasePhysics::shapeDisplacement(), differentiate_wrt(u), acceleration_,
189 *parameters_[parameter_indices].state...);
193 auto block_J = contact_.jacobianFunction(J_.release());
194 block_J->owns_blocks = false;
195 J_ = std::unique_ptr<mfem::HypreParMatrix>(static_cast<mfem::HypreParMatrix*>(&block_J->GetBlock(0, 0)));
197 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
199 J_operator_ = J_.get();
214 const std::set<int>& bdry_attr_surf2,
ContactOptions contact_opts)
216 SLIC_ERROR_ROOT_IF(!is_quasistatic_,
"Contact can only be applied to quasistatic problems.");
217 SLIC_ERROR_ROOT_IF(order > 1,
"Contact can only be applied to linear (order = 1) meshes.");
218 contact_.addContactInteraction(interaction_id, bdry_attr_surf1, bdry_attr_surf2, contact_opts);
229 contact_.update(cycle_, time_, dt);
231 SolidMechanicsBase::completeSetup();
239 mfem::HypreParVector
pressure()
const {
return contact_.mergedPressures(); }
248 warmStartDisplacementContact(dt);
255 mfem::Vector augmented_solution(displacement_.Size() + contact_.numPressureDofs());
256 augmented_solution.SetVector(displacement_, 0);
257 augmented_solution.SetVector(contact_.mergedPressures(), displacement_.Size());
260 nonlin_solver_->solve(augmented_solution);
261 displacement_.Set(1.0, mfem::Vector(augmented_solution, 0, displacement_.Size()));
262 contact_.setPressures(mfem::Vector(augmented_solution, displacement_.Size(), contact_.numPressureDofs()));
263 contact_.update(cycle_, time_, dt);
264 forces_.SetVector(contact_.forces(), 0);
305 for (
auto& bc : bcs_.essentials()) {
307 bc.setDofs(du_, time_ + dt);
310 auto& constrained_dofs = bcs_.allEssentialTrueDofs();
311 for (
int i = 0; i < constrained_dofs.Size(); i++) {
312 int j = constrained_dofs[i];
313 du_[j] -= displacement_(j);
316 if (use_warm_start_) {
318 mfem::Vector augmented_residual(displacement_.Size() + contact_.numPressureDofs());
319 augmented_residual = 0.0;
321 *parameters_[parameter_indices].state...);
323 contact_.setPressures(mfem::Vector(augmented_residual, displacement_.Size(), contact_.numPressureDofs()));
324 contact_.update(cycle_, time_, dt);
325 mfem::Vector r_blk(augmented_residual, 0, displacement_.space().TrueVSize());
328 mfem::Vector augmented_solution(displacement_.space().TrueVSize() + contact_.numPressureDofs());
329 augmented_solution = 0.0;
330 mfem::Vector du(augmented_solution, 0, displacement_.space().TrueVSize());
334 r_blk.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
338 acceleration_, *parameters_[parameter_indices].previous_state...);
342 contact_.update(cycle_, time_, dt);
343 if (contact_.haveLagrangeMultipliers()) {
344 J_offsets_ = mfem::Array<int>({0, displacement_.Size(), displacement_.Size() + contact_.numPressureDofs()});
345 J_constraint_ = contact_.jacobianFunction(J_.release());
348 J_constraint_->owns_blocks =
false;
349 J_ = std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(0, 0)));
351 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(0, 1)));
353 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(1, 0)));
355 std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&J_constraint_->GetBlock(1, 1)));
358 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
359 J_e_21_ = std::unique_ptr<mfem::HypreParMatrix>(J_21_->EliminateCols(bcs_.allEssentialTrueDofs()));
360 J_12_->EliminateRows(bcs_.allEssentialTrueDofs());
362 J_operator_ = J_constraint_.get();
365 auto block_J = contact_.jacobianFunction(J_.release());
366 block_J->owns_blocks =
false;
367 J_ = std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&block_J->GetBlock(0, 0)));
369 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
371 J_operator_ = J_.get();
374 augmented_residual *= -1.0;
377 mfem::EliminateBC(*J_, *J_e_, constrained_dofs, du, r_blk);
378 for (
int i = 0; i < constrained_dofs.Size(); i++) {
379 int j = constrained_dofs[i];
383 auto& lin_solver = nonlin_solver_->linearSolver();
385 lin_solver.SetOperator(*J_operator_);
387 lin_solver.Mult(augmented_residual, augmented_solution);
392 displacement_ += du_;
398 SLIC_ERROR_ROOT_IF(contact_.haveLagrangeMultipliers(),
399 "Lagrange multiplier contact does not currently support sensitivities/adjoints.");
402 mfem::HypreParVector adjoint_essential(displacement_adjoint_load_);
403 adjoint_essential = 0.0;
406 acceleration_, *parameters_[parameter_indices].state...);
407 auto jacobian = assemble(drdu);
409 auto block_J = contact_.jacobianFunction(jacobian.release());
410 block_J->owns_blocks =
false;
411 jacobian = std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&block_J->GetBlock(0, 0)));
412 auto J_T = std::unique_ptr<mfem::HypreParMatrix>(jacobian->Transpose());
414 for (
const auto& bc : bcs_.essentials()) {
415 bc.apply(*J_T, displacement_adjoint_load_, adjoint_essential);
418 auto& lin_solver = nonlin_solver_->linearSolver();
419 lin_solver.SetOperator(*J_T);
420 lin_solver.Mult(displacement_adjoint_load_, adjoint_displacement_);
428 displacement_, acceleration_, *parameters_[parameter_indices].state...));
430 auto drdshape_mat = assemble(drdshape);
432 auto block_J = contact_.jacobianFunction(drdshape_mat.release());
433 block_J->owns_blocks =
false;
434 drdshape_mat = std::unique_ptr<mfem::HypreParMatrix>(
static_cast<mfem::HypreParMatrix*
>(&block_J->GetBlock(0, 0)));
436 drdshape_mat->MultTranspose(adjoint_displacement_, shape_displacement_dual_);
450 using SolidMechanicsBase::acceleration_;
451 using SolidMechanicsBase::adjoint_displacement_;
452 using SolidMechanicsBase::d_residual_d_;
453 using SolidMechanicsBase::DERIVATIVE;
454 using SolidMechanicsBase::displacement_;
455 using SolidMechanicsBase::displacement_adjoint_load_;
456 using SolidMechanicsBase::du_;
457 using SolidMechanicsBase::J_;
458 using SolidMechanicsBase::J_e_;
459 using SolidMechanicsBase::nonlin_solver_;
460 using SolidMechanicsBase::residual_;
461 using SolidMechanicsBase::residual_with_bcs_;
462 using SolidMechanicsBase::shape_displacement_dual_;
463 using SolidMechanicsBase::time_end_step_;
464 using SolidMechanicsBase::use_warm_start_;
465 using SolidMechanicsBase::warmStartDisplacement;
471 std::unique_ptr<mfem::HypreParMatrix>
J_21_;
474 std::unique_ptr<mfem::HypreParMatrix>
J_12_;
477 std::unique_ptr<mfem::HypreParMatrix>
J_22_;
486 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...
The nonlinear solid solver class.
Manages the lifetimes of FEState objects such that restarts are abstracted from physics modules.
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.
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.