20 #ifndef GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_IMPL_HPP_
28 namespace solidMechanicsLagrangianFEMKernels
32 template<
typename SUBREGION_TYPE,
33 typename CONSTITUTIVE_TYPE,
40 SUBREGION_TYPE
const & elementSubRegion,
41 FE_TYPE
const & finiteElementSpace,
42 CONSTITUTIVE_TYPE & inputConstitutiveType,
48 real64 const (&inputGravityVector)[3] ):
55 inputConstitutiveType,
61 m_X( nodeManager.referencePosition()),
62 m_disp( nodeManager.getField< fields::solidMechanics::totalDisplacement >() ),
63 m_uhat( nodeManager.getField< fields::solidMechanics::incrementalDisplacement >() ),
64 m_gravityVector{ inputGravityVector[0], inputGravityVector[1], inputGravityVector[2] },
65 m_density( inputConstitutiveType.getDensity() )
69 template<
typename SUBREGION_TYPE,
70 typename CONSTITUTIVE_TYPE,
78 m_finiteElementSpace.template setup< FE_TYPE >( k, m_meshData, stack.feStack );
80 localIndex const numSupportPoints = m_finiteElementSpace.getNumSupportPoints( stack.feStack );
82 stack.numRows = 3 * numSupportPoints;
83 stack.numCols = stack.numRows;
86 for(
localIndex a = 0; a < numSupportPoints; ++a )
88 localIndex const localNodeIndex = m_elemsToNodes( k, a );
91 for(
int i = 0; i < numDofPerTestSupportPoint; ++i )
93 stack.
xLocal[ a ][ i ] = m_X[ localNodeIndex ][ i ];
94 stack.
u_local[ a ][i] = m_disp[ localNodeIndex ][i];
95 stack.
uhat_local[ a ][i] = m_uhat[ localNodeIndex ][i];
96 stack.localRowDofIndex[a*3+i] = m_dofNumber[localNodeIndex]+i;
97 stack.localColDofIndex[a*3+i] = m_dofNumber[localNodeIndex]+i;
102 real64 const stabilizationScaling = computeStabilizationScaling( k );
103 m_finiteElementSpace.template addGradGradStabilizationMatrix
104 < FE_TYPE, numDofPerTrialSupportPoint,
true >( stack.feStack,
106 -stabilizationScaling );
109 template<
typename SUBREGION_TYPE,
110 typename CONSTITUTIVE_TYPE,
112 template<
typename STRESS_MODIFIER >
118 STRESS_MODIFIER && stressModifier )
const
120 real64 dNdX[ numNodesPerElem ][ 3 ];
121 real64 const detJxW = FE_TYPE::calcGradN( q, stack.xLocal, stack.feStack, dNdX );
123 real64 strainInc[6] = {0};
126 typename CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps stiffness;
128 finiteElement::feOps::symmetricGradient( dNdX, stack.uhat_local, strainInc );
130 m_constitutiveUpdate.smallStrainUpdate( k, q, m_dt, strainInc, stress, stiffness );
132 stressModifier( stress );
136 stress[i] *= -detJxW;
139 real64 const gravityForce[3] = { m_gravityVector[0] * m_density( k, q )* detJxW,
140 m_gravityVector[1] * m_density( k, q )* detJxW,
141 m_gravityVector[2] * m_density( k, q )* detJxW };
143 real64 N[numNodesPerElem];
144 FE_TYPE::calcN( q, stack.feStack, N );
145 finiteElement::feOps::plusGradNajAijPlusNaFi( dNdX,
149 reinterpret_cast< real64 (&)[numNodesPerElem][3]
>(stack.localResidual) );
150 real64 const stabilizationScaling = computeStabilizationScaling( k );
151 m_finiteElementSpace.template addEvaluatedGradGradStabilizationVector< FE_TYPE, numDofPerTrialSupportPoint >( stack.feStack,
153 reinterpret_cast< real64 (&)[numNodesPerElem][3]
>(stack.localResidual),
154 -stabilizationScaling );
155 #if !defined( GEOS_USE_HIP )
156 stiffness.template upperBTDB< numNodesPerElem >( dNdX, -detJxW, stack.localJacobian );
158 stiffness.template BTDB< numNodesPerElem >( dNdX, -detJxW, stack.localJacobian );
163 template<
typename SUBREGION_TYPE,
164 typename CONSTITUTIVE_TYPE,
174 #if !defined( GEOS_USE_HIP )
176 CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps::template fillLowerBTDB< numNodesPerElem >( stack.localJacobian );
178 localIndex const numSupportPoints = m_finiteElementSpace.getNumSupportPoints( stack.feStack );
181 for(
int localNode = 0; localNode < numSupportPoints; ++localNode )
184 for(
int dim = 0; dim < numDofPerTestSupportPoint; ++dim )
186 localIndex const dof = LvArray::integerConversion< localIndex >( stack.localRowDofIndex[ numDofPerTestSupportPoint * localNode + dim ] - m_dofRankOffset );
187 if( dof < 0 || dof >= m_matrix.numRows() )
189 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
190 stack.localRowDofIndex,
191 stack.localJacobian[ numDofPerTestSupportPoint * localNode + dim ],
194 RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[ dof ], stack.localResidual[ numDofPerTestSupportPoint * localNode + dim ] );
195 maxForce = fmax( maxForce, fabs( stack.localResidual[ numDofPerTestSupportPoint * localNode + dim ] ) );
203 template<
typename SUBREGION_TYPE,
204 typename CONSTITUTIVE_TYPE,
206 template<
typename POLICY,
207 typename KERNEL_TYPE >
211 KERNEL_TYPE
const & kernelComponent )
213 return Base::template kernelLaunch< POLICY, KERNEL_TYPE >( numElems, kernelComponent );
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
This class provides an interface to ObjectManagerBase in order to manage edge data.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
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.
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 void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Kernel variables (dof numbers, jacobian and residual) located on the stack.
Kernel variables allocated on the stack.
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 xLocal[numNodesPerElem][3]
C-array stack storage for element local the nodal positions.