GEOS
SolidMechanicsALMKernels.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_CONTACT_KERNELS_SOLIDMECHANICSALMKERNELS_HPP_
21 #define GEOS_PHYSICSSOLVERS_CONTACT_KERNELS_SOLIDMECHANICSALMKERNELS_HPP_
22 
23 #include "SolidMechanicsConformingContactKernelsBase.hpp"
24 #include "mesh/MeshFields.hpp"
25 
26 
27 namespace geos
28 {
29 
30 namespace solidMechanicsALMKernels
31 {
32 
36 template< typename CONSTITUTIVE_TYPE,
37  typename FE_TYPE >
38 class ALM :
40  FE_TYPE >
41 {
42 public:
45  FE_TYPE >;
46 
50 
52  static constexpr int numQuadraturePointsPerElem = FE_TYPE::numQuadraturePoints;
53 
55  static constexpr int numUdofs = Base::numUdofs;
56 
58  static constexpr int numBdofs = Base::numBdofs;
59 
61  static constexpr int numTdofs = Base::numTdofs;
62 
64  using Base::m_faceToNodes;
67  using Base::m_dofNumber;
68  using Base::m_bDofNumber;
70  using Base::m_X;
72  using Base::m_dispJump;
73  using Base::m_oldDispJump;
74  using Base::m_matrix;
75  using Base::m_rhs;
76 
81  ALM( NodeManager const & nodeManager,
82  EdgeManager const & edgeManager,
83  FaceManager const & faceManager,
84  localIndex const targetRegionIndex,
85  FaceElementSubRegion & elementSubRegion,
86  FE_TYPE const & finiteElementSpace,
87  CONSTITUTIVE_TYPE & inputConstitutiveType,
88  arrayView1d< globalIndex const > const uDofNumber,
89  arrayView1d< globalIndex const > const bDofNumber,
90  globalIndex const rankOffset,
92  arrayView1d< real64 > const inputRhs,
93  real64 const inputDt,
94  arrayView1d< localIndex const > const & faceElementList,
95  bool const isSymmetric ):
96  Base( nodeManager,
97  edgeManager,
98  faceManager,
99  targetRegionIndex,
100  elementSubRegion,
101  finiteElementSpace,
102  inputConstitutiveType,
103  uDofNumber,
104  bDofNumber,
105  rankOffset,
106  inputMatrix,
107  inputRhs,
108  inputDt,
109  faceElementList ),
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() )
114  {}
115 
116  //***************************************************************************
117 
122  {
123 
124 public:
125 
127  StackVariables():
130  dispColIndices{},
131  bEqnRowIndices{},
132  bColIndices{},
133  localRu{},
134  localRb{},
135  localAutAtu{ {} },
136  localAbtAtb{ {} },
137  localAbtAtu{ {} },
138  localAutAtb{ {} },
139  tLocal{}
140  {}
141 
144 
147 
150 
153 
156 
159 
162 
165 
168 
171 
174 
175  };
176 
177  //***************************************************************************
178 
179  //START_kernelLauncher
180  template< typename POLICY,
181  typename KERNEL_TYPE >
182  static
183  real64
184  kernelLaunch( localIndex const numElems,
185  KERNEL_TYPE const & kernelComponent )
186  {
187  return Base::template kernelLaunch< POLICY, KERNEL_TYPE >( numElems, kernelComponent );
188  }
189  //END_kernelLauncher
190 
196  inline
197  void setup( localIndex const k,
198  StackVariables & stack ) const
199  {
200  constexpr int shift = numNodesPerElem * 3;
201 
202  int permutation[numNodesPerElem];
203  m_finiteElementSpace.getPermutation( permutation );
204 
205  localIndex const kf0 = m_elemsToFaces[k][0];
206  localIndex const kf1 = m_elemsToFaces[k][1];
207  for( localIndex a=0; a<numNodesPerElem; ++a )
208  {
209  localIndex const kn0 = m_faceToNodes( kf0, a );
210  localIndex const kn1 = m_faceToNodes( kf1, a );
211 
212  for( int i=0; i<3; ++i )
213  {
214  stack.dispEqnRowIndices[a*3+i] = m_dofNumber[kn0]+i-m_dofRankOffset;
215  stack.dispEqnRowIndices[shift + a*3+i] = m_dofNumber[kn1]+i-m_dofRankOffset;
216  stack.dispColIndices[a*3+i] = m_dofNumber[kn0]+i;
217  stack.dispColIndices[shift + a*3+i] = m_dofNumber[kn1]+i;
218  stack.X[ a ][ i ] = m_X[ m_faceToNodes( kf0, permutation[ a ] ) ][ i ];
219  }
220  }
221 
222  for( int j=0; j<3; ++j )
223  {
224  for( int i=0; i<3; ++i )
225  {
226  stack.localRotationMatrix[ i ][ j ] = m_rotationMatrix( k, i, j );
227  }
228  }
229 
230  for( int i=0; i<numTdofs; ++i )
231  {
232  stack.tLocal[i] = m_traction( k, i );
233  stack.dispJumpLocal[i] = m_dispJump( k, i );
234  stack.oldDispJumpLocal[i] = m_oldDispJump( k, i );
235  }
236 
237  for( int i=0; i<3; ++i )
238  {
239  // need to grab the index.
240  stack.bEqnRowIndices[i] = m_bDofNumber[kf0] + i - m_dofRankOffset;
241  stack.bEqnRowIndices[3+i] = m_bDofNumber[kf1] + i - m_dofRankOffset;
242  stack.bColIndices[i] = m_bDofNumber[kf0] + i;
243  stack.bColIndices[3+i] = m_bDofNumber[kf1] + i;
244  }
245  }
246 
248  inline
249  real64 complete( localIndex const k,
250  StackVariables & stack ) const
251  {
252  GEOS_UNUSED_VAR( k );
253  constexpr real64 zero = LvArray::NumericLimits< real64 >::epsilon;
254 
255  real64 matRRtAtu[3][numUdofs], matDRtAtu[3][numUdofs];
256  real64 matRRtAtb[3][numBdofs], matDRtAtb[3][numBdofs];
257 
258  real64 tractionR[numUdofs];
259  real64 tractionRb[numBdofs];
260 
261  real64 tractionNew[3];
262 
263  integer fractureState;
264  m_constitutiveUpdate.updateTraction( m_oldDispJump[k],
265  m_dispJump[k],
266  m_penalty[k],
267  m_traction[k],
268  m_faceArea[k],
269  m_symmetric,
270  m_symmetric,
271  zero,
272  zero,
273  stack.localPenalty,
274  tractionNew,
275  fractureState );
276 
277  // transp(R) * Atu
278  LvArray::tensorOps::Rij_eq_AkiBkj< 3, numUdofs, 3 >( matRRtAtu, stack.localRotationMatrix,
279  stack.localAtu );
280  // transp(R) * Atb
281  LvArray::tensorOps::Rij_eq_AkiBkj< 3, numBdofs, 3 >( matRRtAtb, stack.localRotationMatrix,
282  stack.localAtb );
283 
284  // Compute the traction contribute of the local residuals
285  LvArray::tensorOps::Ri_eq_AjiBj< numUdofs, 3 >( tractionR, matRRtAtu, tractionNew );
286  LvArray::tensorOps::Ri_eq_AjiBj< numBdofs, 3 >( tractionRb, matRRtAtb, tractionNew );
287 
288  // D*RtAtu
289  LvArray::tensorOps::Rij_eq_AikBkj< 3, numUdofs, 3 >( matDRtAtu, stack.localPenalty,
290  matRRtAtu );
291  // D*RtAtb
292  LvArray::tensorOps::Rij_eq_AikBkj< 3, numBdofs, 3 >( matDRtAtb, stack.localPenalty,
293  matRRtAtb );
294 
295  // R*RtAtu
296  LvArray::tensorOps::Rij_eq_AikBkj< 3, numUdofs, 3 >( matRRtAtu, stack.localRotationMatrix,
297  matDRtAtu );
298  // R*RtAtb
299  LvArray::tensorOps::Rij_eq_AikBkj< 3, numBdofs, 3 >( matRRtAtb, stack.localRotationMatrix,
300  matDRtAtb );
301 
302  // transp(Atu)*RRtAtu
303  LvArray::tensorOps::Rij_eq_AkiBkj< numUdofs, numUdofs, 3 >( stack.localAutAtu, stack.localAtu,
304  matRRtAtu );
305  // transp(Atb)*RRtAtb
306  LvArray::tensorOps::Rij_eq_AkiBkj< numBdofs, numBdofs, 3 >( stack.localAbtAtb, stack.localAtb,
307  matRRtAtb );
308 
309  // transp(Atb)*RRtAtu
310  LvArray::tensorOps::Rij_eq_AkiBkj< numBdofs, numUdofs, 3 >( stack.localAbtAtu, stack.localAtb,
311  matRRtAtu );
312 
313  // transp(Atu)*RRtAtb
314  LvArray::tensorOps::Rij_eq_AkiBkj< numUdofs, numBdofs, 3 >( stack.localAutAtb, stack.localAtu,
315  matRRtAtb );
316 
317  // Compute the local residuals
318  LvArray::tensorOps::scaledAdd< numUdofs >( stack.localRu, tractionR, -1 );
319 
320  LvArray::tensorOps::scaledAdd< numBdofs >( stack.localRb, tractionRb, -1 );
321 
322  for( localIndex i=0; i < numUdofs; ++i )
323  {
324  localIndex const dof = LvArray::integerConversion< localIndex >( stack.dispEqnRowIndices[ i ] );
325 
326  if( dof < 0 || dof >= m_matrix.numRows() ) continue;
327 
328  // Is it necessary? Each row should be indepenedent
329  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[dof], stack.localRu[i] );
330 
331  // Fill in matrix
332  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
333  stack.dispColIndices,
334  stack.localAutAtu[i],
335  numUdofs );
336 
337  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
338  stack.bColIndices,
339  stack.localAutAtb[i],
340  numBdofs );
341  }
342 
343  for( localIndex i=0; i < numBdofs; ++i )
344  {
345  localIndex const dof = LvArray::integerConversion< localIndex >( stack.bEqnRowIndices[ i ] );
346 
347  if( dof < 0 || dof >= m_matrix.numRows() ) continue;
348 
349  // Is it necessary? Each row should be indepenedent
350  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[dof], stack.localRb[i] );
351 
352  // Fill in matrix
353  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
354  stack.bColIndices,
355  stack.localAbtAtb[i],
356  numBdofs );
357 
358  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
359  stack.dispColIndices,
360  stack.localAbtAtu[i],
361  numUdofs );
362  }
363 
364  return 0.0;
365  }
366 
367 protected:
368 
369  arrayView2d< real64 const > const m_traction;
370 
371  bool const m_symmetric;
372 
375 
376  arrayView1d< real64 const > const m_faceArea;
377 };
378 
383  globalIndex const,
385  arrayView1d< real64 > const,
386  real64 const,
388  bool const >;
389 
394 {
395 
407  template< typename POLICY, typename CONTACT_WRAPPER >
408  static void
409  launch( localIndex const size,
410  CONTACT_WRAPPER const & contactWrapper,
411  arrayView2d< real64 const > const & penalty,
412  arrayView2d< real64 const > const & traction,
413  arrayView2d< real64 const > const & dispJump,
414  arrayView2d< real64 const > const & deltaDispJump,
415  arrayView1d< real64 const > const & faceArea,
416  arrayView2d< real64 > const & tractionNew )
417  {
418 
419  forAll< POLICY >( size, [=] GEOS_HOST_DEVICE ( localIndex const k )
420  {
421 
422  contactWrapper.updateTractionOnly( dispJump[k], deltaDispJump[k],
423  penalty[k], traction[k], faceArea[k], tractionNew[k] );
424 
425  } );
426  }
427 };
428 
429 } // namespace SolidMechanicsALMKernels
430 
431 } // namespace geos
432 
433 
434 #endif /* GEOS_PHYSICSSOLVERS_CONTACT_KERNELS_SOLIDMECHANICSALMKERNELS_HPP_ */
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
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
Used to forward arguments to a class that implements the InterfaceKernelBase interface.
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.
arrayView2d< real64 > const m_dispJump
The array containing the displacement jump.
ArrayOfArraysView< localIndex const > const m_faceToNodes
The array of array containing the face to node map.
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.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
arrayView3d< real64 const > const m_rotationMatrix
The array containing the rotation matrix for each 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.
arrayView2d< real64 > const m_dispJump
The array containing the displacement jump.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
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
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
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.
real64 oldDispJumpLocal[numTdofs]
Stack storage for the element local old displacement jump vector.