GEOS
SinglePhasePoromechanicsDamage_impl.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_MULTIPHYSICS_SINGLEPHASEPOROMECHANICSDAMAGE_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_SINGLEPHASEPOROMECHANICSDAMAGE_IMPL_HPP_
22 
23 #include "constitutive/fluid/singlefluid/SingleFluidBase.hpp"
26 #include "physicsSolvers/multiphysics/poromechanicsKernels/SinglePhasePoromechanicsDamage.hpp"
27 
28 namespace geos
29 {
30 
31 namespace poromechanicsDamageKernels
32 {
33 
34 template< typename SUBREGION_TYPE,
35  typename CONSTITUTIVE_TYPE,
36  typename FE_TYPE >
38 SinglePhasePoromechanicsDamage( 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 inputDispDofNumber,
46  globalIndex const rankOffset,
48  arrayView1d< real64 > const inputRhs,
49  real64 const inputDt,
50  real64 const (&gravityVector)[3],
51  string const inputFlowDofKey,
52  integer const performStressInitialization,
53  string const fluidModelKey ):
54  Base( nodeManager,
55  edgeManager,
56  faceManager,
57  targetRegionIndex,
58  elementSubRegion,
59  finiteElementSpace,
60  inputConstitutiveType,
61  inputDispDofNumber,
62  rankOffset,
63  inputMatrix,
64  inputRhs,
65  inputDt,
66  gravityVector,
67  inputFlowDofKey,
68  performStressInitialization,
69  fluidModelKey ),
70  m_fluidPressureGradient( elementSubRegion.template getReference< array2d< real64 > >( "pressureGradient" ) )
71 {}
72 
73 template< typename SUBREGION_TYPE,
74  typename CONSTITUTIVE_TYPE,
75  typename FE_TYPE >
80  localIndex const q,
81  StackVariables & stack ) const
82 {
83  real64 porosity = 0.0;
84  real64 porosity_n = 0.0;
85  real64 dPorosity_dVolStrain = 0.0;
86  real64 dPorosity_dPressure = 0.0;
87  real64 dPorosity_dTemperature = 0.0;
88  real64 dSolidDensity_dPressure = 0.0;
89 
90  real64 fluidPressureGradient[3]{};
91 
92  for( integer i=0; i<3; ++i )
93  {
94  fluidPressureGradient[i] = m_fluidPressureGradient( k, i );
95  }
96 
97  // Step 1: call the constitutive model to evaluate the total stress and compute porosity
98  m_constitutiveUpdate.smallStrainUpdatePoromechanics( k, q,
99  m_dt,
100  m_pressure_n[k],
101  m_pressure[k],
102  stack.temperature,
104  stack.strainIncrement,
105  stack.totalStress,
108  stack.stiffness,
109  m_performStressInitialization,
110  porosity,
111  porosity_n,
112  dPorosity_dVolStrain,
113  dPorosity_dPressure,
114  dPorosity_dTemperature,
115  dSolidDensity_dPressure );
116 
117  // Step 2: compute fracture flow term and its derivative w.r.t pressure
118  m_constitutiveUpdate.computeFractureFlowTerm( k, q,
119  fluidPressureGradient,
120  stack.fractureFlowTerm,
122 
123  // Step 3: compute the body force
124  Base::computeBodyForce( k, q,
125  porosity,
126  dPorosity_dVolStrain,
127  dPorosity_dPressure,
128  dPorosity_dTemperature,
129  dSolidDensity_dPressure,
130  stack );
131 
132  // Step 3: compute fluid mass increment
133  Base::computeFluidIncrement( k, q,
134  porosity,
135  porosity_n,
136  dPorosity_dVolStrain,
137  dPorosity_dPressure,
138  dPorosity_dTemperature,
139  stack );
140 }
141 
142 template< typename SUBREGION_TYPE,
143  typename CONSTITUTIVE_TYPE,
144  typename FE_TYPE >
148 assembleMomentumBalanceTerms( real64 const ( &N )[numNodesPerElem],
149  real64 const ( &dNdX )[numNodesPerElem][3],
150  real64 const & detJxW,
151  StackVariables & stack ) const
152 {
153  using namespace PDEUtilities;
154 
155  constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
156  constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
157  constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
158 
159  // Step 1: compute local linear momentum balance residual
160  LinearFormUtilities::compute< displacementTestSpace,
161  DifferentialOperator::SymmetricGradient >
162  (
163  stack.localResidualMomentum,
164  dNdX,
165  stack.totalStress,
166  -detJxW );
167 
168  if( m_gravityAcceleration > 0.0 )
169  {
170  LinearFormUtilities::compute< displacementTestSpace,
171  DifferentialOperator::Identity >
172  (
173  stack.localResidualMomentum,
174  N,
175  stack.bodyForce,
176  detJxW );
177  }
178 
179  LinearFormUtilities::compute< displacementTestSpace,
180  DifferentialOperator::Identity >
181  (
182  stack.localResidualMomentum,
183  N,
184  stack.fractureFlowTerm,
185  detJxW );
186 
187  // Step 2: compute local linear momentum balance residual derivatives with respect to displacement
188  BilinearFormUtilities::compute< displacementTestSpace,
189  displacementTrialSpace,
190  DifferentialOperator::SymmetricGradient,
191  DifferentialOperator::SymmetricGradient >
192  (
194  dNdX,
195  stack.stiffness, // fourth-order tensor handled via DiscretizationOps
196  dNdX,
197  -detJxW );
198 
199  if( m_gravityAcceleration > 0.0 )
200  {
201  BilinearFormUtilities::compute< displacementTestSpace,
202  displacementTrialSpace,
203  DifferentialOperator::Identity,
204  DifferentialOperator::Divergence >
205  (
207  N,
209  dNdX,
210  detJxW );
211  }
212 
213  // Step 3: compute local linear momentum balance residual derivatives with respect to pressure
214  BilinearFormUtilities::compute< displacementTestSpace,
215  pressureTrialSpace,
216  DifferentialOperator::SymmetricGradient,
217  DifferentialOperator::Identity >
218  (
220  dNdX,
222  1.0,
223  -detJxW );
224 
225  if( m_gravityAcceleration > 0.0 )
226  {
227  BilinearFormUtilities::compute< displacementTestSpace,
228  pressureTrialSpace,
229  DifferentialOperator::Identity,
230  DifferentialOperator::Identity >
231  (
233  N,
234  stack.dBodyForce_dPressure,
235  1.0,
236  detJxW );
237  }
238 
239  BilinearFormUtilities::compute< displacementTestSpace,
240  pressureTrialSpace,
241  DifferentialOperator::Identity,
242  DifferentialOperator::Identity >
243  (
245  N,
247  1.0,
248  detJxW );
249 }
250 
251 template< typename SUBREGION_TYPE,
252  typename CONSTITUTIVE_TYPE,
253  typename FE_TYPE >
258  localIndex const q,
259  StackVariables & stack ) const
260 {
261  // Step 1: compute displacement finite element basis functions (N), basis function derivatives (dNdX), and
262  // determinant of the Jacobian transformation matrix times the quadrature weight (detJxW)
263  real64 N[numNodesPerElem]{};
264  real64 dNdX[numNodesPerElem][3]{};
265  FE_TYPE::calcN( q, stack.feStack, N );
266  real64 const detJxW = m_finiteElementSpace.template getGradN< FE_TYPE >( k, q, stack.xLocal,
267  stack.feStack, dNdX );
268 
269  // Step 2: compute strain increment
270  LvArray::tensorOps::fill< 6 >( stack.strainIncrement, 0.0 );
271  FE_TYPE::symmetricGradient( dNdX, stack.uhat_local, stack.strainIncrement );
272 
273  // Step 3: compute 1) the total stress, 2) the body force terms, and 3) the fluidMassIncrement
274  // using quantities returned by the PorousSolid constitutive model.
275  // This function also computes the derivatives of these three quantities wrt primary variables
276  smallStrainUpdate( k, q, stack );
277 
278  // Step 4: use the total stress and the body force to increment the local momentum balance residual
279  // This function also fills the local Jacobian rows corresponding to the momentum balance.
280  assembleMomentumBalanceTerms( N, dNdX, detJxW, stack );
281 
282  // Step 5: use the fluid mass increment to increment the local mass balance residual
283  // This function also fills the local Jacobian rows corresponding to the mass balance.
284  Base::assembleElementBasedFlowTerms( dNdX, detJxW, stack );
285 }
286 
287 template< typename SUBREGION_TYPE,
288  typename CONSTITUTIVE_TYPE,
289  typename FE_TYPE >
293 complete( localIndex const k,
294  StackVariables & stack ) const
295 {
296  real64 const maxForce = Base::complete( k, stack );
297 
298  return maxForce;
299 }
300 
301 template< typename SUBREGION_TYPE,
302  typename CONSTITUTIVE_TYPE,
303  typename FE_TYPE >
304 template< typename POLICY,
305  typename KERNEL_TYPE >
307 kernelLaunch( localIndex const numElems,
308  KERNEL_TYPE const & kernelComponent )
309 {
311 
312  // Define a RAJA reduction variable to get the maximum residual contribution.
313  RAJA::ReduceMax< ReducePolicy< POLICY >, real64 > maxResidual( 0 );
314 
315  forAll< POLICY >( numElems,
316  [=] GEOS_HOST_DEVICE ( localIndex const k )
317  {
318  typename KERNEL_TYPE::StackVariables stack;
319 
320  kernelComponent.setup( k, stack );
321  for( integer q=0; q<KERNEL_TYPE::numQuadraturePointsPerElem; ++q )
322  {
323  kernelComponent.quadraturePointKernel( k, q, stack );
324  }
325  maxResidual.max( kernelComponent.complete( k, stack ) );
326  } );
327  return maxResidual.get();
328 }
329 
330 } // namespace poromechanicsDamageKernels
331 
332 } // namespace geos
333 
334 #endif // GEOS_PHYSICSSOLVERS_MULTIPHYSICS_SINGLEPHASEPOROMECHANICSDAMAGE_IMPL_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:51
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
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
Implements kernels for solving quasi-static single-phase poromechanics with phase-field damage.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
GEOS_HOST_DEVICE void smallStrainUpdate(localIndex const k, localIndex const q, StackVariables &stack) const
Helper function to compute 1) the total stress, 2) the body force term, and 3) the fluidMassIncrement...
GEOS_HOST_DEVICE void assembleMomentumBalanceTerms(real64 const (&N)[numNodesPerElem], real64 const (&dNdX)[numNodesPerElem][3], real64 const &detJxW, StackVariables &stack) const
Assemble the local linear momentum balance residual and derivatives using total stress and body force...
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
SinglePhasePoromechanicsDamage(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 inputDispDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, real64 const (&gravityVector)[3], string const inputFlowDofKey, integer const performStressInitialization, string const fluidModelKey)
Constructor.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
Definition: DataTypes.hpp:192
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
real64 dFractureFlowTerm_dPressure[3]
Derivative of the fracture flow term wrt pressure.
real64(& dLocalResidualMomentum_dDisplacement)[numDispDofPerElem][numDispDofPerElem]
Derivative of linear momentum balance residual wrt displacement.
real64 dBodyForce_dPressure[3]
Derivative of body force wrt pressure.
real64 dTotalStress_dTemperature[6]
Derivative of total stress wrt temperature.
CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps stiffness
Derivative of total stress wrt displacement.
real64 dLocalResidualMomentum_dPressure[numDispDofPerElem][1]
Derivative of linear momentum balance residual wrt pressure.
real64 dTotalStress_dPressure[6]
Derivative of total stress wrt pressure.
real64 dBodyForce_dVolStrainIncrement[3]
Derivative of body force wrt volumetric strain increment.
real64(& localResidualMomentum)[numDispDofPerElem]
Linear momentum balance residual.
real64 deltaTemperatureFromLastStep
Delta temperature since last time step.