GEOS
SolidMechanicsMPM.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 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 
20 #ifndef GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_MPM_HPP_
21 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_MPM_HPP_
22 
24 #include "common/TimingMacros.hpp"
26 #include "kernels/ExplicitMPM.hpp"
28 #include "mesh/mpiCommunications/CommunicationTools.hpp"
29 #include "mesh/mpiCommunications/MPI_iCommData.hpp"
32 #include "MPMSolverFields.hpp"
33 
34 namespace geos
35 {
36 
37 class SpatialPartition;
38 
39 
46 {
47 public:
48 
55  {
56  QuasiStatic,
59  };
60 
67  {
68  OUTFLOW,
69  SYMMETRY
70  };
71 
77  SolidMechanicsMPM( const string & name,
78  Group * const parent );
79 
80 
81  SolidMechanicsMPM( SolidMechanicsMPM const & ) = delete;
82  SolidMechanicsMPM( SolidMechanicsMPM && ) = default;
83 
84  SolidMechanicsMPM & operator=( SolidMechanicsMPM const & ) = delete;
85  SolidMechanicsMPM & operator=( SolidMechanicsMPM && ) = delete;
86 
90  virtual ~SolidMechanicsMPM() override;
91 
95  static string catalogName() { return "SolidMechanics_MPM"; }
99  string getCatalogName() const override { return catalogName(); }
100 
101  virtual void initializePreSubGroups() override;
102 
103  virtual void registerDataOnMesh( Group & meshBodies ) override final;
104 
111  virtual
112  real64 solverStep( real64 const & time_n,
113  real64 const & dt,
114  integer const cycleNumber,
115  DomainPartition & domain ) override;
116 
117  virtual
118  real64 explicitStep( real64 const & time_n,
119  real64 const & dt,
120  integer const cycleNumber,
121  DomainPartition & domain ) override;
122 
123  virtual void updateState( DomainPartition & domain ) override final
124  {
125  // There should be nothing to update
126  GEOS_UNUSED_VAR( domain );
127  };
128 
134  virtual bool execute( real64 const time_n,
135  real64 const dt,
136  integer const cycleNumber,
137  integer const eventCounter,
138  real64 const eventProgress,
139  DomainPartition & domain ) override;
140 
141 
142  template< typename CONSTITUTIVE_BASE,
143  typename KERNEL_WRAPPER,
144  typename ... PARAMS >
145  void assemblyLaunch( DomainPartition & domain,
146  DofManager const & dofManager,
147  CRSMatrixView< real64, globalIndex const > const & localMatrix,
148  arrayView1d< real64 > const & localRhs,
149  PARAMS && ... params );
150 
151 
152  template< typename ... PARAMS >
153  real64 explicitKernelDispatch( MeshLevel & mesh,
154  arrayView1d< string const > const & targetRegions,
155  string const & finiteElementName,
156  real64 const dt,
157  std::string const & elementListName );
158 
170  {
171  static constexpr char const * cflFactorString() { return "cflFactor"; }
172  static constexpr char const * timeIntegrationOptionString() { return "timeIntegrationOption"; }
173  static constexpr char const * solidMaterialNamesString() { return "solidMaterialNames"; }
174  static constexpr char const * forceExternalString() { return "externalForce"; }
175  static constexpr char const * forceInternalString() { return "internalForce"; }
176  static constexpr char const * massString() { return "mass"; }
177  static constexpr char const * velocityString() { return "velocity"; }
178  static constexpr char const * momentumString() { return "momentum"; }
179  static constexpr char const * accelerationString() { return "acceleration"; }
180  static constexpr char const * forceContactString() { return "contactForce"; }
181  static constexpr char const * damageString() { return "damage"; }
182  static constexpr char const * damageGradientString() { return "damageGradient"; }
183  static constexpr char const * maxDamageString() { return "maxDamage"; }
184  static constexpr char const * surfaceNormalString() { return "surfaceNormal"; }
185  static constexpr char const * materialPositionString() { return "materialPosition"; }
186 
187  static constexpr char const * boundaryNodesString() { return "boundaryNodes"; }
188  static constexpr char const * bufferNodesString() { return "bufferNodes"; }
189 
190  dataRepository::ViewKey timeIntegrationOption = { timeIntegrationOptionString() };
191  } solidMechanicsViewKeys;
192 
193  void initialize( NodeManager & nodeManager,
194  ParticleManager & particleManager,
195  SpatialPartition & partition );
196 
197  void resizeGrid( SpatialPartition & partition,
198  NodeManager & nodeManager,
199  real64 const dt );
200 
201  void syncGridFields( std::vector< std::string > const & fieldNames,
202  DomainPartition & domain,
203  NodeManager & nodeManager,
204  MeshLevel & mesh,
205  MPI_Op op );
206 
207  void singleFaceVectorFieldSymmetryBC( const int face,
208  arrayView3d< real64 > const & vectorMultiField,
210  Group & nodeSets );
211 
212  void enforceGridVectorFieldSymmetryBC( arrayView3d< real64 > const & vectorMultiField,
214  Group & nodeSets );
215 
216  void applyEssentialBCs( const real64 dt,
217  const real64 time_n,
218  NodeManager & nodeManager );
219 
220  void computeGridSurfaceNormals( ParticleManager & particleManager,
221  NodeManager & nodeManager );
222 
223  void normalizeGridSurfaceNormals( arrayView2d< real64 const > const & gridMass,
224  arrayView3d< real64 > const & gridSurfaceNormal );
225 
226  void computeContactForces( real64 const dt,
227  arrayView2d< real64 const > const & gridMass,
228  arrayView2d< real64 const > const & gridDamage,
229  arrayView2d< real64 const > const & gridMaxDamage,
230  arrayView3d< real64 const > const & gridVelocity,
231  arrayView3d< real64 const > const & gridMomentum,
232  arrayView3d< real64 const > const & gridSurfaceNormal,
233  arrayView3d< real64 const > const & gridMaterialPosition,
234  arrayView3d< real64 > const & gridContactForce );
235 
236  void computePairwiseNodalContactForce( int const & separable,
237  real64 const & dt,
238  real64 const & mA,
239  real64 const & mB,
246  arraySlice1d< real64 const > const xA, // Position of field A
247  arraySlice1d< real64 const > const xB, // Position of field B
248  arraySlice1d< real64 > const fA,
249  arraySlice1d< real64 > const fB );
250 
251  void computeOrthonormalBasis( const real64 * e1, // input "normal" unit vector.
252  real64 * e2, // output "tangential" unit vector.
253  real64 * e3 ); // output "tangential" unit vector.
254 
255  void setGridFieldLabels( NodeManager & nodeManager );
256 
257  void solverProfiling( std::string label );
258 
259  void solverProfilingIf( std::string label, bool condition );
260 
261  real64 computeNeighborList( ParticleManager & particleManager );
262 
263  void optimizeBinSort( ParticleManager & particleManager );
264 
265  real64 kernel( real64 const & r ); // distance from particle to query point
266 
267  void kernelGradient( arraySlice1d< real64 const > const x, // query point
268  std::vector< real64 > & xp, // particle location
269  real64 const & r, // distance from particle to query point
270  real64 * result );
271 
272  real64 computeKernelField( arraySlice1d< real64 const > const x, // query point
273  arrayView2d< real64 const > const xp, // List of neighbor particle locations.
274  arrayView1d< real64 const > const Vp, // List of neighbor particle volumes.
275  arrayView1d< real64 const > const fp ); // scalar field values (e.g. damage) at neighbor particles
276 
277  void computeKernelFieldGradient( arraySlice1d< real64 const > const x, // query point
278  std::vector< std::vector< real64 > > & xp, // List of neighbor particle locations.
279  std::vector< real64 > & Vp, // List of neighbor particle volumes.
280  std::vector< real64 > & fp, // scalar field values (e.g. damage) at neighbor particles
281  arraySlice1d< real64 > const result );
282 
283  void computeKernelVectorGradient( arraySlice1d< real64 const > const x, // query point
284  std::vector< std::vector< real64 > > & xp, // List of neighbor particle locations.
285  std::vector< real64 > & Vp, // List of neighbor particle volumes.
286  std::vector< std::vector< real64 > > & fp, // vector field values (e.g. velocity) at neighbor particles
287  arraySlice2d< real64 > const result );
288 
289  void computeDamageFieldGradient( ParticleManager & particleManager );
290 
291  void updateSurfaceFlagOverload( ParticleManager & particleManager );
292 
293  void projectDamageFieldGradientToGrid( ParticleManager & particleManager,
294  NodeManager & nodeManager );
295 
296  void updateDeformationGradient( real64 dt,
297  ParticleManager & particleManager );
298 
299  void updateConstitutiveModelDependencies( ParticleManager & particleManager );
300 
301  void updateStress( real64 dt,
302  ParticleManager & particleManager );
303 
304  void particleKinematicUpdate( ParticleManager & particleManager );
305 
306  void computeAndWriteBoxAverage( const real64 dt,
307  const real64 time_n,
308  ParticleManager & particleManager );
309 
310  void initializeGridFields( NodeManager & nodeManager );
311 
312  void boundaryConditionUpdate( real64 dt, real64 time_n );
313 
314  void particleToGrid( ParticleManager & particleManager,
315  NodeManager & nodeManager );
316 
317  void gridTrialUpdate( real64 dt,
318  NodeManager & nodeManager );
319 
320  void enforceContact( real64 dt,
321  DomainPartition & domain,
322  ParticleManager & particleManager,
323  NodeManager & nodeManager,
324  MeshLevel & mesh );
325 
326  void interpolateFTable( real64 dt, real64 time_n );
327 
328  void gridToParticle( real64 dt,
329  ParticleManager & particleManager,
330  NodeManager & nodeManager );
331 
332  void updateSolverDependencies( ParticleManager & particleManager );
333 
334  real64 getStableTimeStep( ParticleManager & particleManager );
335 
336  void deleteBadParticles( ParticleManager & particleManager );
337 
338  void printProfilingResults();
339 
340  void computeSurfaceFlags( ParticleManager & particleManager );
341 
342  void computeSphF( ParticleManager & particleManager );
343 
344  // void directionalOverlapCorrection( real64 dt, ParticleManager & particleManager );
345 
346  int evaluateSeparabilityCriterion( localIndex const & A,
347  localIndex const & B,
348  real64 const & damageA,
349  real64 const & damageB,
350  real64 const & maxDamageA,
351  real64 const & maxDamageB );
352 
353  void flagOutOfRangeParticles( ParticleManager & particleManager );
354 
355  void computeRVectors( ParticleManager & particleManager );
356 
357  void cpdiDomainScaling( ParticleManager & particleManager );
358 
359  void resizeMappingArrays( ParticleManager & particleManager );
360 
361  void populateMappingArrays( ParticleManager & particleManager,
362  NodeManager & nodeManager );
363 
364 protected:
365  virtual void postInputInitialization() override final;
366 
367  std::vector< array2d< localIndex > > m_mappedNodes; // mappedNodes[subregion index][particle index][node index]. dims = {# of subregions,
368  // # of particles, # of nodes a particle on the subregion maps to}
369  std::vector< array2d< real64 > > m_shapeFunctionValues; // mappedNodes[subregion][particle][nodal shape function value]. dims = {# of
370  // subregions, # of particles, # of nodes a particle on the subregion maps to}
371  std::vector< array3d< real64 > > m_shapeFunctionGradientValues; // mappedNodes[subregion][particle][nodal shape function gradient
372  // value][direction]. dims = {# of subregions, # of particles, # of nodes
373  // a particle on the subregion maps to, 3}
374 
375  int m_solverProfiling;
376  std::vector< real64 > m_profilingTimes;
377  std::vector< std::string > m_profilingLabels;
378 
379  TimeIntegrationOption m_timeIntegrationOption;
380 
381  int m_prescribedBcTable;
382  array1d< int > m_boundaryConditionTypes; // TODO: Surely there's a way to have just one variable here
383  array2d< real64 > m_bcTable;
384 
385  int m_prescribedBoundaryFTable;
386  Path m_fTablePath;
387  int m_fTableInterpType;
388  array2d< real64 > m_fTable;
389  array1d< real64 > m_domainF;
390  array1d< real64 > m_domainL;
391 
392  int m_boxAverageHistory;
393  int m_reactionHistory;
394 
395  int m_needsNeighborList;
396  real64 m_neighborRadius;
397  int m_binSizeMultiplier;
398 
399  int m_useDamageAsSurfaceFlag;
400 
401  int m_cpdiDomainScaling;
402 
403  real64 m_smallMass;
404 
405  int m_numContactGroups, m_numContactFlags, m_numVelocityFields;
406  real64 m_separabilityMinDamage;
407  int m_treatFullyDamagedAsSingleField;
408  int m_surfaceDetection;
409  int m_damageFieldPartitioning;
410  int m_contactGapCorrection;
411  // int m_directionalOverlapCorrection;
412  real64 m_frictionCoefficient;
413 
414  int m_planeStrain;
415  int m_numDims;
416 
417  array1d< real64 > m_hEl; // Grid spacing in x-y-z
418  array1d< real64 > m_xLocalMin; // Minimum local grid coordinate including ghost nodes
419  array1d< real64 > m_xLocalMax; // Maximum local grid coordinate including ghost nodes
420  array1d< real64 > m_xLocalMinNoGhost; // Minimum local grid coordinate EXCLUDING ghost nodes
421  array1d< real64 > m_xLocalMaxNoGhost; // Maximum local grid coordinate EXCLUDING ghost nodes
422  array1d< real64 > m_xGlobalMin; // Minimum global grid coordinate excluding buffer nodes
423  array1d< real64 > m_xGlobalMax; // Maximum global grid coordinate excluding buffer nodes
424  array1d< real64 > m_partitionExtent; // Length of each edge of partition including buffer and ghost cells
425  array1d< real64 > m_domainExtent; // Length of each edge of global domain excluding buffer cells
426  array1d< int > m_nEl; // Number of elements in each grid direction including buffer and ghost cells
427  array3d< int > m_ijkMap; // Map from indices in each spatial dimension to local node ID
428 
429 private:
430  struct BinKey
431  {
432  localIndex regionIndex;
433  localIndex subRegionIndex;
434  localIndex binIndex;
435 
436  bool operator==( BinKey const & other ) const
437  {
438  return (regionIndex == other.regionIndex && subRegionIndex == other.subRegionIndex && binIndex == other.binIndex);
439  }
440  };
441 
442  struct BinKeyHash
443  {
444  std::size_t operator()( BinKey const & k ) const
445  {
446  using std::size_t;
447  using std::hash;
448 
449  // Compute individual hash values for first,
450  // second and third and combine them using XOR
451  // and bit shifting:
452  return ((std::hash< localIndex >()( k.regionIndex )
453  ^ (std::hash< localIndex >()( k.subRegionIndex ) << 1)) >> 1)
454  ^ (std::hash< localIndex >()( k.binIndex ) << 1);
455  }
456  };
457 
458  void setParticlesConstitutiveNames( ParticleSubRegionBase & subRegion ) const;
459 };
460 
462  "QuasiStatic",
463  "ImplicitDynamic",
464  "ExplicitDynamic" );
465 
467  "Outflow",
468  "Symmetry" );
469 
470 //**********************************************************************************************************************
471 //**********************************************************************************************************************
472 //**********************************************************************************************************************
473 
474 
475 } /* namespace geos */
476 
477 #endif /* GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_SOLIDMECHANICSLAGRANGIANFEM_HPP_ */
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_UNUSED_PARAM(X)
Mark an unused argument and silence compiler warnings.
Definition: GeosxMacros.hpp:72
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...
Class facilitating the representation of a multi-level discretization of a MeshBody.
Definition: MeshLevel.hpp:42
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
The ParticleManager class provides an interface to ObjectManagerBase in order to manage ParticleRegio...
Class describing a file Path.
Definition: Path.hpp:33
Base class for all physics solvers.
virtual void postInputInitialization() override final
SolidMechanicsMPM(const string &name, Group *const parent)
string getCatalogName() const override
virtual ~SolidMechanicsMPM() override
virtual void registerDataOnMesh(Group &meshBodies) override final
Register wrappers that contain data on the mesh objects.
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
virtual bool execute(real64 const time_n, real64 const dt, integer const cycleNumber, integer const eventCounter, real64 const eventProgress, DomainPartition &domain) override
void initialize()
Run initialization functions on this and all subgroups.
virtual void updateState(DomainPartition &domain) override final
Recompute all dependent quantities from primary variables (including constitutive models)
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 real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
entry function to perform a solver step
Group::wrapperMap::KeyIndex ViewKey
Type alias for KeyIndexT type used for wrapper lookups.
Definition: Group.hpp:1662
bool operator==(InputFlags const left, InputFlags const right)
Comparison operator for InputFlags enumeration.
Definition: InputFlags.hpp:147
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
Definition: DataTypes.hpp:192
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
Array< T, 3, PERMUTATION > array3d
Alias for 3D array.
Definition: DataTypes.hpp:208
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
std::string string
String type.
Definition: DataTypes.hpp:91
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:200
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:184
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
std::size_t size_t
Unsigned size type.
Definition: DataTypes.hpp:79
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:212
Structure to hold scoped key names.