20 #ifndef GEOS_PHYSICSSOLVERS_CONTACT_KERNELS_SOLIDMECHANICSALMKERNELS_HPP_
21 #define GEOS_PHYSICSSOLVERS_CONTACT_KERNELS_SOLIDMECHANICSALMKERNELS_HPP_
23 #include "SolidMechanicsConformingContactKernelsBase.hpp"
30 namespace solidMechanicsALMKernels
36 template<
typename CONSTITUTIVE_TYPE,
86 FE_TYPE
const & finiteElementSpace,
87 CONSTITUTIVE_TYPE & inputConstitutiveType,
95 bool const isSymmetric ):
102 inputConstitutiveType,
110 m_traction( elementSubRegion.getField< fields::contact::traction >().toViewConst()),
111 m_symmetric( isSymmetric ),
112 m_penalty( elementSubRegion.getField< fields::contact::iterativePenalty >().toViewConst() ),
113 m_faceArea( elementSubRegion.getField< fields::elementArea >().toViewConst() )
180 template<
typename POLICY,
181 typename KERNEL_TYPE >
185 KERNEL_TYPE
const & kernelComponent )
187 return Base::template kernelLaunch< POLICY, KERNEL_TYPE >( numElems, kernelComponent );
212 for(
int i=0; i<3; ++i )
222 for(
int j=0; j<3; ++j )
224 for(
int i=0; i<3; ++i )
232 stack.
tLocal[i] = m_traction( k, i );
237 for(
int i=0; i<3; ++i )
250 StackVariables & stack )
const
253 constexpr
real64 zero = LvArray::NumericLimits< real64 >::epsilon;
278 LvArray::tensorOps::Rij_eq_AkiBkj< 3, numUdofs, 3 >( matRRtAtu, stack.localRotationMatrix,
281 LvArray::tensorOps::Rij_eq_AkiBkj< 3, numBdofs, 3 >( matRRtAtb, stack.localRotationMatrix,
285 LvArray::tensorOps::Ri_eq_AjiBj< numUdofs, 3 >( tractionR, matRRtAtu, tractionNew );
286 LvArray::tensorOps::Ri_eq_AjiBj< numBdofs, 3 >( tractionRb, matRRtAtb, tractionNew );
289 LvArray::tensorOps::Rij_eq_AikBkj< 3, numUdofs, 3 >( matDRtAtu, stack.localPenalty,
292 LvArray::tensorOps::Rij_eq_AikBkj< 3, numBdofs, 3 >( matDRtAtb, stack.localPenalty,
296 LvArray::tensorOps::Rij_eq_AikBkj< 3, numUdofs, 3 >( matRRtAtu, stack.localRotationMatrix,
299 LvArray::tensorOps::Rij_eq_AikBkj< 3, numBdofs, 3 >( matRRtAtb, stack.localRotationMatrix,
303 LvArray::tensorOps::Rij_eq_AkiBkj< numUdofs, numUdofs, 3 >( stack.localAutAtu, stack.localAtu,
306 LvArray::tensorOps::Rij_eq_AkiBkj< numBdofs, numBdofs, 3 >( stack.localAbtAtb, stack.localAtb,
310 LvArray::tensorOps::Rij_eq_AkiBkj< numBdofs, numUdofs, 3 >( stack.localAbtAtu, stack.localAtb,
314 LvArray::tensorOps::Rij_eq_AkiBkj< numUdofs, numBdofs, 3 >( stack.localAutAtb, stack.localAtu,
318 LvArray::tensorOps::scaledAdd< numUdofs >( stack.localRu, tractionR, -1 );
320 LvArray::tensorOps::scaledAdd< numBdofs >( stack.localRb, tractionRb, -1 );
324 localIndex const dof = LvArray::integerConversion< localIndex >( stack.dispEqnRowIndices[ i ] );
326 if( dof < 0 || dof >=
m_matrix.numRows() )
continue;
329 RAJA::atomicAdd< parallelDeviceAtomic >( &
m_rhs[dof], stack.localRu[i] );
332 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
333 stack.dispColIndices,
334 stack.localAutAtu[i],
337 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
339 stack.localAutAtb[i],
345 localIndex const dof = LvArray::integerConversion< localIndex >( stack.bEqnRowIndices[ i ] );
347 if( dof < 0 || dof >=
m_matrix.numRows() )
continue;
350 RAJA::atomicAdd< parallelDeviceAtomic >( &
m_rhs[dof], stack.localRb[i] );
353 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
355 stack.localAbtAtb[i],
358 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
359 stack.dispColIndices,
360 stack.localAbtAtu[i],
369 arrayView2d< real64 const >
const m_traction;
371 bool const m_symmetric;
407 template<
typename POLICY,
typename CONTACT_WRAPPER >
410 CONTACT_WRAPPER
const & contactWrapper,
422 contactWrapper.updateTractionOnly( dispJump[k], deltaDispJump[k],
423 penalty[k], traction[k], faceArea[k], tractionNew[k] );
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
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.
globalIndex const m_dofRankOffset
The global rank offset.
FE_TYPE const m_finiteElementSpace
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
static constexpr int maxNumTestSupportPointsPerElem
arrayView1d< real64 > const m_rhs
The global residaul vector.
Used to forward arguments to a class that implements the InterfaceKernelBase interface.
CONSTITUTIVE_TYPE::KernelWrapper const m_constitutiveUpdate
Define the base interface for implicit finite element kernels.
static constexpr int numQuadraturePointsPerElem
Compile time value for the number of quadrature points per element.
ArrayOfArraysView< localIndex const > const m_faceToNodes
The array of array containing the face to node map.
arrayView2d< real64 const > const m_penalty
The array containing the penalty coefficients for each element.
static constexpr int numTdofs
The number of lagrange multiplier dofs per element.
arrayView2d< real64 const > const m_oldDispJump
The array containing the displacement jump of previus time step.
ALM(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, FaceElementSubRegion &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const uDofNumber, arrayView1d< globalIndex const > const bDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, arrayView1d< localIndex const > const &faceElementList, bool const isSymmetric)
Constructor.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
arrayView3d< real64 const > const m_rotationMatrix
The array containing the rotation matrix for each element.
static constexpr int numUdofs
The number of displacement dofs per element.
arrayView1d< globalIndex const > const m_bDofNumber
The global degree of freedom number of bubble.
arrayView2d< localIndex const > const m_elemsToFaces
The array of array containing the element to face map.
static constexpr int numBdofs
The number of bubble dofs per element.
static constexpr int numNodesPerElem
arrayView2d< real64 > const m_dispJump
The array containing the displacement jump.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix 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).
std::int32_t integer
Signed integer type.
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Kernel variables allocated on the stack.
globalIndex dispColIndices[numUdofs]
C-array storage for the element local column degrees of freedom.
globalIndex bColIndices[numBdofs]
C-array storage for the element local column degrees of freedom.
real64 tLocal[numTdofs]
Stack storage for the element local lagrange multiplier vector.
real64 localAbtAtb[numBdofs][numBdofs]
C-array storage for the element local AbtAtb matrix.
globalIndex dispEqnRowIndices[numUdofs]
C-array storage for the element local row degrees of freedom.
real64 localAbtAtu[numBdofs][numUdofs]
C-array storage for the element local AbtAtu matrix.
real64 localAutAtu[numUdofs][numUdofs]
C-array storage for the element local AutAtu matrix.
real64 localAutAtb[numUdofs][numBdofs]
C-array storage for the element local AbtAtu matrix.
real64 localRb[numBdofs]
C-array storage for the element local Rb residual vector.
real64 localRu[numUdofs]
C-array storage for the element local Ru residual vector.
globalIndex bEqnRowIndices[numBdofs]
C-array storage for the element local row degrees of freedom.
A struct to compute the traction after nonlinear solve.
static void launch(localIndex const size, CONTACT_WRAPPER const &contactWrapper, arrayView2d< real64 const > const &penalty, arrayView2d< real64 const > const &traction, arrayView2d< real64 const > const &dispJump, arrayView2d< real64 const > const &deltaDispJump, arrayView1d< real64 const > const &faceArea, arrayView2d< real64 > const &tractionNew)
Launch the kernel function to compute the traction.