GEOS
SolidMechanicsLagrangianFEM.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 
20 #ifndef GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_SOLIDMECHANICSLAGRANGIANFEM_HPP_
21 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_SOLIDMECHANICSLAGRANGIANFEM_HPP_
22 
24 #include "common/TimingMacros.hpp"
26 #include "kernels/StrainHelper.hpp"
27 #include "mesh/mpiCommunications/CommunicationTools.hpp"
28 #include "mesh/mpiCommunications/MPI_iCommData.hpp"
30 
32 
33 namespace geos
34 {
35 
42 {
43 public:
44 
46  static string coupledSolverAttributePrefix() { return "solid"; }
47 
54  {
55  QuasiStatic,
58  };
59 
65  SolidMechanicsLagrangianFEM( const string & name,
66  Group * const parent );
67 
71  static string catalogName() { return "SolidMechanicsLagrangianFEM"; }
75  string getCatalogName() const override { return catalogName(); }
76 
77  virtual void initializePreSubGroups() override;
78 
79  virtual void registerDataOnMesh( Group & meshBodies ) override;
80 
87  virtual
88  real64 solverStep( real64 const & time_n,
89  real64 const & dt,
90  integer const cycleNumber,
91  DomainPartition & domain ) override;
92 
93  virtual
94  real64 explicitStep( real64 const & time_n,
95  real64 const & dt,
96  integer const cycleNumber,
97  DomainPartition & domain ) override;
98 
99  virtual void
100  implicitStepSetup( real64 const & time_n,
101  real64 const & dt,
102  DomainPartition & domain ) override;
103 
104  virtual void
105  setupDofs( DomainPartition const & domain,
106  DofManager & dofManager ) const override;
107 
108  virtual void
110  DofManager & dofManager,
111  CRSMatrix< real64, globalIndex > & localMatrix,
112  ParallelVector & rhs,
113  ParallelVector & solution,
114  bool const setSparsity = false ) override;
115 
116  virtual std::unique_ptr< PreconditionerBase< LAInterface > >
117  createPreconditioner( DomainPartition & domain ) const override;
118 
119  virtual void
120  assembleSystem( real64 const time,
121  real64 const dt,
122  DomainPartition & domain,
123  DofManager const & dofManager,
124  CRSMatrixView< real64, globalIndex const > const & localMatrix,
125  arrayView1d< real64 > const & localRhs ) override;
126 
127  virtual void
128  applySystemSolution( DofManager const & dofManager,
129  arrayView1d< real64 const > const & localSolution,
130  real64 const scalingFactor,
131  real64 const dt,
132  DomainPartition & domain ) override;
133 
134  virtual void updateState( DomainPartition & domain ) override
135  {
136  // There should be nothing to update
137  GEOS_UNUSED_VAR( domain );
138  };
139 
140  virtual void applyBoundaryConditions( real64 const time,
141  real64 const dt,
142  DomainPartition & domain,
143  DofManager const & dofManager,
144  CRSMatrixView< real64, globalIndex const > const & localMatrix,
145  arrayView1d< real64 > const & localRhs ) override;
146 
147  virtual real64
148  calculateResidualNorm( real64 const & time_n,
149  real64 const & dt,
150  DomainPartition const & domain,
151  DofManager const & dofManager,
152  arrayView1d< real64 const > const & localRhs ) override;
153 
154  virtual void resetStateToBeginningOfStep( DomainPartition & domain ) override;
155 
156  virtual void implicitStepComplete( real64 const & time,
157  real64 const & dt,
158  DomainPartition & domain ) override;
159 
163  template< typename TYPE_LIST,
164  typename KERNEL_WRAPPER,
165  typename ... PARAMS >
166  real64 assemblyLaunch( MeshLevel & mesh,
167  DofManager const & dofManager,
168  string_array const & regionNames,
169  string const & materialNamesString,
170  CRSMatrixView< real64, globalIndex const > const & localMatrix,
171  arrayView1d< real64 > const & localRhs,
172  real64 const dt,
173  PARAMS && ... params );
174 
175  real64 explicitKernelDispatch( MeshLevel & mesh,
176  string_array const & targetRegions,
177  string const & finiteElementName,
178  real64 const dt,
179  std::string const & elementListName );
180 
191  DofManager const & dofManager,
192  DomainPartition & domain,
193  CRSMatrixView< real64, globalIndex const > const & localMatrix,
194  arrayView1d< real64 > const & localRhs );
195 
196  void applyTractionBC( real64 const time,
197  DofManager const & dofManager,
198  DomainPartition & domain,
199  arrayView1d< real64 > const & localRhs );
200 
201  void applyChomboPressure( DofManager const & dofManager,
202  DomainPartition & domain,
203  arrayView1d< real64 > const & localRhs );
204 
205 
206  void applyContactConstraint( DofManager const & dofManager,
207  DomainPartition & domain,
208  CRSMatrixView< real64, globalIndex const > const & localMatrix,
209  arrayView1d< real64 > const & localRhs );
210 
211  virtual real64
213  DofManager const & dofManager,
214  arrayView1d< real64 const > const & localSolution ) override;
215 
216  void enableFixedStressPoromechanicsUpdate();
217 
218  virtual void saveSequentialIterationState( DomainPartition & domain ) override;
219 
221  {
222  static constexpr char const * newmarkGammaString() { return "newmarkGamma"; }
223  static constexpr char const * newmarkBetaString() { return "newmarkBeta"; }
224  static constexpr char const * massDampingString() { return "massDamping"; }
225  static constexpr char const * stiffnessDampingString() { return "stiffnessDamping"; }
226  static constexpr char const * timeIntegrationOptionString() { return "timeIntegrationOption"; }
227  static constexpr char const * maxNumResolvesString() { return "maxNumResolves"; }
228  static constexpr char const * strainTheoryString() { return "strainTheory"; }
229  static constexpr char const * solidMaterialNamesString() { return "solidMaterialNames"; }
230  static constexpr char const * contactRelationNameString() { return "contactRelationName"; }
231  static constexpr char const * noContactRelationNameString() { return "NOCONTACT"; }
232  static constexpr char const * maxForceString() { return "maxForce"; }
233  static constexpr char const * elemsAttachedToSendOrReceiveNodesString() { return "elemsAttachedToSendOrReceiveNodes"; }
234  static constexpr char const * elemsNotAttachedToSendOrReceiveNodesString() { return "elemsNotAttachedToSendOrReceiveNodes"; }
235  static constexpr char const * surfaceGeneratorNameString() { return "surfaceGeneratorName"; }
236 
237  static constexpr char const * sendOrReceiveNodesString() { return "sendOrReceiveNodes";}
238  static constexpr char const * nonSendOrReceiveNodesString() { return "nonSendOrReceiveNodes";}
239  static constexpr char const * targetNodesString() { return "targetNodes";}
240  static constexpr char const * forceString() { return "Force";}
241 
242  static constexpr char const * contactPenaltyStiffnessString() { return "contactPenaltyStiffness"; }
243 
244  };
245 
246  SortedArray< localIndex > & getElemsAttachedToSendOrReceiveNodes( ElementSubRegionBase & subRegion )
247  {
248  return subRegion.getReference< SortedArray< localIndex > >( viewKeyStruct::elemsAttachedToSendOrReceiveNodesString() );
249  }
250 
251  SortedArray< localIndex > & getElemsNotAttachedToSendOrReceiveNodes( ElementSubRegionBase & subRegion )
252  {
253  return subRegion.getReference< SortedArray< localIndex > >( viewKeyStruct::elemsNotAttachedToSendOrReceiveNodesString() );
254  }
255 
256  real64 & getMaxForce() { return m_maxForce; }
257  real64 const & getMaxForce() const { return m_maxForce; }
258 
259  void computeRigidBodyModes( DomainPartition & domain ) const;
260 
261  arrayView1d< ParallelVector > const & getRigidBodyModes( DomainPartition & domain ) const
262  {
263  computeRigidBodyModes( domain );
264  return m_rigidBodyModes;
265  }
266 
267  /*
268  * @brief Utility function to set the stress initialization flag
269  * @param[in] performStressInitialization true if the solver has to initialize stress, false otherwise
270  */
271  void setStressInitialization( bool const performStressInitialization )
272  {
273  m_performStressInitialization = performStressInitialization;
274  }
275 
276  TimeIntegrationOption timeIntegrationOption() const { return m_timeIntegrationOption; }
277 
278 protected:
279  virtual void postInputInitialization() override;
280 
282 
283  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override;
284 
285  real64 m_newmarkGamma;
286  real64 m_newmarkBeta;
287  real64 m_massDamping;
288  real64 m_stiffnessDamping;
289  TimeIntegrationOption m_timeIntegrationOption;
290  real64 m_maxForce = 0.0;
291  integer m_maxNumResolves;
292  integer m_strainTheory;
293 
298 
301 
302  real64 m_contactPenaltyStiffness;
303 
304 private:
305 
306  string m_contactRelationName;
307 
308  PhysicsSolverBase *m_surfaceGenerator;
309  string m_surfaceGeneratorName;
310 };
311 
313  "QuasiStatic",
314  "ImplicitDynamic",
315  "ExplicitDynamic" );
316 
317 //**********************************************************************************************************************
318 //**********************************************************************************************************************
319 //**********************************************************************************************************************
320 
321 
322 template< typename TYPE_LIST,
323  typename KERNEL_WRAPPER,
324  typename ... PARAMS >
325 real64 SolidMechanicsLagrangianFEM::assemblyLaunch( MeshLevel & mesh,
326  DofManager const & dofManager,
327  string_array const & regionNames,
328  string const & materialNamesString,
329  CRSMatrixView< real64, globalIndex const > const & localMatrix,
330  arrayView1d< real64 > const & localRhs,
331  real64 const dt,
332  PARAMS && ... params )
333 {
335 
336  NodeManager const & nodeManager = mesh.getNodeManager();
337 
338  string const dofKey = dofManager.getKey( fields::solidMechanics::totalDisplacement::key() );
339  arrayView1d< globalIndex const > const & dofNumber = nodeManager.getReference< globalIndex_array >( dofKey );
340 
341  real64 const gravityVectorData[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( gravityVector() );
342 
343  KERNEL_WRAPPER kernelWrapper( dofNumber,
344  dofManager.rankOffset(),
345  localMatrix,
346  localRhs,
347  dt,
348  gravityVectorData,
349  std::forward< PARAMS >( params )... );
350 
351  return finiteElement::
352  regionBasedKernelApplication< parallelDevicePolicy< >,
353  TYPE_LIST >( mesh,
354  regionNames,
355  this->getDiscretizationName(),
356  materialNamesString,
357  kernelWrapper );
358 
359 }
360 
361 } /* namespace geos */
362 
363 #endif /* GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_SOLIDMECHANICSLAGRANGIANFEM_HPP_ */
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:45
globalIndex rankOffset(string const &fieldName) const
string const & getKey(string const &fieldName) const
Return the key used to record the field in the DofManager.
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
Class facilitating the representation of a multi-level discretization of a MeshBody.
Definition: MeshLevel.hpp:42
NodeManager const & getNodeManager() const
Get the node manager.
Definition: MeshLevel.hpp:155
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
Base class for all physics solvers.
string getDiscretizationName() const
return the name of the discretization object
R1Tensor const gravityVector() const
return the value of the gravity vector specified in PhysicsSolverManager
bool m_performStressInitialization
Flag to indicate that the solver is going to perform stress initialization.
virtual void registerDataOnMesh(Group &meshBodies) override
Register wrappers that contain data on the mesh objects.
bool m_isFixedStressPoromechanicsUpdate
Flag to indicate that the solver is running in fixed stress (sequential) mode.
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
virtual void saveSequentialIterationState(DomainPartition &domain) override
Save the state of the solver for sequential iteration.
SolidMechanicsLagrangianFEM(const string &name, Group *const parent)
virtual void initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
void applyDisplacementBCImplicit(real64 const time, DofManager const &dofManager, DomainPartition &domain, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
array1d< ParallelVector > m_rigidBodyModes
Rigid body modes; TODO remove mutable hack.
virtual void setConstitutiveNamesCallSuper(ElementSubRegionBase &subRegion) const override
This function sets constitutive name fields on an ElementSubRegionBase, and calls the base function i...
virtual real64 scalingForSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution) override
Function to determine if the solution vector should be scaled back in order to maintain a known const...
virtual void postInputInitialization() override
static string coupledSolverAttributePrefix()
String used to form the solverName used to register single-physics solvers in CoupledSolver.
GEOS_DECLTYPE_AUTO_RETURN getReference(LOOKUP_TYPE const &lookup) const
Look up a wrapper and get reference to wrapped object.
Definition: Group.hpp:1275
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 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 implicitStepComplete(real64 const &time, real64 const &dt, DomainPartition &domain) override
perform cleanup for implicit timestep
virtual real64 explicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
Entry function for an explicit time integration step.
virtual std::unique_ptr< PreconditionerBase< LAInterface > > createPreconditioner(DomainPartition &domain) const override
Create a preconditioner for this solver's linear system.
virtual void setupSystem(DomainPartition &domain, DofManager &dofManager, CRSMatrix< real64, globalIndex > &localMatrix, ParallelVector &rhs, ParallelVector &solution, bool const setSparsity=false) override
Set up the linear system (DOF indices and sparsity patterns)
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 applyBoundaryConditions(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
apply boundary condition to system
virtual void resetStateToBeginningOfStep(DomainPartition &domain) override
reset state of physics back to the beginning of the step.
virtual void updateState(DomainPartition &domain) override
Recompute all dependent quantities from primary variables (including constitutive models)
virtual real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
entry function to perform a solver step
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
array1d< globalIndex > globalIndex_array
A 1-dimensional array of geos::globalIndex types.
Definition: DataTypes.hpp:370
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
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
LvArray::SortedArray< T, localIndex, LvArray::ChaiBuffer > SortedArray
A sorted array of local indices.
Definition: DataTypes.hpp:266
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "richardson", "preconditioner")
Declare strings associated with enumeration values.
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.
Structure to hold scoped key names.