GEOSX
FixedStressThermoPoromechanics_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_FIXEDSTRESSTHERMOPOROMECHANICS_IMPL_HPP_
20 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_FIXEDSTRESSTHERMOPOROMECHANICS_IMPL_HPP_
21 
23 #include "physicsSolvers/fluidFlow/FlowSolverBaseFields.hpp"
24 #include "physicsSolvers/multiphysics/PoromechanicsFields.hpp"
26 
27 namespace geos
28 {
29 
30 namespace solidMechanicsLagrangianFEMKernels
31 {
32 
33 template< typename SUBREGION_TYPE,
34  typename CONSTITUTIVE_TYPE,
35  typename FE_TYPE >
36 
38 FixedStressThermoPoromechanics( NodeManager const & nodeManager,
39  EdgeManager const & edgeManager,
40  FaceManager const & faceManager,
41  localIndex const targetRegionIndex,
42  SUBREGION_TYPE const & elementSubRegion,
43  FE_TYPE const & finiteElementSpace,
44  CONSTITUTIVE_TYPE & inputConstitutiveType,
45  arrayView1d< globalIndex const > const inputDofNumber,
46  globalIndex const rankOffset,
48  arrayView1d< real64 > const inputRhs,
49  real64 const inputDt,
50  real64 const (&inputGravityVector)[3] ):
51  Base( nodeManager,
52  edgeManager,
53  faceManager,
54  targetRegionIndex,
55  elementSubRegion,
56  finiteElementSpace,
57  inputConstitutiveType,
58  inputDofNumber,
59  rankOffset,
60  inputMatrix,
61  inputRhs,
62  inputDt ),
63  m_X( nodeManager.referencePosition()),
64  m_disp( nodeManager.getField< fields::solidMechanics::totalDisplacement >() ),
65  m_uhat( nodeManager.getField< fields::solidMechanics::incrementalDisplacement >() ),
66  m_gravityVector{ inputGravityVector[0], inputGravityVector[1], inputGravityVector[2] },
67  m_bulkDensity( elementSubRegion.template getField< fields::poromechanics::bulkDensity >() ),
68  m_pressure( elementSubRegion.template getField< fields::flow::pressure >() ),
69  m_pressure_n( elementSubRegion.template getField< fields::flow::pressure_n >() ),
70  m_initialTemperature( elementSubRegion.template getField< fields::flow::initialTemperature >() ),
71  m_temperature( elementSubRegion.template getField< fields::flow::temperature >() ),
72  m_temperature_n( elementSubRegion.template getField< fields::flow::temperature_n >() )
73 {}
74 
75 template< typename SUBREGION_TYPE,
76  typename CONSTITUTIVE_TYPE,
77  typename FE_TYPE >
81 setup( localIndex const k,
82  StackVariables & stack ) const
83 {
84  m_finiteElementSpace.template setup< FE_TYPE >( k, m_meshData, stack.feStack );
85  localIndex const numSupportPoints =
86  m_finiteElementSpace.template numSupportPoints< FE_TYPE >( stack.feStack );
87  stack.numRows = 3 * numSupportPoints;
88  stack.numCols = stack.numRows;
89 
90  for( localIndex a = 0; a < numSupportPoints; ++a )
91  {
92  localIndex const localNodeIndex = m_elemsToNodes( k, a );
93 
94  for( int i = 0; i < 3; ++i )
95  {
96 #if defined(CALC_FEM_SHAPE_IN_KERNEL)
97  stack.xLocal[ a ][ i ] = m_X[ localNodeIndex ][ i ];
98 #endif
99  stack.u_local[ a ][i] = m_disp[ localNodeIndex ][i];
100  stack.uhat_local[ a ][i] = m_uhat[ localNodeIndex ][i];
101  stack.localRowDofIndex[a*3+i] = m_dofNumber[localNodeIndex]+i;
102  stack.localColDofIndex[a*3+i] = m_dofNumber[localNodeIndex]+i;
103  }
104  }
105 
106  // Add stabilization to block diagonal parts of the local jacobian
107  // (this is a no-operation with FEM classes)
108  real64 const stabilizationScaling = computeStabilizationScaling( k );
109  m_finiteElementSpace.template addGradGradStabilizationMatrix
110  < FE_TYPE, numDofPerTrialSupportPoint, true >( stack.feStack,
111  stack.localJacobian,
112  -stabilizationScaling );
113 }
114 
115 template< typename SUBREGION_TYPE,
116  typename CONSTITUTIVE_TYPE,
117  typename FE_TYPE >
122  localIndex const q,
123  StackVariables & stack ) const
124 {
125  real64 dNdX[ numNodesPerElem ][ 3 ];
126  real64 const detJxW = m_finiteElementSpace.template getGradN< FE_TYPE >( k, q, stack.xLocal,
127  stack.feStack, dNdX );
128 
129  real64 strainInc[6] = {0};
130  real64 totalStress[6] = {0};
131 
132  typename CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps stiffness;
133 
134  FE_TYPE::symmetricGradient( dNdX, stack.uhat_local, strainInc );
135 
136  // Evaluate total stress and its derivatives
137  // TODO: allow for a customization of the kernel to pass the average pressure to the small strain update (to account for cap pressure
138  // later)
139  m_constitutiveUpdate.smallStrainUpdatePoromechanicsFixedStress( k, q,
140  m_dt,
141  m_pressure[k],
142  m_pressure_n[k],
143  m_temperature[k],
144  m_temperature_n[k],
145  strainInc,
146  totalStress,
147  stiffness );
148 
149  for( localIndex i=0; i<6; ++i )
150  {
151  totalStress[i] *= -detJxW;
152  }
153 
154  // Here we consider the bodyForce is purely from the solid
155  // Warning: here, we lag (in iteration) the displacement dependence of bulkDensity
156  real64 const gravityForce[3] = { m_gravityVector[0] * m_bulkDensity( k, q )* detJxW,
157  m_gravityVector[1] * m_bulkDensity( k, q )* detJxW,
158  m_gravityVector[2] * m_bulkDensity( k, q )* detJxW };
159 
160  real64 N[numNodesPerElem];
161  FE_TYPE::calcN( q, stack.feStack, N );
162  FE_TYPE::plusGradNajAijPlusNaFi( dNdX,
163  totalStress,
164  N,
165  gravityForce,
166  reinterpret_cast< real64 (&)[numNodesPerElem][3] >(stack.localResidual) );
167  real64 const stabilizationScaling = computeStabilizationScaling( k );
168  m_finiteElementSpace.template
169  addEvaluatedGradGradStabilizationVector< FE_TYPE,
170  numDofPerTrialSupportPoint >( stack.feStack,
171  stack.uhat_local,
172  reinterpret_cast< real64 (&)[numNodesPerElem][3] >(stack.localResidual),
173  -stabilizationScaling );
174  stiffness.template upperBTDB< numNodesPerElem >( dNdX, -detJxW, stack.localJacobian );
175 }
176 
177 template< typename SUBREGION_TYPE,
178  typename CONSTITUTIVE_TYPE,
179  typename FE_TYPE >
183 complete( localIndex const k,
184  StackVariables & stack ) const
185 {
186  GEOS_UNUSED_VAR( k );
187  real64 maxForce = 0;
188 
189  // TODO: Does this work if BTDB is non-symmetric?
190  CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps::template fillLowerBTDB< numNodesPerElem >( stack.localJacobian );
191  localIndex const numSupportPoints =
192  m_finiteElementSpace.template numSupportPoints< FE_TYPE >( stack.feStack );
193  for( int localNode = 0; localNode < numSupportPoints; ++localNode )
194  {
195  for( int dim = 0; dim < numDofPerTestSupportPoint; ++dim )
196  {
197  localIndex const dof =
198  LvArray::integerConversion< localIndex >( stack.localRowDofIndex[ numDofPerTestSupportPoint * localNode + dim ] - m_dofRankOffset );
199  if( dof < 0 || dof >= m_matrix.numRows() )
200  continue;
201  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
202  stack.localRowDofIndex,
203  stack.localJacobian[ numDofPerTestSupportPoint * localNode + dim ],
204  stack.numRows );
205 
206  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[ dof ], stack.localResidual[ numDofPerTestSupportPoint * localNode + dim ] );
207  maxForce = fmax( maxForce, fabs( stack.localResidual[ numDofPerTestSupportPoint * localNode + dim ] ) );
208  }
209  }
210  return maxForce;
211 }
212 
213 template< typename SUBREGION_TYPE,
214  typename CONSTITUTIVE_TYPE,
215  typename FE_TYPE >
216 template< typename POLICY,
217  typename KERNEL_TYPE >
218 real64
220  KERNEL_TYPE const & kernelComponent )
221 {
222  return Base::template kernelLaunch< POLICY, KERNEL_TYPE >( numElems, kernelComponent );
223 }
224 
225 } // namespace solidMechanicsLagrangianFEMKernels
226 
227 } // namespace geos
228 
229 #endif // GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_FIXEDSTRESSTHERMOPOROMECHANICS_IMPL_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:48
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:83
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:50
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
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
GEOS_HOST_DEVICE GEOS_FORCE_INLINE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
FixedStressThermoPoromechanics(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 void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack) const
Performs a state update at a quadrature point.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:350
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
GEOSX_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
GEOSX_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
real64 u_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal displacement.
real64 uhat_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal incremental displacement.