GEOSX
ImplicitSmallStrainQuasiStatic.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 TotalEnergies
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
19 #ifndef GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_HPP_
20 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_HPP_
21 
24 
25 namespace geos
26 {
27 
28 namespace solidMechanicsLagrangianFEMKernels
29 {
30 
55 template< typename SUBREGION_TYPE,
56  typename CONSTITUTIVE_TYPE,
57  typename FE_TYPE >
59  public finiteElement::ImplicitKernelBase< SUBREGION_TYPE,
60  CONSTITUTIVE_TYPE,
61  FE_TYPE,
62  3,
63  3 >
64 {
65 public:
67  using Base = finiteElement::ImplicitKernelBase< SUBREGION_TYPE,
68  CONSTITUTIVE_TYPE,
69  FE_TYPE,
70  3,
71  3 >;
72 
79  using Base::m_dofNumber;
81  using Base::m_matrix;
82  using Base::m_rhs;
86  using Base::m_meshData;
87  using Base::m_dt;
88 
94  ImplicitSmallStrainQuasiStatic( NodeManager const & nodeManager,
95  EdgeManager const & edgeManager,
96  FaceManager const & faceManager,
97  localIndex const targetRegionIndex,
98  SUBREGION_TYPE const & elementSubRegion,
99  FE_TYPE const & finiteElementSpace,
100  CONSTITUTIVE_TYPE & inputConstitutiveType,
101  arrayView1d< globalIndex const > const inputDofNumber,
102  globalIndex const rankOffset,
104  arrayView1d< real64 > const inputRhs,
105  real64 const inputDt,
106  real64 const (&inputGravityVector)[3] );
107 
108  //*****************************************************************************
117  {
118 public:
119 
123  Base::StackVariables(),
124  xLocal(),
125  u_local(),
126  uhat_local(),
128  {}
129 
130 #if !defined(CALC_FEM_SHAPE_IN_KERNEL)
132  int xLocal;
133 #else
135  real64 xLocal[ numNodesPerElem ][ 3 ];
136 #endif
137 
140 
143 
146  };
147  //*****************************************************************************
148 
158  void setup( localIndex const k,
159  StackVariables & stack ) const;
160 
161 
168  {
175  GEOS_HOST_DEVICE inline constexpr
176  void operator() ( localIndex const a, localIndex const b )
177  {
178  GEOS_UNUSED_VAR( a );
179  GEOS_UNUSED_VAR( b );
180  }
181 
187  GEOS_HOST_DEVICE inline constexpr
188  void operator() ( real64 (& stress)[6] )
189  {
190  GEOS_UNUSED_VAR( stress );
191  }
192  };
193 
204  template< typename STRESS_MODIFIER = NoOpFunctors >
207  localIndex const q,
208  StackVariables & stack,
209  STRESS_MODIFIER && stressModifier = NoOpFunctors{} ) const;
210 
215  inline
216  real64 complete( localIndex const k,
217  StackVariables & stack ) const;
218 
222  template< typename POLICY,
223  typename KERNEL_TYPE >
224  static real64
225  kernelLaunch( localIndex const numElems,
226  KERNEL_TYPE const & kernelComponent );
227 
228 protected:
231 
234 
237 
240 
243 
251  inline
253  {
254  // TODO: generalize this to other constitutive models (currently we assume linear elasticity).
255  return 2.0 * m_constitutiveUpdate.getShearModulus( k );
256  }
257 };
258 
262  globalIndex,
264  arrayView1d< real64 > const,
265  real64 const,
266  real64 const (&)[3] >;
267 
268 } // namespace solidMechanicsLagrangianFEMKernels
269 
270 } // namespace geos
271 
273 
274 #endif // GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:48
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:83
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:42
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:43
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:45
Define the base interface for implicit finite element kernels.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
Definition: KernelBase.hpp:256
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
FE_TYPE::template MeshData< SUBREGION_TYPE > m_meshData
Data structure containing mesh data used to setup the finite element.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
CONSTITUTIVE_TYPE::KernelWrapper const m_constitutiveUpdate
Definition: KernelBase.hpp:263
Used to forward arguments to a class that implements the KernelBase interface.
Definition: KernelBase.hpp:281
arrayView2d< real64 const, nodes::INCR_DISPLACEMENT_USD > const m_uhat
The rank-global incremental displacement array.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
GEOS_HOST_DEVICE void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack, STRESS_MODIFIER &&stressModifier=NoOpFunctors{}) const
Performs a state update at a quadrature point.
arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > const m_disp
The rank-global displacement array.
ImplicitSmallStrainQuasiStatic(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, real64 const (&inputGravityVector)[3])
Constructor.
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
GEOS_HOST_DEVICE real64 computeStabilizationScaling(localIndex const k) const
Get a parameter representative of the stiffness, used as physical scaling for the stabilization matri...
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:350
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
GEOSX_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:236
GEOSX_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
Kernel variables allocated on the stack.
Definition: KernelBase.hpp:136
Internal struct to provide no-op defaults used in the inclusion of lambda functions into kernel compo...
constexpr GEOS_HOST_DEVICE void operator()(localIndex const a, localIndex const b)
operator() no-op used for adding an additional dynamics term inside the jacobian assembly loop.
real64 uhat_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal incremental displacement.
real64 u_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal displacement.
real64 constitutiveStiffness[6][6]
Stack storage for the constitutive stiffness at a quadrature point.