Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
contact_constraint.hpp
Go to the documentation of this file.
1 // Copyright Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
14 #pragma once
15 
16 #include <vector>
17 
18 #include "smith/smith_config.hpp"
19 
20 #ifdef SMITH_USE_TRIBOL
21 
22 #include "tribol/interface/tribol.hpp"
23 
28 
29 namespace mfem {
30 class Vector;
31 class HypreParMatrix;
32 } // namespace mfem
33 
34 namespace smith {
35 
39 enum ContactFields
40 {
41  SHAPE,
42  DISP,
43 };
44 
53 static std::unique_ptr<mfem::HypreParMatrix> obtainBlock(mfem::BlockOperator* block_operator, int iblock, int jblock)
54 {
55  SLIC_ERROR_IF(iblock < 0 || jblock < 0, "block indicies must be non-negative");
56  SLIC_ERROR_IF(iblock > block_operator->NumRowBlocks() || jblock > block_operator->NumColBlocks(),
57  "one or more block indicies are too large and the requested block does not exist");
58  SLIC_ERROR_IF(block_operator->IsZeroBlock(iblock, jblock), "attempting to extract a null block");
59  auto Ablk = dynamic_cast<mfem::HypreParMatrix*>(&block_operator->GetBlock(iblock, jblock));
60  SLIC_ERROR_IF(!Ablk, "failed cast block to mfem::HypreParMatrix");
61  // deep copy --> unique_ptr
62  auto Ablk_unique = std::make_unique<mfem::HypreParMatrix>(*Ablk);
63  return Ablk_unique;
64 };
65 
66 class FiniteElementState;
67 
76 class ContactConstraint : public Constraint {
77  public:
88  ContactConstraint(int interaction_id, const mfem::ParMesh& mesh, const std::set<int>& bdry_attr_surf1,
89  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts,
90  const std::string& name = "contact_constraint")
91  : Constraint(name), contact_(mesh), contact_opts_{contact_opts}
92  {
93  contact_opts_.enforcement = ContactEnforcement::LagrangeMultiplier;
94  contact_.addContactInteraction(interaction_id, bdry_attr_surf1, bdry_attr_surf2, contact_opts_);
95  interaction_id_ = interaction_id;
96  }
97 
99  virtual ~ContactConstraint() {}
100 
110  mfem::Vector evaluate(double time, double dt, const std::vector<ConstFieldPtr>& fields,
111  bool update_fields = true) const override
112  {
113  if (update_fields) {
114  // note: Tribol does not use cycle.
115  int cycle = 0;
116  contact_.updateGaps(cycle, time, dt, *fields[ContactFields::SHAPE], *fields[ContactFields::DISP]);
117  auto gaps_hpv = contact_.mergedGaps(false);
118  // Note: this copy is needed to prevent the HypreParVector pointer from going out of scope. see
119  // https://github.com/mfem/mfem/issues/5029
120  cached_gap_.SetSize(gaps_hpv.Size());
121  cached_gap_.Set(1.0, gaps_hpv);
122  }
123  return cached_gap_;
124  };
125 
137  std::unique_ptr<mfem::HypreParMatrix> jacobian(double time, double dt, const std::vector<ConstFieldPtr>& fields,
138  int direction, bool update_fields = true,
139  [[maybe_unused]] bool fresh_derivative = true) const override
140  {
141  SLIC_ERROR_IF(direction != ContactFields::DISP, "requesting a non displacement-field derivative");
142  // if the field has been updated or we are requesting a fresh derivative
143  // then re-compute the gap Jacobian
144  // otherwise use previously cached Jacobian
145  if (update_fields || fresh_derivative) {
146  int cycle = 0;
147  if (update_fields) {
148  contact_.updateGaps(cycle, time, dt, *fields[ContactFields::SHAPE], *fields[ContactFields::DISP], true);
149  } else {
150  contact_.updateGaps(cycle, time, dt, std::nullopt, std::nullopt, true);
151  }
152  J_contact_ = contact_.mergedJacobian();
153  }
154  // obtain (1, 0) block entry from the 2 x 2 block contact linear system
155  auto dgdu = obtainBlock(J_contact_.get(), 1, 0);
156  return dgdu;
157  };
158 
172  mfem::Vector residual_contribution(double time, double dt, const std::vector<ConstFieldPtr>& fields,
173  const mfem::Vector& multipliers, int direction, bool update_fields = true,
174  bool fresh_derivative = true) const override
175  {
176  SLIC_ERROR_IF(direction != ContactFields::DISP, "requesting a non displacement-field derivative");
177  int cycle = 0;
178  if (update_fields || fresh_derivative) {
179  if (update_fields) {
180  contact_.update(cycle, time, dt, *fields[ContactFields::SHAPE], *fields[ContactFields::DISP], multipliers);
181  } else {
182  contact_.update(cycle, time, dt, std::nullopt, std::nullopt, multipliers);
183  }
184  }
185  return contact_.forces();
186  };
187 
200  std::unique_ptr<mfem::HypreParMatrix> residual_contribution_jacobian(double time, double dt,
201  const std::vector<ConstFieldPtr>& fields,
202  const mfem::Vector& multipliers, int direction,
203  bool update_fields = true,
204  bool fresh_derivative = true) const override
205  {
206  SLIC_ERROR_IF(direction != ContactFields::DISP, "requesting a non displacement-field derivative");
207 
208  int cycle = 0;
209  if (update_fields || fresh_derivative) {
210  if (update_fields) {
211  contact_.update(cycle, time, dt, *fields[ContactFields::SHAPE], *fields[ContactFields::DISP], multipliers);
212  } else {
213  contact_.update(cycle, time, dt, std::nullopt, std::nullopt, multipliers);
214  }
215  J_contact_ = contact_.mergedJacobian();
216  }
217  // obtain (0, 0) block entry from the 2 x 2 block contact linear system
218  auto Hessian = obtainBlock(J_contact_.get(), 0, 0);
219  return Hessian;
220  };
221 
232  std::unique_ptr<mfem::HypreParMatrix> jacobian_tilde(double time, double dt, const std::vector<ConstFieldPtr>& fields,
233  int direction, bool update_fields = true,
234  [[maybe_unused]] bool fresh_derivative = true) const override
235  {
236  SLIC_ERROR_IF(direction != ContactFields::DISP, "requesting a non displacement-field derivative");
237 
238  int cycle = 0;
239  if (update_fields || fresh_derivative) {
240  mfem::Vector p = contact_.mergedPressures();
241  if (update_fields) {
242  contact_.update(cycle, time, dt, *fields[ContactFields::SHAPE], *fields[ContactFields::DISP], p);
243  } else {
244  contact_.update(cycle, time, dt, std::nullopt, std::nullopt, p);
245  }
246  J_contact_ = contact_.mergedJacobian();
247  }
248  // obtain (0, 1) block entry from the 2 x 2 block contact linear system
249  auto dgduT = obtainBlock(J_contact_.get(), 0, 1);
250  std::unique_ptr<mfem::HypreParMatrix> dgdu(dgduT->Transpose());
251  return dgdu;
252  };
253 
254  int numPressureDofs() const { return contact_.numPressureDofs(); }
255 
256  protected:
260  mutable ContactData contact_;
261 
265  ContactOptions contact_opts_;
266 
270  int interaction_id_;
271 
275  mutable std::unique_ptr<mfem::BlockOperator> J_contact_;
276 
280  mutable mfem::Vector cached_gap_;
281 };
282 
283 } // namespace smith
284 
285 #endif
Specifies interface for evaluating distributed constriants from fields as well as their Jacobians and...
This file contains enumerations and record types for contact configuration.
Class for storing contact data.
Defines common types and helper functions for using the residual and scalar_objective classes.
Accelerator functionality.
Definition: smith.cpp:36