GEOS
SolidMechanicsLagrangeContactBubbleStab.hpp
Go to the documentation of this file.
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 
21 #ifndef GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSLAGRANGECONTACTBUBBLESTAB_HPP_
22 #define GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSLAGRANGECONTACTBUBBLESTAB_HPP_
23 
24 #include "physicsSolvers/solidMechanics/contact/ContactSolverBase.hpp"
25 
26 namespace geos
27 {
28 
29 class NumericalMethodsManager;
30 
32 {
33 public:
34 
35  SolidMechanicsLagrangeContactBubbleStab( const string & name,
36  Group * const parent );
37 
39 
44  static string catalogName()
45  {
46  return "SolidMechanicsLagrangeContactBubbleStab";
47  }
51  string getCatalogName() const override { return catalogName(); }
52 
53  virtual void registerDataOnMesh( Group & MeshBodies ) override final;
54 
55  real64 solverStep( real64 const & time_n,
56  real64 const & dt,
57  integer const cycleNumber,
58  DomainPartition & domain ) override final;
59 
60  virtual void
61  setupDofs( DomainPartition const & domain,
62  DofManager & dofManager ) const override;
63 
64  virtual void
66  DofManager & dofManager,
68  ParallelVector & rhs,
69  ParallelVector & solution,
70  bool const setSparsity = true ) override final;
71 
72  virtual void setSparsityPattern( DomainPartition & domain,
73  DofManager & dofManager,
74  CRSMatrix< real64, globalIndex > & localMatrix,
75  SparsityPattern< globalIndex > & pattern ) override final;
76 
77  virtual void
78  implicitStepSetup( real64 const & time_n,
79  real64 const & dt,
80  DomainPartition & domain ) override final;
81 
82  virtual void
83  implicitStepComplete( real64 const & time_n,
84  real64 const & dt,
85  DomainPartition & domain ) override final;
86 
87  virtual void
88  assembleSystem( real64 const time,
89  real64 const dt,
90  DomainPartition & domain,
91  DofManager const & dofManager,
92  CRSMatrixView< real64, globalIndex const > const & localMatrix,
93  arrayView1d< real64 > const & localRhs ) override;
94 
95  virtual real64
97  real64 const & dt,
98  DomainPartition const & domain,
99  DofManager const & dofManager,
100  arrayView1d< real64 const > const & localRhs ) override;
101 
102  virtual void
103  applySystemSolution( DofManager const & dofManager,
104  arrayView1d< real64 const > const & localSolution,
105  real64 const scalingFactor,
106  real64 const dt,
107  DomainPartition & domain ) override;
108 
109  void assembleContact( real64 const dt,
110  DomainPartition & domain,
111  DofManager const & dofManager,
112  CRSMatrixView< real64, globalIndex const > const & localMatrix,
113  arrayView1d< real64 > const & localRhs );
114 
115  void assembleStabilization( real64 const dt,
116  DomainPartition & domain,
117  DofManager const & dofManager,
118  CRSMatrixView< real64, globalIndex const > const & localMatrix,
119  arrayView1d< real64 > const & localRhs );
120 
121  real64 calculateContactResidualNorm( DomainPartition const & domain,
122  DofManager const & dofManager,
123  arrayView1d< real64 const > const & localRhs );
124 
132  template< typename LAMBDA >
133  void forFiniteElementOnFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
134  {
135  std::map< string,
136  array1d< localIndex > > const & faceTypesToFaceElements = m_faceTypesToFaceElements.at( meshName );
137 
138  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
139  {
140  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
141 
142  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
143 
144  lambda( finiteElementName, subRegionFE, faceElemList );
145  }
146  }
147 
155  template< typename LAMBDA >
156  void forFiniteElementOnStickFractureSubRegions( string const & meshName, LAMBDA && lambda ) const
157  {
158  bool const isStickState = true;
159 
160  std::map< string, array1d< localIndex > > const &
161  faceTypesToFaceElements = m_faceTypesToFaceElementsStick.at( meshName );
162 
163  for( const auto & [finiteElementName, faceElementList] : faceTypesToFaceElements )
164  {
165  arrayView1d< localIndex const > const faceElemList = faceElementList.toViewConst();
166 
167  finiteElement::FiniteElementBase const & subRegionFE = *(m_faceTypeToFiniteElements.at( finiteElementName ));
168 
169  lambda( finiteElementName, subRegionFE, faceElemList, isStickState );
170  }
171  }
172 
179  void updateStickSlipList( DomainPartition const & domain );
180 
186  void createFaceTypeList( DomainPartition const & domain );
187 
193  void createBubbleCellList( DomainPartition & domain ) const;
194 
199  void computeRotationMatrices( DomainPartition & domain ) const;
200 
201 
202 private:
210  void addCouplingNumNonzeros( DomainPartition & domain,
211  DofManager & dofManager,
212  arrayView1d< localIndex > const & rowLengths ) const;
213 
220  void addCouplingSparsityPattern( DomainPartition const & domain,
221  DofManager const & dofManager,
222  SparsityPatternView< globalIndex > const & pattern ) const;
223 
225  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElements;
226 
228  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElementsStick;
229 
231  std::map< string, std::map< string, array1d< localIndex > > > m_faceTypesToFaceElementsSlip;
232 
234  std::map< string, std::unique_ptr< geos::finiteElement::FiniteElementBase > > m_faceTypeToFiniteElements;
235 
236  struct viewKeyStruct : ContactSolverBase::viewKeyStruct
237  {
238  constexpr static char const * rotationMatrixString() { return "rotationMatrix"; }
239  };
240 
241 };
242 
243 } /* namespace geos */
244 
245 #endif /* GEOS_PHYSICSSOLVERS_CONTACT_SOLIDMECHANICSLAGRANGECONTACTBUBBLESTAB_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 registerDataOnMesh(Group &MeshBodies) override final
Register wrappers that contain data on the mesh objects.
void createBubbleCellList(DomainPartition &domain) const
Create the list of elements belonging to CellElementSubRegion that are enriched with the bubble basis...
void updateStickSlipList(DomainPartition const &domain)
Create the list of finite elements of the same type for each FaceElementSubRegion (Triangle or Quadri...
void createFaceTypeList(DomainPartition const &domain)
Create the list of finite elements of the same type for each FaceElementSubRegion (Triangle or Quadri...
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 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 implicitStepComplete(real64 const &time_n, real64 const &dt, DomainPartition &domain) override final
perform cleanup for implicit timestep
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
static string catalogName()
name of the node manager in the object catalog
void forFiniteElementOnStickFractureSubRegions(string const &meshName, LAMBDA &&lambda) const
Loop over the finite element type on the stick fracture subregions of meshName and apply callback.
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
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)
virtual real64 calculateResidualNorm(real64 const &time, real64 const &dt, DomainPartition const &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localRhs) override
calculate the norm of the global system residual
void forFiniteElementOnFractureSubRegions(string const &meshName, LAMBDA &&lambda) const
Loop over the finite element type on the fracture subregions of meshName and apply callback.
void computeRotationMatrices(DomainPartition &domain) const
Compute rotation matrices and unit normal vectors for Face elements.
real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override final
entry function to perform a solver step
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override final
function to perform setup for implicit timestep
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
mapBase< TKEY, TVAL, std::integral_constant< bool, true > > map
Ordered map type.
Definition: DataTypes.hpp:339
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.