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 TotalEnergies
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/solidMechanics/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 setSparsityPattern( DomainPartition & domain,
63  DofManager & dofManager,
64  CRSMatrix< real64, globalIndex > & localMatrix,
65  SparsityPattern< globalIndex > & pattern ) override final;
66 
67  virtual void implicitStepSetup( real64 const & time_n,
68  real64 const & dt,
69  DomainPartition & domain ) override final;
70 
71  virtual void implicitStepComplete( real64 const & time_n,
72  real64 const & dt,
73  DomainPartition & domain ) override final;
74 
75  virtual void assembleSystem( real64 const time,
76  real64 const dt,
77  DomainPartition & domain,
78  DofManager const & dofManager,
79  CRSMatrixView< real64, globalIndex const > const & localMatrix,
80  arrayView1d< real64 > const & localRhs ) override;
81 
82  void assembleContact( real64 const time,
83  real64 const dt,
84  DomainPartition & domain,
85  DofManager const & dofManager,
86  CRSMatrixView< real64, globalIndex const > const & localMatrix,
87  arrayView1d< real64 > const & localRhs );
88 
89  void assembleForceResidualPressureContribution( DomainPartition & domain,
90  real64 const & dt,
91  DofManager const & dofManager,
92  CRSMatrixView< real64, globalIndex const > const & localMatrix,
93  arrayView1d< real64 > const & localRhs );
94 
95  virtual real64 calculateResidualNorm( real64 const & time_n,
96  real64 const & dt,
97  DomainPartition const & domain,
98  DofManager const & dofManager,
99  arrayView1d< real64 const > const & localRhs ) override;
100 
101  virtual void applySystemSolution( DofManager const & dofManager,
102  arrayView1d< real64 const > const & localSolution,
103  real64 const scalingFactor,
104  real64 const dt,
105  DomainPartition & domain ) override;
106 
107  void updateState( DomainPartition & domain ) override final;
108 
109  virtual bool updateConfiguration( DomainPartition & domain,
110  integer configurationLoopIter ) override final;
111 
112 
120  template< typename LAMBDA >
121  void forFiniteElementOnFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
122  {
123 
124  stdMap< string,
125  array1d< localIndex > > const & faceTypesToFaceElements = m_faceTypesToFaceElements.at( meshName );
126 
127  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
128  {
129  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
130 
131  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
132 
133  lambda( finiteElementName, subRegionFE, faceElemList );
134  }
135 
136  }
137 
145  template< typename LAMBDA >
146  void forFiniteElementOnStickFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
147  {
148 
149  bool const isStickState = true;
150 
152  faceTypesToFaceElements = m_faceTypesToFaceElementsStick.at( meshName );
153 
154  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
155  {
156  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
157 
158  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
159 
160  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
161  }
162 
163  }
164 
172  template< typename LAMBDA >
173  void forFiniteElementOnSlipFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
174  {
175 
176  bool const isStickState = false;
177 
179  faceTypesToFaceElements = m_faceTypesToFaceElementsSlip.at( meshName );
180 
181  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
182  {
183  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
184 
185  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
186 
187  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
188  }
189 
190  }
191 
198  void updateStickSlipList( DomainPartition const & domain );
199 
205  void createFaceTypeList( DomainPartition const & domain );
206 
212  void createBubbleCellList( DomainPartition & domain ) const;
213 
214 private:
215 
223  void addCouplingNumNonzeros( DomainPartition & domain,
224  DofManager & dofManager,
225  arrayView1d< localIndex > const & rowLengths ) const;
226 
233  void addCouplingSparsityPattern( DomainPartition const & domain,
234  DofManager const & dofManager,
235  SparsityPatternView< globalIndex > const & pattern ) const;
236 
237  void computeTolerances( DomainPartition & domain ) const;
238 
240  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElements;
241 
243  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElementsStick;
244 
246  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElementsSlip;
247 
250 
251  struct viewKeyStruct : ContactSolverBase::viewKeyStruct
252  {
253 
254  constexpr static char const * normalDisplacementToleranceString() { return "normalDisplacementTolerance"; }
255 
256  constexpr static char const * normalTractionToleranceString() { return "normalTractionTolerance"; }
257 
258  constexpr static char const * slidingToleranceString() { return "slidingTolerance"; }
259 
260  constexpr static char const * dispJumpUpdPenaltyString() { return "dispJumpUpdPenalty"; }
261 
262  constexpr static char const * simultaneousString() { return "simultaneous"; }
263 
264  constexpr static char const * symmetricString() { return "symmetric"; }
265 
266  constexpr static char const * iterativePenaltyNFacString() { return "iterPenaltyN"; }
267 
268  constexpr static char const * iterativePenaltyTFacString() { return "iterPenaltyT"; }
269 
270  constexpr static char const * tolJumpDispNFacString() { return "tolJumpN"; }
271 
272  constexpr static char const * tolJumpDispTFacString() { return "tolJumpT"; }
273 
274  constexpr static char const * tolNormalTracFacString() { return "tolNormalTrac"; }
275 
276  constexpr static char const * tolTauLimitString() { return "tolTauLimit"; }
277 
278  };
279 
282  real64 m_slidingCheckTolerance = 5.e-02;
283 
285  int m_simultaneous = 1;
286 
289  int m_symmetric = 1;
290 
292  real64 m_iterPenaltyNFac = 10.0;
293 
295  real64 m_iterPenaltyTFac = 0.1;
296 
298  real64 m_tolJumpDispNFac = 1.e-07;
299 
301  real64 m_tolJumpDispTFac = 1.e-05;
302 
304  real64 m_tolNormalTracFac = 0.5;
305 
306 };
307 
308 } /* namespace geos */
309 
310 #endif /* GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSAUGMENTEDLAGRANGIANCONTACT_HPP_ */
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:45
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 void setSparsityPattern(DomainPartition &domain, DofManager &dofManager, CRSMatrix< real64, globalIndex > &localMatrix, SparsityPattern< globalIndex > &pattern) override final
Set the sparsity pattern of the linear system matrix.
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, integer configurationLoopIter) 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:179
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
LvArray::CRSMatrix< T, COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
Definition: DataTypes.hpp:305
std::string string
String type.
Definition: DataTypes.hpp:90
LvArray::SparsityPatternView< COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > SparsityPatternView
Alias for Sparsity pattern View.
Definition: DataTypes.hpp:301
LvArray::SparsityPattern< COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > SparsityPattern
Alias for Sparsity pattern class.
Definition: DataTypes.hpp:297
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
internal::StdMapWrapper< std::map< Key, T, Compare, Allocator >, USE_STD_CONTAINER_BOUNDS_CHECKING > stdMap
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.