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 #include <format>
11 #include <memory>
12 
13 #include "axom/slic.hpp"
14 #include "mpi.h"
15 
17 
18 #ifdef SMITH_USE_TRIBOL
19 #include "tribol/interface/tribol.hpp"
20 #include "tribol/interface/mfem_tribol.hpp"
21 #include "tribol/mesh/CouplingScheme.hpp"
22 #endif
23 
24 namespace smith {
25 
26 #ifdef SMITH_USE_TRIBOL
27 
28 ContactData::ContactData(const mfem::ParMesh& mesh)
29  : mesh_{mesh},
30  reference_nodes_{static_cast<const mfem::ParGridFunction*>(mesh.GetNodes())},
31  current_coords_{*reference_nodes_},
32  have_lagrange_multipliers_{false},
33  num_pressure_dofs_{0},
34  offsets_up_to_date_{false}
35 {
36 }
37 
38 ContactData::~ContactData() { tribol::finalize(); }
39 
40 void ContactData::addContactInteraction(int interaction_id, const std::set<int>& bdry_attr_surf1,
41  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts)
42 {
43  // Disallow duplicate ids globally: Tribol coupling schemes are keyed by cs_id (interaction_id) globally.
44  auto* cs = tribol::CouplingSchemeManager::getInstance().findData(static_cast<tribol::IndexT>(interaction_id));
45  SLIC_ERROR_ROOT_IF(cs != nullptr,
46  std::format("Contact interaction id {} is already registered with Tribol.", interaction_id));
47 
48  interactions_.emplace_back(interaction_id, mesh_, bdry_attr_surf1, bdry_attr_surf2, current_coords_, contact_opts);
49  if (contact_opts.enforcement == ContactEnforcement::LagrangeMultiplier) {
50  have_lagrange_multipliers_ = true;
51  num_pressure_dofs_ += interactions_.back().numPressureDofs();
52  offsets_up_to_date_ = false;
53  }
54  // specify all contact boundaries
55  mfem::Array<int> contact_bdry_attribs;
56  contact_bdry_attribs.SetSize(mesh_.bdr_attributes.Max());
57  contact_bdry_attribs = 0;
58  // attributes start at 1,
59  // shift by -1 to account for zero-based array indexing
60  for (const auto& bdry_attr : bdry_attr_surf1) {
61  contact_bdry_attribs[bdry_attr - 1] = 1;
62  }
63  for (const auto& bdry_attr : bdry_attr_surf2) {
64  contact_bdry_attribs[bdry_attr - 1] = 1;
65  }
66  // dofs for the current contact interaction
67  mfem::Array<int> contact_interaction_dofs;
68  reference_nodes_->ParFESpace()->GetEssentialTrueDofs(contact_bdry_attribs, contact_interaction_dofs);
69  // add dofs for current contact interaction call to all contact_dofs_
70  contact_dofs_.Append(contact_interaction_dofs);
71  // sort and delete duplicates
72  contact_dofs_.Sort();
73  contact_dofs_.Unique();
74 }
75 
76 void ContactData::reset()
77 {
78  for (auto& interaction : interactions_) {
79  FiniteElementState zero = interaction.pressure();
80  zero = 0.0;
81  interaction.setPressure(zero);
82  }
83 }
84 
85 void ContactData::updateGaps(int cycle, double time, double& dt,
86  std::optional<std::reference_wrapper<const mfem::Vector>> u_shape,
87  std::optional<std::reference_wrapper<const mfem::Vector>> u, bool eval_jacobian)
88 {
89  cycle_ = cycle;
90  time_ = time;
91  dt_ = dt;
92 
93  if (u_shape && u) {
94  setDisplacements(u_shape->get(), u->get());
95  }
96 
97  for (auto& interaction : interactions_) {
98  interaction.evalJacobian(eval_jacobian);
99  }
100  // This updates the redecomposed surface mesh based on the current displacement, then transfers field quantities to
101  // the updated mesh.
102  if (u_shape && u) {
103  tribol::updateMfemParallelDecomposition();
104  }
105  // This function computes gaps (and optionally geometric Jacobian blocks) based on the current mesh.
106  tribol::update(cycle, time, dt);
107 }
108 
109 void ContactData::update(int cycle, double time, double& dt,
110  std::optional<std::reference_wrapper<const mfem::Vector>> u_shape,
111  std::optional<std::reference_wrapper<const mfem::Vector>> u,
112  std::optional<std::reference_wrapper<const mfem::Vector>> p)
113 {
114  // First pass: update gaps if coordinates are provided
115  if (u_shape && u) {
116  updateGaps(cycle, time, dt, u_shape, u, false);
117  } else {
118  // Ensure internal timing is updated even if coordinates are not
119  cycle_ = cycle;
120  time_ = time;
121  dt_ = dt;
122  }
123 
124  // second pass: update pressures and compute forces/Jacobians if p is provided
125  if (p) {
126  // with updated gaps, we can update pressure for contact interactions (active set detection and penalty)
127  setPressures(p->get());
128 
129  for (auto& interaction : interactions_) {
130  interaction.evalJacobian(true);
131  }
132  // This second call is required to synchronize the updated pressures to Tribol's internal redecomposed surface mesh
133  // and to ensure Tribol's internal state is correctly reset for the second pass.
134  tribol::updateMfemParallelDecomposition();
135  tribol::update(cycle, time, dt);
136  }
137 }
138 
139 FiniteElementDual ContactData::forces() const
140 {
141  FiniteElementDual f(*reference_nodes_->ParFESpace(), "contact force");
142  for (const auto& interaction : interactions_) {
143  f += interaction.forces();
144  }
145  return f;
146 }
147 
148 mfem::HypreParVector ContactData::mergedPressures() const
149 {
150  updateDofOffsets();
151  mfem::HypreParVector merged_p(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
152  global_pressure_dof_offsets_.GetData());
153  for (size_t i{0}; i < interactions_.size(); ++i) {
154  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
155  mfem::Vector p_interaction;
156  p_interaction.MakeRef(
157  merged_p, pressure_dof_offsets_[static_cast<int>(i)],
158  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
159  p_interaction.Set(1.0, interactions_[i].pressure());
160  }
161  }
162  return merged_p;
163 }
164 
165 mfem::HypreParVector ContactData::mergedGaps(bool zero_inactive) const
166 {
167  updateDofOffsets();
168  mfem::HypreParVector merged_g(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
169  global_pressure_dof_offsets_.GetData());
170  for (size_t i{0}; i < interactions_.size(); ++i) {
171  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
172  auto g = interactions_[i].gaps();
173  if (zero_inactive) {
174  for (auto dof : interactions_[i].inactiveDofs()) {
175  g[dof] = 0.0;
176  }
177  }
178  mfem::Vector g_interaction(
179  merged_g, pressure_dof_offsets_[static_cast<int>(i)],
180  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
181  g_interaction.Set(1.0, g);
182  }
183  }
184  return merged_g;
185 }
186 
187 std::unique_ptr<mfem::BlockOperator> ContactData::mergedJacobian() const
188 {
189  updateDofOffsets();
190  // this is the BlockOperator we are returning with the following blocks:
191  // | df_(contact)/dx df_(contact)/dp |
192  // | dg/dx I_(inactive) |
193  // where I_(inactive) is a matrix with ones on the diagonal of inactive pressure true degrees of freedom
194  auto block_J = std::make_unique<mfem::BlockOperator>(jacobian_offsets_);
195  block_J->owns_blocks = true;
196  // rather than returning different blocks for each contact interaction with Lagrange multipliers, merge them all into
197  // a single block
198  mfem::Array2D<const mfem::HypreParMatrix*> dgdu_blocks(static_cast<int>(interactions_.size()), 1);
199  mfem::Array2D<const mfem::HypreParMatrix*> dfdp_blocks(1, static_cast<int>(interactions_.size()));
200  for (size_t i{0}; i < interactions_.size(); ++i) {
201  dgdu_blocks(static_cast<int>(i), 0) = nullptr;
202  dfdp_blocks(0, static_cast<int>(i)) = nullptr;
203  }
204 
205  for (size_t i{0}; i < interactions_.size(); ++i) {
206  // this is the BlockOperator for one of the contact interactions, post-processed for Smith's assembly conventions
207  auto interaction_J = interactions_[i].jacobianContribution();
208  interaction_J->owns_blocks = false; // we'll manage the ownership of the blocks on our own...
209 
210  // add the contact interaction's contribution to df_(contact)/dx (the 0, 0 block)
211  if (!interaction_J->IsZeroBlock(0, 0)) {
212  SLIC_ERROR_ROOT_IF(!dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(0, 0)),
213  "Only HypreParMatrix constraint matrix blocks are currently supported.");
214  if (block_J->IsZeroBlock(0, 0)) {
215  block_J->SetBlock(0, 0, &interaction_J->GetBlock(0, 0));
216  } else {
217  block_J->SetBlock(0, 0,
218  mfem::Add(1.0, static_cast<mfem::HypreParMatrix&>(block_J->GetBlock(0, 0)), 1.0,
219  static_cast<mfem::HypreParMatrix&>(interaction_J->GetBlock(0, 0))));
220  delete &interaction_J->GetBlock(0, 0);
221  }
222  }
223 
224  // add the contact interaction's contribution to df_(contact)/dp and dg/dx (for Lagrange multipliers)
225  if (!interaction_J->IsZeroBlock(1, 0) && !interaction_J->IsZeroBlock(0, 1)) {
226  auto dgdu = dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(1, 0));
227  auto dfdp = dynamic_cast<mfem::HypreParMatrix*>(&interaction_J->GetBlock(0, 1));
228  SLIC_ERROR_ROOT_IF(!dgdu, "Only HypreParMatrix constraint matrix blocks are currently supported.");
229  SLIC_ERROR_ROOT_IF(!dfdp, "Only HypreParMatrix constraint matrix blocks are currently supported.");
230 
231  dgdu_blocks(static_cast<int>(i), 0) = dgdu;
232  dfdp_blocks(0, static_cast<int>(i)) = dfdp;
233  }
234  }
235  if (haveLagrangeMultipliers()) {
236  // merge all of the contributions from all of the contact interactions
237  block_J->SetBlock(1, 0, mfem::HypreParMatrixFromBlocks(dgdu_blocks));
238  // store the transpose explicitly (rather than as a TransposeOperator) for solvers that need HypreParMatrixs
239  block_J->SetBlock(0, 1, mfem::HypreParMatrixFromBlocks(dfdp_blocks));
240  // explicitly delete the blocks
241  for (size_t i{0}; i < interactions_.size(); ++i) {
242  delete dgdu_blocks(static_cast<int>(i), 0);
243  delete dfdp_blocks(0, static_cast<int>(i));
244  }
245  // build I_(inactive): a diagonal matrix with ones on inactive dofs and zeros elsewhere
246  mfem::Array<const mfem::Array<int>*> inactive_tdofs_vector(static_cast<int>(interactions_.size()));
247  int inactive_tdofs_ct = 0;
248  for (int i{0}; i < inactive_tdofs_vector.Size(); ++i) {
249  inactive_tdofs_vector[i] = &interactions_[static_cast<size_t>(i)].inactiveDofs();
250  inactive_tdofs_ct += inactive_tdofs_vector[i]->Size();
251  }
252  mfem::Array<int> inactive_tdofs(inactive_tdofs_ct);
253  inactive_tdofs_ct = 0;
254  for (int i{0}; i < inactive_tdofs_vector.Size(); ++i) {
255  if (inactive_tdofs_vector[i]) {
256  for (int d{0}; d < inactive_tdofs_vector[i]->Size(); ++d) {
257  inactive_tdofs[d + inactive_tdofs_ct] = (*inactive_tdofs_vector[i])[d] + pressure_dof_offsets_[i];
258  }
259  inactive_tdofs_ct += inactive_tdofs_vector[i]->Size();
260  }
261  }
262  mfem::Array<int> rows(numPressureDofs() + 1);
263  rows = 0;
264  inactive_tdofs_ct = 0;
265  for (int i{0}; i < numPressureDofs(); ++i) {
266  if (inactive_tdofs_ct < inactive_tdofs.Size() && inactive_tdofs[inactive_tdofs_ct] == i) {
267  ++inactive_tdofs_ct;
268  }
269  rows[i + 1] = inactive_tdofs_ct;
270  }
271  mfem::Vector ones(inactive_tdofs_ct);
272  ones = 1.0;
273  mfem::SparseMatrix inactive_diag(rows.GetData(), inactive_tdofs.GetData(), ones.GetData(), numPressureDofs(),
274  numPressureDofs(), false, false, true);
275  rows.GetMemory().ClearOwnerFlags();
276  inactive_tdofs.GetMemory().ClearOwnerFlags();
277  ones.GetMemory().ClearOwnerFlags();
278  inactive_diag.GetMemoryI().ClearOwnerFlags();
279  inactive_diag.GetMemoryJ().ClearOwnerFlags();
280  inactive_diag.GetMemoryData().ClearOwnerFlags();
281  auto block_1_1 =
282  new mfem::HypreParMatrix(mesh_.GetComm(), global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1],
283  global_pressure_dof_offsets_, &inactive_diag);
284  constexpr int mfem_owned_host_flag = 3;
285  block_1_1->SetOwnerFlags(mfem_owned_host_flag, block_1_1->OwnsOffd(), block_1_1->OwnsColMap());
286  block_J->SetBlock(1, 1, block_1_1);
287  // end building I_(inactive)
288  }
289  return block_J;
290 }
291 
292 void ContactData::residualFunction(const mfem::Vector& u_shape, const mfem::Vector& u, mfem::Vector& r)
293 {
294  const int disp_size = reference_nodes_->ParFESpace()->GetTrueVSize();
295 
296  // u_const should not change in this method; const cast is to create vector views which are copied to Tribol
297  // displacements and pressures and used to compute the (non-contact) residual
298  auto& u_const = const_cast<mfem::Vector&>(u);
299  const mfem::Vector u_blk(u_const, 0, disp_size);
300  const mfem::Vector p_blk(u_const, disp_size, numPressureDofs());
301 
302  mfem::Vector r_blk(r, 0, disp_size);
303  mfem::Vector g_blk(r, disp_size, numPressureDofs());
304 
305  update(cycle_, time_, dt_, u_shape, u_blk, p_blk);
306 
307  r_blk += forces();
308  // calling mergedGaps() with true will zero out gap on inactive dofs (so the residual converges and the linearized
309  // system makes sense)
310  g_blk.Set(1.0, mergedGaps(true));
311 }
312 
313 std::unique_ptr<mfem::BlockOperator> ContactData::jacobianFunction(std::unique_ptr<mfem::HypreParMatrix> orig_J) const
314 {
315  auto J_contact = mergedJacobian();
316  if (J_contact->IsZeroBlock(0, 0)) {
317  J_contact->SetBlock(0, 0, orig_J.release());
318  } else {
319  J_contact->SetBlock(0, 0,
320  mfem::Add(1.0, *orig_J, 1.0, static_cast<mfem::HypreParMatrix&>(J_contact->GetBlock(0, 0))));
321  }
322 
323  return J_contact;
324 }
325 
326 void ContactData::setPressures(const mfem::Vector& merged_pressures) const
327 {
328  updateDofOffsets();
329  for (size_t i{0}; i < interactions_.size(); ++i) {
330  FiniteElementState p_interaction(interactions_[i].pressureSpace());
331  if (interactions_[i].getContactOptions().enforcement == ContactEnforcement::LagrangeMultiplier) {
332  // merged_pressures_const should not change; const cast is to create a vector view for copying to tribol pressures
333  auto& merged_pressures_const = const_cast<mfem::Vector&>(merged_pressures);
334  const mfem::Vector p_interaction_ref(
335  merged_pressures_const, pressure_dof_offsets_[static_cast<int>(i)],
336  pressure_dof_offsets_[static_cast<int>(i) + 1] - pressure_dof_offsets_[static_cast<int>(i)]);
337  p_interaction.Set(1.0, p_interaction_ref);
338  } else // enforcement == ContactEnforcement::Penalty
339  {
340  p_interaction.Set(interactions_[i].getContactOptions().penalty, interactions_[i].gaps());
341  }
342  for (auto dof : interactions_[i].inactiveDofs()) {
343  p_interaction[dof] = 0.0;
344  }
345  interactions_[i].setPressure(p_interaction);
346  }
347 }
348 
349 void ContactData::setDisplacements(const mfem::Vector& shape_u, const mfem::Vector& u)
350 {
351  mfem::ParGridFunction prolonged_shape_disp{current_coords_};
352  reference_nodes_->ParFESpace()->GetProlongationMatrix()->Mult(u, current_coords_);
353  reference_nodes_->ParFESpace()->GetProlongationMatrix()->Mult(shape_u, prolonged_shape_disp);
354 
355  current_coords_ += *reference_nodes_;
356  current_coords_ += prolonged_shape_disp;
357 }
358 
359 void ContactData::updateDofOffsets() const
360 {
361  if (offsets_up_to_date_) {
362  return;
363  }
364  jacobian_offsets_ = mfem::Array<int>({0, reference_nodes_->ParFESpace()->GetTrueVSize(),
365  numPressureDofs() + reference_nodes_->ParFESpace()->GetTrueVSize()});
366  pressure_dof_offsets_.SetSize(static_cast<int>(interactions_.size()) + 1);
367  pressure_dof_offsets_ = 0;
368  for (size_t i{0}; i < interactions_.size(); ++i) {
369  pressure_dof_offsets_[static_cast<int>(i + 1)] =
370  pressure_dof_offsets_[static_cast<int>(i)] + interactions_[i].numPressureDofs();
371  }
372  global_pressure_dof_offsets_.SetSize(mesh_.GetNRanks() + 1);
373  global_pressure_dof_offsets_ = 0;
374  global_pressure_dof_offsets_[mesh_.GetMyRank() + 1] = numPressureDofs();
375  MPI_Allreduce(MPI_IN_PLACE, global_pressure_dof_offsets_.GetData(), global_pressure_dof_offsets_.Size(), MPI_INT,
376  MPI_SUM, mesh_.GetComm());
377  for (int i{1}; i < mesh_.GetNRanks(); ++i) {
378  global_pressure_dof_offsets_[i + 1] += global_pressure_dof_offsets_[i];
379  }
380  if (HYPRE_AssumedPartitionCheck()) {
381  auto total_dofs = global_pressure_dof_offsets_[global_pressure_dof_offsets_.Size() - 1];
382  // If the number of ranks is less than 2, ensure the size of global_pressure_dof_offsets_ is large enough
383  if (mesh_.GetNRanks() < 2) {
384  global_pressure_dof_offsets_.SetSize(3);
385  }
386  global_pressure_dof_offsets_[0] = global_pressure_dof_offsets_[mesh_.GetMyRank()];
387  global_pressure_dof_offsets_[1] = global_pressure_dof_offsets_[mesh_.GetMyRank() + 1];
388  global_pressure_dof_offsets_[2] = total_dofs;
389  // If the number of ranks is greater than 2, shrink the size of global_pressure_dof_offsets_
390  if (mesh_.GetNRanks() > 2) {
391  global_pressure_dof_offsets_.SetSize(3);
392  }
393  }
394  offsets_up_to_date_ = true;
395 }
396 
397 std::unique_ptr<mfem::HypreParMatrix> ContactData::contactSubspaceTransferOperator()
398 {
399  const MPI_Comm comm = reference_nodes_->ParFESpace()->GetComm();
400  HYPRE_BigInt* col_offsets = reference_nodes_->ParFESpace()->GetTrueDofOffsets();
401  HYPRE_BigInt ncols_glb = reference_nodes_->ParFESpace()->GlobalTrueVSize();
402 
403  // number of rows of the restriction
404  // operator owned by the local MPI process
405  int nrows_loc = contact_dofs_.Size();
406  // should nrows_glb be of type HYPRE_BigInt?
407  // global number of rows of restriction
408  // operator
409  int nrows_glb = 0;
410  MPI_Allreduce(&nrows_loc, &nrows_glb, 1, MPI_INT, MPI_SUM, comm);
411  // determine rows offsets of the restriction operator
412  int row_offset = 0;
413  MPI_Scan(&nrows_loc, &row_offset, 1, MPI_INT, MPI_SUM, comm);
414  row_offset -= nrows_loc;
415  HYPRE_BigInt row_offsets[2];
416  row_offsets[0] = row_offset;
417  row_offsets[1] = row_offset + nrows_loc;
418 
419  // create mfem::SparseMatrix restriction matrix
420  // restriction from displacement dofs to
421  // contact dofs
422  // one nonzero (unit) entry per row
423  mfem::SparseMatrix Rsparse(nrows_loc, ncols_glb);
424  mfem::Array<int> col(1);
425  col = 0;
426  mfem::Vector entry(1);
427  entry = 1.0;
428  HYPRE_BigInt col_offset = col_offsets[0];
429  for (int k = 0; k < nrows_loc; k++) {
430  // local per process contact dof
431  // to global column number
432  col[0] = col_offset + contact_dofs_[k];
433  Rsparse.SetRow(k, col, entry);
434  }
435  Rsparse.Finalize();
436 
437  // convert local sparse restriction matrix
438  // to distributed mfem::HypreParMatrix
439  int* I = Rsparse.GetI();
440  HYPRE_BigInt* J = Rsparse.GetJ();
441  double* data = Rsparse.GetData();
442 
443  std::unique_ptr<mfem::HypreParMatrix> restriction_operator = std::make_unique<mfem::HypreParMatrix>(
444  comm, nrows_loc, nrows_glb, ncols_glb, I, J, data, row_offsets, col_offsets);
445  // convert restriction operator
446  // to a contact dof to displacement
447  // dof transfer operator
448  std::unique_ptr<mfem::HypreParMatrix> transfer_operator;
449  transfer_operator.reset(restriction_operator->Transpose());
450  return transfer_operator;
451 }
452 
453 #else
454 
455 ContactData::ContactData([[maybe_unused]] const mfem::ParMesh& mesh)
456  : have_lagrange_multipliers_{false}, num_pressure_dofs_{0}
457 {
458 }
459 
461 
462 void ContactData::addContactInteraction([[maybe_unused]] int interaction_id,
463  [[maybe_unused]] const std::set<int>& bdry_attr_surf1,
464  [[maybe_unused]] const std::set<int>& bdry_attr_surf2,
465  [[maybe_unused]] ContactOptions contact_opts)
466 {
467  SLIC_WARNING_ROOT("Smith built without Tribol support. No contact interaction will be added.");
468 }
469 
470 void ContactData::updateGaps([[maybe_unused]] int cycle, [[maybe_unused]] double time, [[maybe_unused]] double& dt,
471  [[maybe_unused]] std::optional<std::reference_wrapper<const mfem::Vector>> u_shape,
472  [[maybe_unused]] std::optional<std::reference_wrapper<const mfem::Vector>> u,
473  [[maybe_unused]] bool eval_jacobian)
474 {
475 }
476 
477 void ContactData::update([[maybe_unused]] int cycle, [[maybe_unused]] double time, [[maybe_unused]] double& dt,
478  [[maybe_unused]] std::optional<std::reference_wrapper<const mfem::Vector>> u_shape,
479  [[maybe_unused]] std::optional<std::reference_wrapper<const mfem::Vector>> u,
480  [[maybe_unused]] std::optional<std::reference_wrapper<const mfem::Vector>> p)
481 {
482 }
483 
485 {
486  FiniteElementDual f(*reference_nodes_->ParFESpace(), "contact force");
487  f = 0.0;
488  return f;
489 }
490 
491 mfem::HypreParVector ContactData::mergedPressures() const { return mfem::HypreParVector(); }
492 
493 mfem::HypreParVector ContactData::mergedGaps([[maybe_unused]] bool zero_inactive) const
494 {
495  return mfem::HypreParVector();
496 }
497 
498 std::unique_ptr<mfem::BlockOperator> ContactData::mergedJacobian() const
499 {
500  jacobian_offsets_ = mfem::Array<int>(
501  {0, reference_nodes_->ParFESpace()->GetTrueVSize(), reference_nodes_->ParFESpace()->GetTrueVSize()});
502  return std::make_unique<mfem::BlockOperator>(jacobian_offsets_);
503 }
504 
505 void ContactData::residualFunction([[maybe_unused]] const mfem::Vector& u_shape, [[maybe_unused]] const mfem::Vector& u,
506  [[maybe_unused]] mfem::Vector& r)
507 {
508 }
509 
510 std::unique_ptr<mfem::BlockOperator> ContactData::jacobianFunction(std::unique_ptr<mfem::HypreParMatrix> orig_J) const
511 {
512  auto J_contact = mergedJacobian();
513  if (J_contact->IsZeroBlock(0, 0)) {
514  J_contact->SetBlock(0, 0, orig_J.release());
515  } else {
516  J_contact->SetBlock(0, 0,
517  mfem::Add(1.0, *orig_J, 1.0, static_cast<mfem::HypreParMatrix&>(J_contact->GetBlock(0, 0))));
518  }
519 
520  return J_contact;
521 }
522 
523 void ContactData::setPressures([[maybe_unused]] const mfem::Vector& true_pressures) const {}
524 
525 void ContactData::setDisplacements([[maybe_unused]] const mfem::Vector& u_shape,
526  [[maybe_unused]] const mfem::Vector& true_displacement)
527 {
528 }
529 
530 std::unique_ptr<mfem::HypreParMatrix> ContactData::contactSubspaceTransferOperator()
531 {
532  std::unique_ptr<mfem::HypreParMatrix> transfer_operator = nullptr;
533  return transfer_operator;
534 }
535 
536 #endif
537 
538 } // namespace smith
void residualFunction(const mfem::Vector &u_shape, const mfem::Vector &u, mfem::Vector &r)
Computes the residual including contact terms.
void update(int cycle, double time, double &dt, std::optional< std::reference_wrapper< const mfem::Vector >> u_shape=std::nullopt, std::optional< std::reference_wrapper< const mfem::Vector >> u=std::nullopt, std::optional< std::reference_wrapper< const mfem::Vector >> p=std::nullopt)
Updates the positions, forces, and Jacobian contributions associated with contact.
mfem::HypreParVector mergedGaps(bool zero_inactive=false) const
Returns nodal gaps from all contact interactions on the contact surface true degrees of freedom.
std::unique_ptr< mfem::BlockOperator > jacobianFunction(std::unique_ptr< mfem::HypreParMatrix > orig_J) const
Computes the Jacobian including contact terms, given the non-contact Jacobian terms.
ContactData(const mfem::ParMesh &mesh)
The constructor.
void setPressures(const mfem::Vector &merged_pressures) const
Set the pressure field.
~ContactData()
Destructor to finalize Tribol.
void updateGaps(int cycle, double time, double &dt, std::optional< std::reference_wrapper< const mfem::Vector >> u_shape=std::nullopt, std::optional< std::reference_wrapper< const mfem::Vector >> u=std::nullopt, bool eval_jacobian=false)
Updates only the gap contributions associated with contact.
mfem::HypreParVector mergedPressures() const
Returns pressures from all contact interactions on the contact surface true degrees of freedom.
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.
std::unique_ptr< mfem::HypreParMatrix > contactSubspaceTransferOperator()
Computes the subspace transfer operator.
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.