Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
contact_data.hpp
Go to the documentation of this file.
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 
13 #pragma once
14 
15 #include <functional>
16 #include <memory>
17 #include <optional>
18 #include <set>
19 #include <vector>
20 
21 #include "mfem.hpp"
22 
23 #include "smith/smith_config.hpp"
26 #ifdef SMITH_USE_TRIBOL
28 #endif
29 
30 namespace smith {
31 
32 namespace contact {
33 
39  .enforcement = ContactEnforcement::Penalty,
41  .penalty = 1.0e3};
42 
43 } // namespace contact
44 
49 class ContactData {
50  public:
56  ContactData(const mfem::ParMesh& mesh);
57 
61  ~ContactData();
62 
74  void addContactInteraction(int interaction_id, const std::set<int>& bdry_attr_surf1,
75  const std::set<int>& bdry_attr_surf2, ContactOptions contact_opts);
76 
87  void updateGaps(int cycle, double time, double& dt,
88  std::optional<std::reference_wrapper<const mfem::Vector>> u_shape = std::nullopt,
89  std::optional<std::reference_wrapper<const mfem::Vector>> u = std::nullopt,
90  bool eval_jacobian = false);
91 
102  void update(int cycle, double time, double& dt,
103  std::optional<std::reference_wrapper<const mfem::Vector>> u_shape = std::nullopt,
104  std::optional<std::reference_wrapper<const mfem::Vector>> u = std::nullopt,
105  std::optional<std::reference_wrapper<const mfem::Vector>> p = std::nullopt);
106 
111  void reset();
112 
120  FiniteElementDual forces() const;
121 
131  mfem::HypreParVector mergedPressures() const;
132 
145  mfem::HypreParVector mergedGaps(bool zero_inactive = false) const;
146 
159  std::unique_ptr<mfem::BlockOperator> mergedJacobian() const;
160 
173  void residualFunction(const mfem::Vector& u_shape, const mfem::Vector& u, mfem::Vector& r);
174 
183  std::unique_ptr<mfem::BlockOperator> jacobianFunction(std::unique_ptr<mfem::HypreParMatrix> orig_J) const;
184 
190  std::unique_ptr<mfem::HypreParMatrix> contactSubspaceTransferOperator();
191 
199  {
200 #ifdef SMITH_USE_TRIBOL
201  return !interactions_.empty();
202 #else
203  return false;
204 #endif
205  }
206 
207 #ifdef SMITH_USE_TRIBOL
213  const std::vector<ContactInteraction>& getContactInteractions() const { return interactions_; }
214 #endif
215 
222  bool haveLagrangeMultipliers() const { return have_lagrange_multipliers_; }
223 
229  int numPressureDofs() const { return num_pressure_dofs_; };
230 
231  protected:
243  void setPressures(const mfem::Vector& merged_pressures) const;
244 
251  void setDisplacements(const mfem::Vector& u_shape, const mfem::Vector& u);
252 
253  private:
254 #ifdef SMITH_USE_TRIBOL
260  void updateDofOffsets() const;
261 
265  const mfem::ParMesh& mesh_;
266 #endif
267 
271  const mfem::ParGridFunction* reference_nodes_;
272 
273 #ifdef SMITH_USE_TRIBOL
277  mfem::ParGridFunction current_coords_;
278 
284  std::vector<ContactInteraction> interactions_;
285 #endif
286 
291  bool have_lagrange_multipliers_;
292 
296  int num_pressure_dofs_;
297 
304  mutable bool offsets_up_to_date_;
305 
313  mutable mfem::Array<int> jacobian_offsets_;
314 
323  mutable mfem::Array<int> pressure_dof_offsets_;
324 
333  mutable mfem::Array<HYPRE_BigInt> global_pressure_dof_offsets_;
334 
338  mfem::Array<int> contact_dofs_;
339 
340  int cycle_{0};
341  double time_{0.0};
342  double dt_{1.0};
343 };
344 
345 } // namespace smith
This class stores all ContactInteractions for a problem, calls Tribol functions that act on all conta...
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.
bool haveContactInteractions() const
Have there been contact interactions added?
ContactData(const mfem::ParMesh &mesh)
The constructor.
int numPressureDofs() const
Get the number of Lagrange multiplier true degrees of freedom.
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.
bool haveLagrangeMultipliers() const
Are any contact interactions enforced using Lagrange multipliers?
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.
void reset()
Resets the contact pressures to zero.
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...
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: smith.cpp:36
Stores the options for a contact pair.
ContactMethod method
The contact methodology to be applied.