Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
contact_data.hpp
Go to the documentation of this file.
1 // Copyright (c) 2019-2024, Lawrence Livermore National Security, LLC and
2 // other Serac Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
15 #include "mfem.hpp"
16 
17 #include "serac/serac_config.hpp"
20 #ifdef SERAC_USE_TRIBOL
22 #endif
23 
24 namespace serac {
25 
26 namespace contact {
27 
33  .enforcement = ContactEnforcement::Penalty,
35  .penalty = 1.0e3};
36 
37 } // namespace contact
38 
43 class ContactData {
44 public:
50  ContactData(const mfem::ParMesh& mesh);
51 
55  ~ContactData();
56 
65  void addContactInteraction(int interaction_id, const std::set<int>& bdry_attr_surf1,
66  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts);
67 
75  void update(int cycle, double time, double& dt);
76 
82  FiniteElementDual forces() const;
83 
93  mfem::HypreParVector mergedPressures() const;
94 
105  mfem::HypreParVector mergedGaps(bool zero_inactive = false) const;
106 
119  std::unique_ptr<mfem::BlockOperator> mergedJacobian() const;
120 
128  void residualFunction(const mfem::Vector& u, mfem::Vector& r);
129 
137  std::unique_ptr<mfem::BlockOperator> jacobianFunction(const mfem::Vector& u, mfem::HypreParMatrix* orig_J) const;
138 
150  void setPressures(const mfem::Vector& merged_pressures) const;
151 
157  void setDisplacements(const mfem::Vector& u);
158 
166  {
167 #ifdef SERAC_USE_TRIBOL
168  return !interactions_.empty();
169 #else
170  return false;
171 #endif
172  }
173 
180  bool haveLagrangeMultipliers() const { return have_lagrange_multipliers_; }
181 
187  int numPressureDofs() const { return num_pressure_dofs_; };
188 
189 private:
190 #ifdef SERAC_USE_TRIBOL
196  void updateDofOffsets() const;
197 
201  const mfem::ParMesh& mesh_;
202 #endif
203 
207  const mfem::ParGridFunction* reference_nodes_;
208 
209 #ifdef SERAC_USE_TRIBOL
213  mfem::ParGridFunction current_coords_;
214 
218  std::vector<ContactInteraction> interactions_;
219 #endif
220 
225  bool have_lagrange_multipliers_;
226 
230  int num_pressure_dofs_;
231 
238  mutable bool offsets_up_to_date_;
239 
247  mutable mfem::Array<int> jacobian_offsets_;
248 
257  mutable mfem::Array<int> pressure_dof_offsets_;
258 
267  mutable mfem::Array<HYPRE_BigInt> global_pressure_dof_offsets_;
268 };
269 
270 } // namespace serac
This class stores all ContactInteractions for a problem, calls Tribol functions that act on all conta...
void setDisplacements(const mfem::Vector &u)
Update the current coordinates based on the new displacement field.
void setPressures(const mfem::Vector &merged_pressures) const
Set the pressure field.
mfem::HypreParVector mergedPressures() const
Returns pressures from all contact interactions on the contact surface true degrees of freedom.
int numPressureDofs() const
Get the number of Lagrange multiplier true degrees of freedom.
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.
~ContactData()
Destructor to finalize Tribol.
ContactData(const mfem::ParMesh &mesh)
The constructor.
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.
std::unique_ptr< mfem::BlockOperator > jacobianFunction(const mfem::Vector &u, mfem::HypreParMatrix *orig_J) const
Computes the Jacobian including contact terms, given the non-contact Jacobian terms.
bool haveLagrangeMultipliers() const
Are any contact interactions enforced using Lagrange multipliers?
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 residualFunction(const mfem::Vector &u, mfem::Vector &r)
Computes the residual including contact terms.
bool haveContactInteractions() const
Have there been contact interactions added?
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
This file contains enumerations and record types for contact configuration.
const ContactOptions default_contact_options
Default contact options: frictionless mortar with penalty = 1000 enforcement.
Class for storing a contact interaction and interfacing with Tribol.
This contains a class that represents the dual of a finite element vector space, i....
Accelerator functionality.
Definition: serac.cpp:38
Stores the options for a contact pair.
ContactMethod method
The contact methodology to be applied.