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 postInputInitialization() override;
51 
52  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override final;
53 
55 
56  virtual void setupDofs( DomainPartition const & domain,
57  DofManager & dofManager ) const override;
58 
59  virtual void setupSystem( DomainPartition & domain,
60  DofManager & dofManager,
62  ParallelVector & rhs,
63  ParallelVector & solution,
64  bool const setSparsity = true ) override final;
65 
66  virtual void setSparsityPattern( DomainPartition & domain,
67  DofManager & dofManager,
68  CRSMatrix< real64, globalIndex > & localMatrix,
69  SparsityPattern< globalIndex > & pattern ) override final;
70 
71  virtual void implicitStepSetup( real64 const & time_n,
72  real64 const & dt,
73  DomainPartition & domain ) override final;
74 
75  virtual void implicitStepComplete( real64 const & time_n,
76  real64 const & dt,
77  DomainPartition & domain ) override final;
78 
79  virtual void assembleSystem( real64 const time,
80  real64 const dt,
81  DomainPartition & domain,
82  DofManager const & dofManager,
83  CRSMatrixView< real64, globalIndex const > const & localMatrix,
84  arrayView1d< real64 > const & localRhs ) override;
85 
86  void assembleContact( real64 const time,
87  real64 const dt,
88  DomainPartition & domain,
89  DofManager const & dofManager,
90  CRSMatrixView< real64, globalIndex const > const & localMatrix,
91  arrayView1d< real64 > const & localRhs );
92 
93  void assembleForceResidualPressureContribution( DomainPartition & domain,
94  real64 const & dt,
95  DofManager const & dofManager,
96  CRSMatrixView< real64, globalIndex const > const & localMatrix,
97  arrayView1d< real64 > const & localRhs );
98 
99  virtual real64 calculateResidualNorm( real64 const & time_n,
100  real64 const & dt,
101  DomainPartition const & domain,
102  DofManager const & dofManager,
103  arrayView1d< real64 const > const & localRhs ) override;
104 
105  virtual void applySystemSolution( DofManager const & dofManager,
106  arrayView1d< real64 const > const & localSolution,
107  real64 const scalingFactor,
108  real64 const dt,
109  DomainPartition & domain ) override;
110 
111  void updateState( DomainPartition & domain ) override final;
112 
113  virtual bool updateConfiguration( DomainPartition & domain,
114  integer configurationLoopIter ) override final;
115 
116 
124  template< typename LAMBDA >
125  void forFiniteElementOnFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
126  {
127 
128  stdMap< string,
129  array1d< localIndex > > const & faceTypesToFaceElements = m_faceTypesToFaceElements.at( meshName );
130 
131  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
132  {
133  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
134 
135  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
136 
137  lambda( finiteElementName, subRegionFE, faceElemList );
138  }
139 
140  }
141 
149  template< typename LAMBDA >
150  void forFiniteElementOnStickFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
151  {
152 
153  bool const isStickState = true;
154 
156  faceTypesToFaceElements = m_faceTypesToFaceElementsStick.at( meshName );
157 
158  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
159  {
160  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
161 
162  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
163 
164  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
165  }
166 
167  }
168 
176  template< typename LAMBDA >
177  void forFiniteElementOnSlipFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
178  {
179 
180  bool const isStickState = false;
181 
183  faceTypesToFaceElements = m_faceTypesToFaceElementsSlip.at( meshName );
184 
185  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
186  {
187  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
188 
189  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
190 
191  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
192  }
193 
194  }
195 
202  void updateStickSlipList( DomainPartition const & domain );
203 
209  void createFaceTypeList( DomainPartition const & domain );
210 
216  void createBubbleCellList( DomainPartition & domain ) const;
217 
218 private:
219 
224  void validateTetrahedralQuadrature( Group & meshBodies );
225 
233  void addCouplingNumNonzeros( DomainPartition & domain,
234  DofManager & dofManager,
235  arrayView1d< localIndex > const & rowLengths ) const;
236 
243  void addCouplingSparsityPattern( DomainPartition const & domain,
244  DofManager const & dofManager,
245  SparsityPatternView< globalIndex > const & pattern ) const;
246 
247  void computeTolerances( DomainPartition & domain ) const;
248 
250  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElements;
251 
253  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElementsStick;
254 
256  stdMap< string, stdMap< string, array1d< localIndex > > > m_faceTypesToFaceElementsSlip;
257 
260 
261  struct viewKeyStruct : ContactSolverBase::viewKeyStruct
262  {
263 
264  constexpr static char const * normalDisplacementToleranceString() { return "normalDisplacementTolerance"; }
265 
266  constexpr static char const * normalTractionToleranceString() { return "normalTractionTolerance"; }
267 
268  constexpr static char const * slidingToleranceString() { return "slidingTolerance"; }
269 
270  constexpr static char const * dispJumpUpdPenaltyString() { return "dispJumpUpdPenalty"; }
271 
272  constexpr static char const * simultaneousString() { return "simultaneous"; }
273 
274  constexpr static char const * symmetricString() { return "symmetric"; }
275 
276  constexpr static char const * iterativePenaltyNFacString() { return "iterPenaltyN"; }
277 
278  constexpr static char const * iterativePenaltyTFacString() { return "iterPenaltyT"; }
279 
280  constexpr static char const * tolJumpDispNFacString() { return "tolJumpN"; }
281 
282  constexpr static char const * tolJumpDispTFacString() { return "tolJumpT"; }
283 
284  constexpr static char const * tolNormalTracFacString() { return "tolNormalTrac"; }
285 
286  constexpr static char const * tolTauLimitString() { return "tolTauLimit"; }
287 
288  };
289 
292  real64 m_slidingCheckTolerance = 5.e-02;
293 
295  int m_simultaneous = 1;
296 
299  int m_symmetric = 1;
300 
302  real64 m_iterPenaltyNFac = 10.0;
303 
305  real64 m_iterPenaltyTFac = 0.1;
306 
308  real64 m_tolJumpDispNFac = 1.e-07;
309 
311  real64 m_tolJumpDispTFac = 1.e-05;
312 
314  real64 m_tolNormalTracFac = 0.5;
315 
316 };
317 
318 } /* namespace geos */
319 
320 #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 initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
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.