GEOS
ExplicitFiniteStrain.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_KERNELS_EXPLICITFINITESTRAIN_HPP_
21 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_EXPLICITFINITESTRAIN_HPP_
22 
23 #include "ExplicitSmallStrain.hpp"
24 #include "finiteElement/Kinematics.h"
25 
26 namespace geos
27 {
28 
29 namespace solidMechanicsLagrangianFEMKernels
30 {
31 
37 #define UPDATE_STRESS 2
38 
39 
48 template< typename SUBREGION_TYPE,
49  typename CONSTITUTIVE_TYPE,
50  typename FE_TYPE >
51 class ExplicitFiniteStrain : public ExplicitSmallStrain< SUBREGION_TYPE,
52  CONSTITUTIVE_TYPE,
53  FE_TYPE >
54 {
55 public:
57  using Base = ExplicitSmallStrain< SUBREGION_TYPE,
58  CONSTITUTIVE_TYPE,
59  FE_TYPE >;
60 
68 
69  using Base::m_dt;
70  using Base::m_u;
71  using Base::m_vel;
72  using Base::m_acc;
73 #if !defined(CALCFEMSHAPE)
74  using Base::m_X;
75 #endif
76 
80  ExplicitFiniteStrain( NodeManager & nodeManager,
81  EdgeManager const & edgeManager,
82  FaceManager const & faceManager,
83  localIndex const targetRegionIndex,
84  SUBREGION_TYPE const & elementSubRegion,
85  FE_TYPE const & finiteElementSpace,
86  CONSTITUTIVE_TYPE & inputConstitutiveType,
87  real64 const dt,
88  string const elementListName );
89 
90 
91  //*****************************************************************************
96  {
100 
101 
107  Base::StackVariables(),
108  uLocal{ {0.0} }
109  {}
110 
113  };
114  //*****************************************************************************
115 
116 
121  void setup( localIndex const k,
122  StackVariables & stack ) const;
123 
128  void quadraturePointKernel( localIndex const k,
129  localIndex const q,
130  StackVariables & stack ) const;
131 
132 
136  template< typename POLICY,
137  typename KERNEL_TYPE >
138  static real64
139  kernelLaunch( localIndex const numElems,
140  KERNEL_TYPE const & kernelComponent );
141 
142 
143 
144 };
147  real64,
148  string const >;
149 
150 } // namespace solidMechanicsLagrangianFEMKernels
151 
152 } // namespace geos
153 
154 #endif //GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_EXPLICITFINITESTRAIN_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:43
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:44
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
Definition: KernelBase.hpp:257
arrayView1d< integer const > const m_elemGhostRank
The element ghost rank array.
Definition: KernelBase.hpp:260
Used to forward arguments to a class that implements the KernelBase interface.
Definition: KernelBase.hpp:282
Implements kernels for solving the equations of motion using the explicit Newmark method under the fi...
ExplicitFiniteStrain(NodeManager &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, real64 const dt, string const elementListName)
Constructor.
Implements kernels for solving the equations of motion using the explicit Newmark method under the sm...
arrayView2d< real64 const, nodes::VELOCITY_USD > const m_vel
The array containing the nodal velocity array.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
real64 const m_dt
The time increment for this time integration step.
arrayView2d< real64, nodes::ACCELERATION_USD > const m_acc
arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > const m_u
The array containing the nodal displacement array.
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
real64 uLocal[numNodesPerElem][numDofPerTrialSupportPoint]
Local stack storage for nodal displacements.
real64 varLocal[numNodesPerElem][numDofPerTestSupportPoint]
C-array stack storage for element local primary variable values.
real64 fLocal[numNodesPerElem][numDofPerTrialSupportPoint]
C-array stack storage for the element local force.