GEOS
SolidMechanicsAugmentedLagrangianContact.hpp
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
16 /*
17  * SolidMechanicsAugmentedLagrangianContact.hpp
18  *
19  */
20 
21 #ifndef GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSAUGMENTEDLAGRANGIANCONTACT_HPP_
22 #define GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSAUGMENTEDLAGRANGIANCONTACT_HPP_
23 
24 #include "physicsSolvers/contact/ContactSolverBase.hpp"
25 
26 namespace geos
27 {
28 
30 {
31 public:
32  SolidMechanicsAugmentedLagrangianContact( const string & name,
33  Group * const parent );
34 
36 
41  static string catalogName()
42  {
43  return "SolidMechanicsAugmentedLagrangianContact";
44  }
48  string getCatalogName() const override { return catalogName(); }
49 
50  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override final;
51 
52  virtual void setupDofs( DomainPartition const & domain,
53  DofManager & dofManager ) const override;
54 
55  virtual void setupSystem( DomainPartition & domain,
56  DofManager & dofManager,
58  ParallelVector & rhs,
59  ParallelVector & solution,
60  bool const setSparsity = true ) override final;
61 
62  virtual void implicitStepSetup( real64 const & time_n,
63  real64 const & dt,
64  DomainPartition & domain ) override final;
65 
66  virtual void implicitStepComplete( real64 const & time_n,
67  real64 const & dt,
68  DomainPartition & domain ) override final;
69 
70  virtual void assembleSystem( real64 const time,
71  real64 const dt,
72  DomainPartition & domain,
73  DofManager const & dofManager,
74  CRSMatrixView< real64, globalIndex const > const & localMatrix,
75  arrayView1d< real64 > const & localRhs ) override;
76 
77  virtual real64 calculateResidualNorm( real64 const & time_n,
78  real64 const & dt,
79  DomainPartition const & domain,
80  DofManager const & dofManager,
81  arrayView1d< real64 const > const & localRhs ) override;
82 
83  virtual void applySystemSolution( DofManager const & dofManager,
84  arrayView1d< real64 const > const & localSolution,
85  real64 const scalingFactor,
86  real64 const dt,
87  DomainPartition & domain ) override;
88 
89  void updateState( DomainPartition & domain ) override final;
90 
91  virtual bool updateConfiguration( DomainPartition & domain ) override final;
92 
93 
101  template< typename LAMBDA >
102  void forFiniteElementOnFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
103  {
104 
105  std::map< string,
106  array1d< localIndex > > const & faceTypesToFaceElements = m_faceTypesToFaceElements.at( meshName );
107 
108  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
109  {
110  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
111 
112  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
113 
114  lambda( finiteElementName, subRegionFE, faceElemList );
115  }
116 
117  }
118 
126  template< typename LAMBDA >
127  void forFiniteElementOnStickFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
128  {
129 
130  bool const isStickState = true;
131 
132  std::map< string, array1d< localIndex > > const &
133  faceTypesToFaceElements = m_faceTypesToFaceElementsStick.at( meshName );
134 
135  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
136  {
137  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
138 
139  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
140 
141  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
142  }
143 
144  }
145 
153  template< typename LAMBDA >
154  void forFiniteElementOnSlipFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
155  {
156 
157  bool const isStickState = false;
158 
159  std::map< string, array1d< localIndex > > const &
160  faceTypesToFaceElements = m_faceTypesToFaceElementsSlip.at( meshName );
161 
162  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
163  {
164  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
165 
166  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
167 
168  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
169  }
170 
171  }
172 
179  void updateStickSlipList( DomainPartition const & domain );
180 
186  void createFaceTypeList( DomainPartition const & domain );
187 
193  void createBubbleCellList( DomainPartition & domain ) const;
194 
195 private:
196 
204  void addCouplingNumNonzeros( DomainPartition & domain,
205  DofManager & dofManager,
206  arrayView1d< localIndex > const & rowLengths ) const;
207 
214  void addCouplingSparsityPattern( DomainPartition const & domain,
215  DofManager const & dofManager,
216  SparsityPatternView< globalIndex > const & pattern ) const;
217 
218  void computeTolerances( DomainPartition & domain ) const;
219 
221  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElements;
222 
224  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElementsStick;
225 
227  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElementsSlip;
228 
230  std::map< string, std::unique_ptr< geos::finiteElement::FiniteElementBase > > m_faceTypeToFiniteElements;
231 
232  struct viewKeyStruct : ContactSolverBase::viewKeyStruct
233  {
234 
235  constexpr static char const * normalDisplacementToleranceString() { return "normalDisplacementTolerance"; }
236 
237  constexpr static char const * normalTractionToleranceString() { return "normalTractionTolerance"; }
238 
239  constexpr static char const * slidingToleranceString() { return "slidingTolerance"; }
240 
241  constexpr static char const * dispJumpUpdPenaltyString() { return "dispJumpUpdPenalty"; }
242  };
243 
246  real64 const m_slidingCheckTolerance = 0.05;
247 
249  bool m_simultaneous = true;
250 
253  bool m_symmetric = true;
254 
255 };
256 
257 } /* namespace geos */
258 
259 #endif /* GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSAUGMENTEDLAGRANGIANCONTACT_HPP_ */
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:44
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
virtual void setupSystem(DomainPartition &domain, DofManager &dofManager, CRSMatrix< real64, globalIndex > &localMatrix, ParallelVector &rhs, ParallelVector &solution, bool const setSparsity=true) override final
Set up the linear system (DOF indices and sparsity patterns)
static string catalogName()
name of the node manager in the object catalog
void createFaceTypeList(DomainPartition const &domain)
Create the list of finite elements of the same type for each FaceElementSubRegion (Triangle or Quadri...
void createBubbleCellList(DomainPartition &domain) const
Create the list of elements belonging to CellElementSubRegion that are enriched with the bubble basis...
virtual real64 calculateResidualNorm(real64 const &time_n, real64 const &dt, DomainPartition const &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localRhs) override
calculate the norm of the global system residual
virtual void assembleSystem(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
function to assemble the linear system matrix and rhs
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void implicitStepComplete(real64 const &time_n, real64 const &dt, DomainPartition &domain) override final
perform cleanup for implicit timestep
virtual void registerDataOnMesh(dataRepository::Group &meshBodies) override final
Register wrappers that contain data on the mesh objects.
virtual bool updateConfiguration(DomainPartition &domain) override final
updates the configuration (if needed) based on the state after a converged Newton loop.
void forFiniteElementOnFractureSubRegions(string const &meshName, LAMBDA &&lambda) const
Loop over the finite element type on the fracture subregions of meshName and apply callback.
void forFiniteElementOnSlipFractureSubRegions(string const &meshName, LAMBDA &&lambda) const
Loop over the finite element type on the slip fracture subregions of meshName and apply callback.
void updateState(DomainPartition &domain) override final
Recompute all dependent quantities from primary variables (including constitutive models)
virtual void applySystemSolution(DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor, real64 const dt, DomainPartition &domain) override
Function to apply the solution vector to the state.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override final
function to perform setup for implicit timestep
void updateStickSlipList(DomainPartition const &domain)
Create the list of finite elements of the same type for each FaceElementSubRegion (Triangle or Quadri...
void forFiniteElementOnStickFractureSubRegions(string const &meshName, LAMBDA &&lambda) const
Loop over the finite element type on the stick fracture subregions of meshName and apply callback.
Base class for FEM element implementations.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
std::string string
String type.
Definition: DataTypes.hpp:91
LvArray::SparsityPatternView< COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > SparsityPatternView
Alias for Sparsity pattern View.
Definition: DataTypes.hpp:302
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
LvArray::CRSMatrix< T, COL_INDEX, localIndex, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
Definition: DataTypes.hpp:306
mapBase< TKEY, TVAL, std::integral_constant< bool, true > > map
Ordered map type.
Definition: DataTypes.hpp:369
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.