Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
contact_data.cpp
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 
8 
9 #include <cstddef>
10 
11 #include "axom/slic.hpp"
12 #include "mpi.h"
13 
15 
16 #ifdef SMITH_USE_TRIBOL
17 #include "tribol/interface/tribol.hpp"
18 #include "tribol/interface/mfem_tribol.hpp"
19 #endif
20 
21 namespace smith {
22 
23 #ifdef SMITH_USE_TRIBOL
24 
25 ContactData::ContactData(const mfem::ParMesh& mesh)
26  : mesh_{mesh},
27  reference_nodes_{static_cast<const mfem::ParGridFunction*>(mesh.GetNodes())},
28  current_coords_{*reference_nodes_},
29  have_lagrange_multipliers_{false},
30  num_pressure_dofs_{0},
31  offsets_up_to_date_{false}
32 {
33 }
34 
35 ContactData::~ContactData() { tribol::finalize(); }
36 
37 void ContactData::addContactInteraction(int interaction_id, const std::set<int>& bdry_attr_surf1,
38  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts)
39 {
40  interactions_.emplace_back(interaction_id, mesh_, bdry_attr_surf1, bdry_attr_surf2, current_coords_, contact_opts);
41  if (contact_opts.enforcement == ContactEnforcement::LagrangeMultiplier) {
42  have_lagrange_multipliers_ = true;
43  num_pressure_dofs_ += interactions_.back().numPressureDofs();
44  offsets_up_to_date_ = false;
45  }
46 }
47 
48 void ContactData::reset()
49 {
50  for (auto& interaction : interactions_) {
51  FiniteElementState zero = interaction.pressure();
52  zero = 0.0;
53  interaction.setPressure(zero);
54  }
55 }
56 
57 void ContactData::update(int cycle, double time, double& dt)
58 {
59  cycle_ = cycle;
60  time_ = time;
61  dt_ = dt;
62  // This updates the redecomposed surface mesh based on the current displacement, then transfers field quantities to
63  // the updated mesh.
64  tribol::updateMfemParallelDecomposition();
65  // This function computes forces, gaps, and Jacobian contributions based on the current field quantities. Note the
66  // fields (with the exception of pressure) are stored on the redecomposed surface mesh until transferred by calling
67  // forces(), mergedGaps(), etc.
68  tribol::update(cycle, time, dt);
69 }
70 
71 FiniteElementDual ContactData::forces() const
72 {
73  FiniteElementDual f(*reference_nodes_->ParFESpace(), "contact force");
74  for (const auto& interaction : interactions_) {
75  f += interaction.forces();
76  }
77  return f;
78 }
79 
80 mfem::HypreParVector ContactData::mergedPressures() const
81 {
82  updateDofOffsets();
83  mfem::HypreParVector merged_p(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
84  global_pressure_dof_offsets_.GetData());
85  for (size_t i{0}; i < interactions_.size(); ++i) {
86  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
87  mfem::Vector p_interaction;
88  p_interaction.MakeRef(
89  merged_p, pressure_dof_offsets_[static_cast<int>(i)],
90  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
91  p_interaction.Set(1.0, interactions_[i].pressure());
92  }
93  }
94  return merged_p;
95 }
96 
97 mfem::HypreParVector ContactData::mergedGaps(bool zero_inactive) const
98 {
99  updateDofOffsets();
100  mfem::HypreParVector merged_g(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
101  global_pressure_dof_offsets_.GetData());
102  for (size_t i{0}; i < interactions_.size(); ++i) {
103  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
104  auto g = interactions_[i].gaps();
105  if (zero_inactive) {
106  for (auto dof : interactions_[i].inactiveDofs()) {
107  g[dof] = 0.0;
108  }
109  }
110  mfem::Vector g_interaction(
111  merged_g, pressure_dof_offsets_[static_cast<int>(i)],
112  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
113  g_interaction.Set(1.0, g);
114  }
115  }
116  return merged_g;
117 }
118 
119 std::unique_ptr<mfem::BlockOperator> ContactData::mergedJacobian() const
120 {
121  updateDofOffsets();
122  // this is the BlockOperator we are returning with the following blocks:
123  // | df_(contact)/dx df_(contact)/dp |
124  // | dg/dx I_(inactive) |
125  // where I_(inactive) is a matrix with ones on the diagonal of inactive pressure true degrees of freedom
126  auto block_J = std::make_unique<mfem::BlockOperator>(jacobian_offsets_);
127  block_J->owns_blocks = true;
128  // rather than returning different blocks for each contact interaction with Lagrange multipliers, merge them all into
129  // a single block
130  mfem::Array2D<const mfem::HypreParMatrix*> dgdu_blocks(static_cast<int>(interactions_.size()), 1);
131  mfem::Array2D<const mfem::HypreParMatrix*> dfdp_blocks(1, static_cast<int>(interactions_.size()));
132  for (size_t i{0}; i < interactions_.size(); ++i) {
133  dgdu_blocks(static_cast<int>(i), 0) = nullptr;
134  dfdp_blocks(0, static_cast<int>(i)) = nullptr;
135  }
136 
137  for (size_t i{0}; i < interactions_.size(); ++i) {
138  // this is the BlockOperator for one of the contact interactions
139  auto interaction_J = interactions_[i].jacobian();
140  interaction_J->owns_blocks = false; // we'll manage the ownership of the blocks on our own...
141  // add the contact interaction's contribution to df_(contact)/dx (the 0, 0 block)
142  if (!interaction_J->IsZeroBlock(0, 0)) {
143  SLIC_ERROR_ROOT_IF(!dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(0, 0)),
144  "Only HypreParMatrix constraint matrix blocks are currently supported.");
145  if (block_J->IsZeroBlock(0, 0)) {
146  block_J->SetBlock(0, 0, &interaction_J->GetBlock(0, 0));
147  } else {
148  block_J->SetBlock(0, 0,
149  mfem::Add(1.0, static_cast<mfem::HypreParMatrix&>(block_J->GetBlock(0, 0)), 1.0,
150  static_cast<mfem::HypreParMatrix&>(interaction_J->GetBlock(0, 0))));
151  delete &interaction_J->GetBlock(0, 0);
152  }
153  }
154  // add the contact interaction's (other) contribution to df_(contact)/dx (for penalty) or to df_(contact)/dp and
155  // dg/dx (for Lagrange multipliers)
156  if (!interaction_J->IsZeroBlock(1, 0) && !interaction_J->IsZeroBlock(0, 1)) {
157  auto dgdu = dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(1, 0));
158  auto dfdp = dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(0, 1));
159  SLIC_ERROR_ROOT_IF(!dgdu, "Only HypreParMatrix constraint matrix blocks are currently supported.");
160  SLIC_ERROR_ROOT_IF(!dfdp, "Only HypreParMatrix constraint matrix blocks are currently supported.");
161  // zero out rows and cols not in the active set
162  auto inactive_dofs = interactions_[i].inactiveDofs();
163  dgdu->EliminateRows(inactive_dofs);
164  dfdp->EliminateCols(inactive_dofs);
165  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::Penalty) {
166  // compute contribution to df_(contact)/dx (the 0, 0 block) for penalty
167  std::unique_ptr<mfem::HypreParMatrix> BTB(mfem::ParMult(dfdp, dgdu, true));
168  delete &interaction_J->GetBlock(1, 0);
169  delete &interaction_J->GetBlock(0, 1);
170  if (block_J->IsZeroBlock(0, 0)) {
171  mfem::Vector penalty(reference_nodes_->ParFESpace()->GetTrueVSize());
172  penalty = interactions_[i].getContactOptions().penalty;
173  BTB->ScaleRows(penalty);
174  block_J->SetBlock(0, 0, BTB.release());
175  } else {
176  block_J->SetBlock(0, 0,
177  mfem::Add(1.0, static_cast<mfem::HypreParMatrix&>(block_J->GetBlock(0, 0)),
178  interactions_[i].getContactOptions().penalty, *BTB));
179  }
180  } else // enforcement == ContactEnforcement::LagrangeMultiplier
181  {
182  // compute contribution to off-diagonal blocks for Lagrange multiplier
183  dgdu_blocks(static_cast<int>(i), 0) = dgdu;
184  dfdp_blocks(0, static_cast<int>(i)) = dfdp;
185  }
186  if (!interaction_J->IsZeroBlock(1, 1)) {
187  // we track our own active set, so get rid of the tribol inactive dof block
188  delete &interaction_J->GetBlock(1, 1);
189  }
190  }
191  }
192  if (haveLagrangeMultipliers()) {
193  // merge all of the contributions from all of the contact interactions
194  block_J->SetBlock(1, 0, mfem::HypreParMatrixFromBlocks(dgdu_blocks));
195  // store the transpose explicitly (rather than as a TransposeOperator) for solvers that need HypreParMatrixs
196  block_J->SetBlock(0, 1, mfem::HypreParMatrixFromBlocks(dfdp_blocks));
197  // explicitly delete the blocks
198  for (size_t i{0}; i < interactions_.size(); ++i) {
199  delete dgdu_blocks(static_cast<int>(i), 0);
200  delete dfdp_blocks(0, static_cast<int>(i));
201  }
202  // build I_(inactive): a diagonal matrix with ones on inactive dofs and zeros elsewhere
203  mfem::Array<const mfem::Array<int>*> inactive_tdofs_vector(static_cast<int>(interactions_.size()));
204  int inactive_tdofs_ct = 0;
205  for (int i{0}; i < inactive_tdofs_vector.Size(); ++i) {
206  inactive_tdofs_vector[i] = &interactions_[static_cast<size_t>(i)].inactiveDofs();
207  inactive_tdofs_ct += inactive_tdofs_vector[i]->Size();
208  }
209  mfem::Array<int> inactive_tdofs(inactive_tdofs_ct);
210  inactive_tdofs_ct = 0;
211  for (int i{0}; i < inactive_tdofs_vector.Size(); ++i) {
212  if (inactive_tdofs_vector[i]) {
213  for (int d{0}; d < inactive_tdofs_vector[i]->Size(); ++d) {
214  inactive_tdofs[d + inactive_tdofs_ct] = (*inactive_tdofs_vector[i])[d] + pressure_dof_offsets_[i];
215  }
216  inactive_tdofs_ct += inactive_tdofs_vector[i]->Size();
217  }
218  }
219  inactive_tdofs.GetMemory().SetHostPtrOwner(false);
220  mfem::Array<int> rows(numPressureDofs() + 1);
221  rows = 0;
222  inactive_tdofs_ct = 0;
223  for (int i{0}; i < numPressureDofs(); ++i) {
224  if (inactive_tdofs_ct < inactive_tdofs.Size() && inactive_tdofs[inactive_tdofs_ct] == i) {
225  ++inactive_tdofs_ct;
226  }
227  rows[i + 1] = inactive_tdofs_ct;
228  }
229  rows.GetMemory().SetHostPtrOwner(false);
230  mfem::Vector ones(inactive_tdofs_ct);
231  ones = 1.0;
232  ones.GetMemory().SetHostPtrOwner(false);
233  mfem::SparseMatrix inactive_diag(rows.GetData(), inactive_tdofs.GetData(), ones.GetData(), numPressureDofs(),
234  numPressureDofs(), false, false, true);
235  // if the size of ones is zero, SparseMatrix creates its own memory which it
236  // owns. explicitly prevent this...
237  inactive_diag.SetDataOwner(false);
238  auto block_1_1 =
239  new mfem::HypreParMatrix(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
240  global_pressure_dof_offsets_, &inactive_diag);
241  block_1_1->SetOwnerFlags(3, 3, 1);
242  block_J->SetBlock(1, 1, block_1_1);
243  // end building I_(inactive)
244  }
245  return block_J;
246 }
247 
248 void ContactData::residualFunction(const mfem::Vector& u_shape, const mfem::Vector& u, mfem::Vector& r)
249 {
250  const int disp_size = reference_nodes_->ParFESpace()->GetTrueVSize();
251 
252  // u_const should not change in this method; const cast is to create vector views which are copied to Tribol
253  // displacements and pressures and used to compute the (non-contact) residual
254  auto& u_const = const_cast<mfem::Vector&>(u);
255  const mfem::Vector u_blk(u_const, 0, disp_size);
256  const mfem::Vector p_blk(u_const, disp_size, numPressureDofs());
257 
258  mfem::Vector r_blk(r, 0, disp_size);
259  mfem::Vector g_blk(r, disp_size, numPressureDofs());
260 
261  setDisplacements(u_shape, u_blk);
262 
263  // we need to call update first to update gaps
264  for (auto& interaction : interactions_) {
265  interaction.evalJacobian(false);
266  }
267  update(cycle_, time_, dt_);
268  // with updated gaps, we can update pressure for contact interactions with penalty enforcement
269  setPressures(p_blk);
270  // call update again with the right pressures
271  for (auto& interaction : interactions_) {
272  interaction.evalJacobian(true);
273  }
274  update(cycle_, time_, dt_);
275 
276  r_blk += forces();
277  // calling mergedGaps() with true will zero out gap on inactive dofs (so the residual converges and the linearized
278  // system makes sense)
279  g_blk.Set(1.0, mergedGaps(true));
280 }
281 
282 std::unique_ptr<mfem::BlockOperator> ContactData::jacobianFunction(mfem::HypreParMatrix* orig_J) const
283 {
284  auto J_contact = mergedJacobian();
285  if (J_contact->IsZeroBlock(0, 0)) {
286  J_contact->SetBlock(0, 0, orig_J);
287  } else {
288  J_contact->SetBlock(0, 0,
289  mfem::Add(1.0, *orig_J, 1.0, static_cast<mfem::HypreParMatrix&>(J_contact->GetBlock(0, 0))));
290  }
291 
292  return J_contact;
293 }
294 
295 void ContactData::setPressures(const mfem::Vector& merged_pressures) const
296 {
297  updateDofOffsets();
298  for (size_t i{0}; i < interactions_.size(); ++i) {
299  FiniteElementState p_interaction(interactions_[i].pressureSpace());
300  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
301  // merged_pressures_const should not change; const cast is to create a vector view for copying to tribol pressures
302  auto& merged_pressures_const = const_cast<mfem::Vector&>(merged_pressures);
303  const mfem::Vector p_interaction_ref(
304  merged_pressures_const, pressure_dof_offsets_[static_cast<int>(i)],
305  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
306  p_interaction.Set(1.0, p_interaction_ref);
307  } else // enforcement == ContactEnforcement::Penalty
308  {
309  p_interaction.Set(interactions_[i].getContactOptions().penalty, interactions_[i].gaps());
310  }
311  for (auto dof : interactions_[i].inactiveDofs()) {
312  p_interaction[dof] = 0.0;
313  }
314  interactions_[i].setPressure(p_interaction);
315  }
316 }
317 
318 void ContactData::setDisplacements(const mfem::Vector& shape_u, const mfem::Vector& u)
319 {
320  mfem::ParGridFunction prolonged_shape_disp{current_coords_};
321  reference_nodes_->ParFESpace()->GetProlongationMatrix()->Mult(u, current_coords_);
322  reference_nodes_->ParFESpace()->GetProlongationMatrix()->Mult(shape_u, prolonged_shape_disp);
323 
324  current_coords_ += *reference_nodes_;
325  current_coords_ += prolonged_shape_disp;
326 }
327 
328 void ContactData::updateDofOffsets() const
329 {
330  if (offsets_up_to_date_) {
331  return;
332  }
333  jacobian_offsets_ = mfem::Array<int>({0, reference_nodes_->ParFESpace()->GetTrueVSize(),
334  numPressureDofs() + reference_nodes_->ParFESpace()->GetTrueVSize()});
335  pressure_dof_offsets_.SetSize(static_cast<int>(interactions_.size()) + 1);
336  pressure_dof_offsets_ = 0;
337  for (size_t i{0}; i < interactions_.size(); ++i) {
338  pressure_dof_offsets_[static_cast<int>(i + 1)] =
339  pressure_dof_offsets_[static_cast<int>(i)] + interactions_[i].numPressureDofs();
340  }
341  global_pressure_dof_offsets_.SetSize(mesh_.GetNRanks() + 1);
342  global_pressure_dof_offsets_ = 0;
343  global_pressure_dof_offsets_[mesh_.GetMyRank() + 1] = numPressureDofs();
344  MPI_Allreduce(MPI_IN_PLACE, global_pressure_dof_offsets_.GetData(), global_pressure_dof_offsets_.Size(), MPI_INT,
345  MPI_SUM, mesh_.GetComm());
346  for (int i{1}; i < mesh_.GetNRanks(); ++i) {
347  global_pressure_dof_offsets_[i + 1] += global_pressure_dof_offsets_[i];
348  }
349  if (HYPRE_AssumedPartitionCheck()) {
350  auto total_dofs = global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1];
351  // If the number of ranks is less than 2, ensure the size of global_pressure_dof_offsets_ is large enough
352  if (mesh_.GetNRanks() < 2) {
353  global_pressure_dof_offsets_.SetSize(3);
354  }
355  global_pressure_dof_offsets_[0] = global_pressure_dof_offsets_[mesh_.GetMyRank()];
356  global_pressure_dof_offsets_[1] = global_pressure_dof_offsets_[mesh_.GetMyRank() + 1];
357  global_pressure_dof_offsets_[2] = total_dofs;
358  // If the number of ranks is greater than 2, shrink the size of global_pressure_dof_offsets_
359  if (mesh_.GetNRanks() > 2) {
360  global_pressure_dof_offsets_.SetSize(3);
361  }
362  }
363  offsets_up_to_date_ = true;
364 }
365 
366 #else
367 
368 ContactData::ContactData([[maybe_unused]] const mfem::ParMesh& mesh)
369  : have_lagrange_multipliers_{false}, num_pressure_dofs_{0}
370 {
371 }
372 
374 
375 void ContactData::addContactInteraction([[maybe_unused]] int interaction_id,
376  [[maybe_unused]] const std::set<int>& bdry_attr_surf1,
377  [[maybe_unused]] const std::set<int>& bdry_attr_surf2,
378  [[maybe_unused]] ContactOptions contact_opts)
379 {
380  SLIC_WARNING_ROOT("Smith built without Tribol support. No contact interaction will be added.");
381 }
382 
383 void ContactData::update([[maybe_unused]] int cycle, [[maybe_unused]] double time, [[maybe_unused]] double& dt) {}
384 
386 {
387  FiniteElementDual f(*reference_nodes_->ParFESpace(), "contact force");
388  f = 0.0;
389  return f;
390 }
391 
392 mfem::HypreParVector ContactData::mergedPressures() const { return mfem::HypreParVector(); }
393 
394 mfem::HypreParVector ContactData::mergedGaps([[maybe_unused]] bool zero_inactive) const
395 {
396  return mfem::HypreParVector();
397 }
398 
399 std::unique_ptr<mfem::BlockOperator> ContactData::mergedJacobian() const
400 {
401  jacobian_offsets_ = mfem::Array<int>(
402  {0, reference_nodes_->ParFESpace()->GetTrueVSize(), reference_nodes_->ParFESpace()->GetTrueVSize()});
403  return std::make_unique<mfem::BlockOperator>(jacobian_offsets_);
404 }
405 
406 void ContactData::residualFunction([[maybe_unused]] const mfem::Vector& u_shape, [[maybe_unused]] const mfem::Vector& u,
407  [[maybe_unused]] mfem::Vector& r)
408 {
409 }
410 
411 std::unique_ptr<mfem::BlockOperator> ContactData::jacobianFunction(mfem::HypreParMatrix* orig_J) const
412 {
413  auto J_contact = mergedJacobian();
414  if (J_contact->IsZeroBlock(0, 0)) {
415  J_contact->SetBlock(0, 0, orig_J);
416  } else {
417  J_contact->SetBlock(0, 0,
418  mfem::Add(1.0, *orig_J, 1.0, static_cast<mfem::HypreParMatrix&>(J_contact->GetBlock(0, 0))));
419  }
420 
421  return J_contact;
422 }
423 
424 void ContactData::setPressures([[maybe_unused]] const mfem::Vector& true_pressures) const {}
425 
426 void ContactData::setDisplacements([[maybe_unused]] const mfem::Vector& u_shape,
427  [[maybe_unused]] const mfem::Vector& true_displacement)
428 {
429 }
430 
431 #endif
432 
433 } // namespace smith
void residualFunction(const mfem::Vector &u_shape, const mfem::Vector &u, mfem::Vector &r)
Computes the residual including contact terms.
mfem::HypreParVector mergedGaps(bool zero_inactive=false) const
Returns nodal gaps from all contact interactions on the contact surface true degrees of freedom.
void update(int cycle, double time, double &dt)
Updates the positions, forces, and Jacobian contributions associated with contact.
ContactData(const mfem::ParMesh &mesh)
The constructor.
void setPressures(const mfem::Vector &merged_pressures) const
Set the pressure field.
~ContactData()
Destructor to finalize Tribol.
mfem::HypreParVector mergedPressures() const
Returns pressures from all contact interactions on the contact surface true degrees of freedom.
std::unique_ptr< mfem::BlockOperator > jacobianFunction(mfem::HypreParMatrix *orig_J) const
Computes the Jacobian including contact terms, given the non-contact Jacobian terms.
FiniteElementDual forces() const
Get the contact constraint residual (i.e. nodal forces) from all contact interactions.
std::unique_ptr< mfem::BlockOperator > mergedJacobian() const
Returns a 2x2 block Jacobian on displacement/pressure true degrees of freedom from contact constraint...
void addContactInteraction(int interaction_id, const std::set< int > &bdry_attr_surf1, const std::set< int > &bdry_attr_surf2, ContactOptions contact_opts)
Add another contact interaction.
void setDisplacements(const mfem::Vector &u_shape, const mfem::Vector &u)
Update the current coordinates based on the new displacement field.
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
Class for storing contact data.
This file contains the declaration of structure that manages the MFEM objects that make up the state ...
Accelerator functionality.
Definition: smith.cpp:36
Stores the options for a contact pair.