GEOSX
ExplicitFiniteStrain_impl.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_EXPLICITFINITESTRAIN_IMPL_HPP_
20 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_EXPLICITFINITESTRAIN_IMPL_HPP_
21 
22 #include "constitutive/solid/SolidUtilities.hpp"
23 #include "ExplicitFiniteStrain.hpp"
24 #include "ExplicitSmallStrain_impl.hpp"
25 #include "finiteElement/Kinematics.h"
26 
27 namespace geos
28 {
29 
30 namespace solidMechanicsLagrangianFEMKernels
31 {
32 
33 template< typename SUBREGION_TYPE,
34  typename CONSTITUTIVE_TYPE,
35  typename FE_TYPE >
37  EdgeManager const & edgeManager,
38  FaceManager const & faceManager,
39  localIndex const targetRegionIndex,
40  SUBREGION_TYPE const & elementSubRegion,
41  FE_TYPE const & finiteElementSpace,
42  CONSTITUTIVE_TYPE & inputConstitutiveType,
43  real64 const dt,
44  string const elementListName ):
45  Base( nodeManager,
46  edgeManager,
47  faceManager,
48  targetRegionIndex,
49  elementSubRegion,
50  finiteElementSpace,
51  inputConstitutiveType,
52  dt,
53  elementListName )
54 {}
55 
56 template< typename SUBREGION_TYPE,
57  typename CONSTITUTIVE_TYPE,
58  typename FE_TYPE >
60 inline
62  StackVariables & stack ) const
63 {
64  for( localIndex a=0; a< numNodesPerElem; ++a )
65  {
66  localIndex const nodeIndex = m_elemsToNodes( k, a );
67  for( int i=0; i<numDofPerTrialSupportPoint; ++i )
68  {
69 #if defined(CALC_FEM_SHAPE_IN_KERNEL)
70  stack.xLocal[ a ][ i ] = m_X[ nodeIndex ][ i ];
71 #endif
72  stack.uLocal[ a ][ i ] = m_u[ nodeIndex ][ i ];
73  stack.varLocal[ a ][ i ] = m_vel[ nodeIndex ][ i ];
74  }
75  }
76 }
77 
78 template< typename SUBREGION_TYPE,
79  typename CONSTITUTIVE_TYPE,
80  typename FE_TYPE >
82 inline
84  localIndex const q,
85  StackVariables & stack ) const
86 {
87  real64 dNdX[ numNodesPerElem ][ 3 ];
88  real64 const detJ = m_finiteElementSpace.template getGradN< FE_TYPE >( k, q, stack.xLocal, dNdX );
89 
90  real64 dUhatdX[3][3] = { {0} };
91  real64 dUdX[3][3] = { {0} };
92  real64 F[3][3] = { {0} };
93  real64 Ldt[3][3] = { {0} };
94  real64 fInv[3][3] = { {0} };
95 
96  FE_TYPE::gradient( dNdX, stack.varLocal, dUhatdX );
97  FE_TYPE::gradient( dNdX, stack.uLocal, dUdX );
98 
99  LvArray::tensorOps::scale< 3, 3 >( dUhatdX, m_dt );
100 
101  // calculate du/dX
102  LvArray::tensorOps::scaledCopy< 3, 3 >( F, dUhatdX, 0.5 );
103  LvArray::tensorOps::add< 3, 3 >( F, dUdX );
104  LvArray::tensorOps::addIdentity< 3 >( F, 1.0 );
105  LvArray::tensorOps::invert< 3 >( fInv, F );
106 
107  // chain rule: calculate dv/dx^(n+1/2) = dv/dX * dX/dx^(n+1/2)
108  LvArray::tensorOps::Rij_eq_AikBkj< 3, 3, 3 >( Ldt, dUhatdX, fInv );
109 
110  // calculate gradient (end of step)
111  LvArray::tensorOps::copy< 3, 3 >( F, dUhatdX );
112  LvArray::tensorOps::add< 3, 3 >( F, dUdX );
113  LvArray::tensorOps::addIdentity< 3 >( F, 1.0 );
114  real64 const detF = LvArray::tensorOps::invert< 3 >( fInv, F );
115 
116  real64 Rot[ 3 ][ 3 ]{};
117  real64 Dadt[ 6 ]{};
118  HughesWinget( Rot, Dadt, Ldt );
119 
120  real64 stress[ 6 ]{};
121  constitutive::SolidUtilities::
122  hypoUpdate_StressOnly( m_constitutiveUpdate,
123  k, q, m_dt, Dadt, Rot, stress );
124 
125  real64 P[ 3 ][ 3 ]{};
126  LvArray::tensorOps::Rij_eq_symAikBjk< 3 >( P, stress, fInv );
127  LvArray::tensorOps::scale< 3, 3 >( P, -detJ * detF );
128 
129  FE_TYPE::plusGradNajAij( dNdX, P, stack.fLocal );
130 }
131 
132 
133 template< typename SUBREGION_TYPE,
134  typename CONSTITUTIVE_TYPE,
135  typename FE_TYPE >
136 template< typename POLICY,
137  typename KERNEL_TYPE >
138 real64
140  KERNEL_TYPE const & kernelComponent )
141 {
142  return Base::template kernelLaunch< POLICY, KERNEL_TYPE >( numElems, kernelComponent );
143 }
144 
145 #undef UPDATE_STRESS
146 
147 
148 } // namespace solidMechanicsLagrangianFEMKernels
149 
150 } // namespace geos
151 
152 #endif //GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_EXPLICITFINITESTRAIN_IMPL_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:48
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
GEOS_HOST_DEVICE void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack) const
Performs a state update at a quadrature point.
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.
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Performs the setup phase for the kernel.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
Implements kernels for solving the equations of motion using the explicit Newmark method under the sm...
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
GEOSX_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
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.