GEOS
SinglePhasePoromechanics_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_POROMECHANICSKERNELS_SINGLEPHASEPOROMECHANICS_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSKERNELS_SINGLEPHASEPOROMECHANICS_IMPL_HPP_
22 
23 #include "constitutive/fluid/singlefluid/SingleFluidBase.hpp"
24 #include "finiteElement/BilinearFormUtilities.hpp"
25 #include "finiteElement/elementFormulations/FiniteElementOperators.hpp"
26 #include "finiteElement/LinearFormUtilities.hpp"
28 
29 namespace geos
30 {
31 
32 namespace poromechanicsKernels
33 {
34 
35 template< typename SUBREGION_TYPE,
36  typename CONSTITUTIVE_TYPE,
37  typename FE_TYPE >
39 SinglePhasePoromechanics( NodeManager const & nodeManager,
40  EdgeManager const & edgeManager,
41  FaceManager const & faceManager,
42  localIndex const targetRegionIndex,
43  SUBREGION_TYPE const & elementSubRegion,
44  FE_TYPE const & finiteElementSpace,
45  CONSTITUTIVE_TYPE & inputConstitutiveType,
46  arrayView1d< globalIndex const > const inputDispDofNumber,
47  globalIndex const rankOffset,
49  arrayView1d< real64 > const inputRhs,
50  real64 const inputDt,
51  real64 const (&gravityVector)[3],
52  string const inputFlowDofKey,
53  integer const performStressInitialization,
54  string const fluidModelKey ):
55  Base( nodeManager,
56  edgeManager,
57  faceManager,
58  targetRegionIndex,
59  elementSubRegion,
60  finiteElementSpace,
61  inputConstitutiveType,
62  inputDispDofNumber,
63  rankOffset,
64  inputMatrix,
65  inputRhs,
66  inputDt,
67  gravityVector,
68  inputFlowDofKey,
69  fluidModelKey ),
70  m_fluidDensity( elementSubRegion.template getConstitutiveModel< constitutive::SingleFluidBase >( elementSubRegion.template getReference< string >( fluidModelKey ) ).density() ),
71  m_fluidDensity_n( elementSubRegion.template getConstitutiveModel< constitutive::SingleFluidBase >( elementSubRegion.template getReference< string >( fluidModelKey ) ).density_n() ),
72  m_dFluidDensity( elementSubRegion.template getConstitutiveModel< constitutive::SingleFluidBase >( elementSubRegion.template getReference< string >( fluidModelKey ) ).dDensity() ),
73  m_performStressInitialization( performStressInitialization )
74 {}
75 
76 template< typename SUBREGION_TYPE,
77  typename CONSTITUTIVE_TYPE,
78  typename FE_TYPE >
83  localIndex const q,
84  StackVariables & stack ) const
85 {
86  real64 porosity = 0.0;
87  real64 porosity_n = 0.0;
88  real64 dPorosity_dVolStrain = 0.0;
89  real64 dPorosity_dPressure = 0.0;
90  real64 dPorosity_dTemperature = 0.0;
91  real64 dSolidDensity_dPressure = 0.0;
92 
93  // Step 1: call the constitutive model to evaluate the total stress and compute porosity
94  m_constitutiveUpdate.smallStrainUpdatePoromechanics( k, q,
95  m_dt,
96  m_pressure[k],
97  m_pressure_n[k],
98  stack.deltaTemperature,
99  stack.deltaTemperatureFromLastStep,
100  stack.strainIncrement,
101  stack.totalStress,
102  stack.dTotalStress_dPressure,
103  stack.dTotalStress_dTemperature,
104  stack.stiffness,
105  m_performStressInitialization,
106  porosity,
107  porosity_n,
108  dPorosity_dVolStrain,
109  dPorosity_dPressure,
110  dPorosity_dTemperature,
111  dSolidDensity_dPressure );
112 
113  // Step 2: compute the body force
114  computeBodyForce( k, q,
115  porosity,
116  dPorosity_dVolStrain,
117  dPorosity_dPressure,
118  dPorosity_dTemperature,
119  dSolidDensity_dPressure,
120  stack );
121 
122  // Step 3: compute fluid mass increment
123  computeFluidIncrement( k, q,
124  porosity,
125  porosity_n,
126  dPorosity_dVolStrain,
127  dPorosity_dPressure,
128  dPorosity_dTemperature,
129  stack );
130 }
131 
132 template< typename SUBREGION_TYPE,
133  typename CONSTITUTIVE_TYPE,
134  typename FE_TYPE >
139  localIndex const q,
140  real64 const & porosity,
141  real64 const & dPorosity_dVolStrain,
142  real64 const & dPorosity_dPressure,
143  real64 const & dPorosity_dTemperature,
144  real64 const & dSolidDensity_dPressure,
145  StackVariables & stack ) const
146 {
147  GEOS_UNUSED_VAR( dPorosity_dTemperature );
148 
149  real64 const mixtureDensity = ( 1.0 - porosity ) * m_solidDensity( k, q ) + porosity * m_fluidDensity( k, q );
150  real64 const dMixtureDens_dVolStrainIncrement = dPorosity_dVolStrain * ( -m_solidDensity( k, q ) + m_fluidDensity( k, q ) );
151  real64 const dMixtureDens_dPressure = dPorosity_dPressure * ( -m_solidDensity( k, q ) + m_fluidDensity( k, q ) )
152  + ( 1.0 - porosity ) * dSolidDensity_dPressure
153  + porosity * m_dFluidDensity( k, q, DerivOffset::dP );
154 
155  LvArray::tensorOps::scaledCopy< 3 >( stack.bodyForce, m_gravityVector, mixtureDensity );
156  LvArray::tensorOps::scaledCopy< 3 >( stack.dBodyForce_dVolStrainIncrement, m_gravityVector, dMixtureDens_dVolStrainIncrement );
157  LvArray::tensorOps::scaledCopy< 3 >( stack.dBodyForce_dPressure, m_gravityVector, dMixtureDens_dPressure );
158 }
159 
160 template< typename SUBREGION_TYPE,
161  typename CONSTITUTIVE_TYPE,
162  typename FE_TYPE >
167  localIndex const q,
168  real64 const & porosity,
169  real64 const & porosity_n,
170  real64 const & dPorosity_dVolStrain,
171  real64 const & dPorosity_dPressure,
172  real64 const & dPorosity_dTemperature,
173  StackVariables & stack ) const
174 {
175  GEOS_UNUSED_VAR( dPorosity_dTemperature );
176 
177  stack.fluidMassIncrement = porosity * m_fluidDensity( k, q ) - porosity_n * m_fluidDensity_n( k, q );
178  stack.dFluidMassIncrement_dVolStrainIncrement = dPorosity_dVolStrain * m_fluidDensity( k, q );
179  stack.dFluidMassIncrement_dPressure = dPorosity_dPressure * m_fluidDensity( k, q ) + porosity * m_dFluidDensity( k, q, 0 );
180 }
181 
182 template< typename SUBREGION_TYPE,
183  typename CONSTITUTIVE_TYPE,
184  typename FE_TYPE >
188 assembleMomentumBalanceTerms( real64 const ( &N )[numNodesPerElem],
189  real64 const ( &dNdX )[numNodesPerElem][3],
190  real64 const & detJxW,
191  StackVariables & stack ) const
192 {
193  using namespace PDEUtilities;
194 
195  constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
196  constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
197  constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
198 
199  // Step 1: compute local linear momentum balance residual
200  LinearFormUtilities::compute< displacementTestSpace,
201  DifferentialOperator::SymmetricGradient >
202  (
203  stack.localResidualMomentum,
204  dNdX,
205  stack.totalStress,
206  -detJxW );
207 
208  LinearFormUtilities::compute< displacementTestSpace,
209  DifferentialOperator::Identity >
210  (
211  stack.localResidualMomentum,
212  N,
213  stack.bodyForce,
214  detJxW );
215 
216  // Step 2: compute local linear momentum balance residual derivatives with respect to displacement
217  BilinearFormUtilities::compute< displacementTestSpace,
218  displacementTrialSpace,
219  DifferentialOperator::SymmetricGradient,
220  DifferentialOperator::SymmetricGradient >
221  (
222  stack.dLocalResidualMomentum_dDisplacement,
223  dNdX,
224  stack.stiffness, // fourth-order tensor handled via DiscretizationOps
225  dNdX,
226  -detJxW );
227 
228  BilinearFormUtilities::compute< displacementTestSpace,
229  displacementTrialSpace,
230  DifferentialOperator::Identity,
231  DifferentialOperator::Divergence >
232  (
233  stack.dLocalResidualMomentum_dDisplacement,
234  N,
235  stack.dBodyForce_dVolStrainIncrement,
236  dNdX,
237  detJxW );
238 
239  // Step 3: compute local linear momentum balance residual derivatives with respect to pressure
240  BilinearFormUtilities::compute< displacementTestSpace,
241  pressureTrialSpace,
242  DifferentialOperator::SymmetricGradient,
243  DifferentialOperator::Identity >
244  (
245  stack.dLocalResidualMomentum_dPressure,
246  dNdX,
247  stack.dTotalStress_dPressure,
248  1.0,
249  -detJxW );
250 
251  BilinearFormUtilities::compute< displacementTestSpace,
252  pressureTrialSpace,
253  DifferentialOperator::Identity,
254  DifferentialOperator::Identity >
255  (
256  stack.dLocalResidualMomentum_dPressure,
257  N,
258  stack.dBodyForce_dPressure,
259  1.0,
260  detJxW );
261 }
262 
263 template< typename SUBREGION_TYPE,
264  typename CONSTITUTIVE_TYPE,
265  typename FE_TYPE >
269 assembleElementBasedFlowTerms( real64 const ( &dNdX )[numNodesPerElem][3],
270  real64 const & detJxW,
271  StackVariables & stack ) const
272 {
273  using namespace PDEUtilities;
274 
275  constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
276  constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
277  constexpr FunctionSpace pressureTestSpace = pressureTrialSpace;
278 
279  // Step 1: compute local mass balance residual
280  LinearFormUtilities::compute< pressureTestSpace,
281  DifferentialOperator::Identity >
282  (
283  stack.localResidualMass,
284  1.0,
285  stack.fluidMassIncrement,
286  detJxW );
287 
288  // Step 2: compute local mass balance residual derivatives with respect to displacement
289  BilinearFormUtilities::compute< pressureTestSpace,
290  displacementTrialSpace,
291  DifferentialOperator::Identity,
292  DifferentialOperator::Divergence >
293  (
295  1.0,
297  dNdX,
298  detJxW );
299 
300  // Step 3: compute local mass balance residual derivatives with respect to pressure
301  BilinearFormUtilities::compute< pressureTestSpace,
302  pressureTrialSpace,
303  DifferentialOperator::Identity,
304  DifferentialOperator::Identity >
305  (
307  1.0,
309  1.0,
310  detJxW );
311 }
312 
313 template< typename SUBREGION_TYPE,
314  typename CONSTITUTIVE_TYPE,
315  typename FE_TYPE >
320  localIndex const q,
321  StackVariables & stack ) const
322 {
323  // Governing equations (strong form)
324  // ---------------------------------
325  //
326  // divergence( totalStress ) + bodyForce = 0 (quasi-static linear momentum balance)
327  // dFluidMass_dTime + divergence( fluidMassFlux ) = source (fluid phase mass balance)
328  //
329  // with currently the following dependencies on the strainIncrement tensor and pressure
330  //
331  // totalStress = totalStress( strainIncrement, pressure)
332  // bodyForce = bodyForce( strainIncrement, pressure)
333  // fluidMass = fluidMass( strainIncrement, pressure)
334  // fluidMassFlux = fludiMassFlux( pressure)
335  //
336  // Note that the fluidMassFlux will depend on the strainIncrement if a stress-dependent constitutive
337  // model is assumed. A dependency on pressure can also occur in the source term of the mass
338  // balance equation, e.g. if a Peaceman well model is used.
339  //
340  // In this kernel cell-based contributions to Jacobian matrix and residual
341  // vector are computed. The face-based contributions, namely associated with the fluid
342  // mass flux term are computed in a different kernel. The source term in the mass balance
343  // equation is also treated elsewhere.
344  //
345  // Integration in time is performed using a backward Euler scheme (in the mass balance
346  // equation LHS and RHS are multiplied by the timestep).
347  //
348  // Using a weak formulation of the governing equation the following terms are assembled in this kernel
349  //
350  // Rmom = - \int symmetricGradient( \eta ) : totalStress + \int \eta \cdot bodyForce = 0
351  // Rmas = \int \chi ( fluidMass - fluidMass_n) = 0
352  //
353  // dRmom_dVolStrain = - \int_Omega symmetricGradient( \eta ) : dTotalStress_dVolStrain
354  // + \int \eta \cdot dBodyForce_dVolStrain
355  // dRmom_dPressure = - \int_Omega symmetricGradient( \eta ) : dTotalStress_dPressure
356  // + \int \eta \cdot dBodyForce_dPressure
357  // dRmas_dVolStrain = \int \chi dFluidMass_dVolStrain
358  // dRmas_dPressure = \int \chi dFluidMass_dPressure
359  //
360  // with \eta and \chi test basis functions for the displacement and pressure field, respectively.
361  // A continuous interpolation is used for the displacement, with \eta continuous finite element
362  // basis functions. A piecewise-constant approximation is used for the pressure.
363 
364  // Step 1: compute displacement finite element basis functions (N), basis function derivatives (dNdX), and
365  // determinant of the Jacobian transformation matrix times the quadrature weight (detJxW)
366  real64 N[numNodesPerElem]{};
367  real64 dNdX[numNodesPerElem][3]{};
368  FE_TYPE::calcN( q, stack.feStack, N );
369  real64 const detJxW = FE_TYPE::calcGradN( q, stack.xLocal,
370  stack.feStack, dNdX );
371 
372  // Step 2: compute strain increment
373  LvArray::tensorOps::fill< 6 >( stack.strainIncrement, 0.0 );
374  finiteElement::feOps::symmetricGradient( dNdX, stack.uhat_local, stack.strainIncrement );
375 
376  // Step 3: compute 1) the total stress, 2) the body force terms, and 3) the fluidMassIncrement
377  // using quantities returned by the PorousSolid constitutive model.
378  // This function also computes the derivatives of these three quantities wrt primary variables
379  smallStrainUpdate( k, q, stack );
380 
381  // Step 4: use the total stress and the body force to increment the local momentum balance residual
382  // This function also fills the local Jacobian rows corresponding to the momentum balance.
383  assembleMomentumBalanceTerms( N, dNdX, detJxW, stack );
384 
385  // Step 5: use the fluid mass increment to increment the local mass balance residual
386  // This function also fills the local Jacobian rows corresponding to the mass balance.
387  assembleElementBasedFlowTerms( dNdX, detJxW, stack );
388 }
389 
390 template< typename SUBREGION_TYPE,
391  typename CONSTITUTIVE_TYPE,
392  typename FE_TYPE >
396 complete( localIndex const k,
397  StackVariables & stack ) const
398 {
399  GEOS_UNUSED_VAR( k );
400  real64 maxForce = 0;
401  localIndex const numSupportPoints =
402  m_finiteElementSpace.getNumSupportPoints( stack.feStack );
403  integer numDisplacementDofs = numSupportPoints * numDofPerTestSupportPoint;
404 
405  for( int localNode = 0; localNode < numSupportPoints; ++localNode )
406  {
407  for( int dim = 0; dim < numDofPerTestSupportPoint; ++dim )
408  {
409 
410  localIndex const dof = LvArray::integerConversion< localIndex >( stack.localRowDofIndex[numDofPerTestSupportPoint*localNode + dim] - m_dofRankOffset );
411 
412  // we need this check to filter out ghost nodes in the assembly
413  if( dof < 0 || dof >= m_matrix.numRows() )
414  {
415  continue;
416  }
417  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
418  stack.localRowDofIndex,
419  stack.dLocalResidualMomentum_dDisplacement[numDofPerTestSupportPoint * localNode + dim],
420  numDisplacementDofs );
421 
422  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[dof], stack.localResidualMomentum[numDofPerTestSupportPoint * localNode + dim] );
423  maxForce = fmax( maxForce, fabs( stack.localResidualMomentum[numDofPerTestSupportPoint * localNode + dim] ) );
424  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
425  &stack.localPressureDofIndex,
426  stack.dLocalResidualMomentum_dPressure[numDofPerTestSupportPoint * localNode + dim],
427  1 );
428  }
429  }
430 
431  localIndex const dof = LvArray::integerConversion< localIndex >( stack.localPressureDofIndex - m_dofRankOffset );
432 
433  // we need this check to filter out ghost cells in the assembly
434  if( 0 <= dof && dof < m_matrix.numRows() )
435  {
436  m_matrix.template addToRowBinarySearchUnsorted< serialAtomic >( dof,
437  stack.localRowDofIndex,
439  numDisplacementDofs );
440  m_matrix.template addToRow< serialAtomic >( dof,
441  &stack.localPressureDofIndex,
443  1 );
444  RAJA::atomicAdd< serialAtomic >( &m_rhs[dof], stack.localResidualMass[0] );
445  }
446  return maxForce;
447 }
448 
449 template< typename SUBREGION_TYPE,
450  typename CONSTITUTIVE_TYPE,
451  typename FE_TYPE >
452 template< typename POLICY,
453  typename KERNEL_TYPE >
454 real64
456 kernelLaunch( localIndex const numElems,
457  KERNEL_TYPE const & kernelComponent )
458 {
460 
461  // Define a RAJA reduction variable to get the maximum residual contribution.
462  RAJA::ReduceMax< ReducePolicy< POLICY >, real64 > maxResidual( 0 );
463 
464  forAll< POLICY >( numElems,
465  [=] GEOS_HOST_DEVICE ( localIndex const k )
466  {
467  typename KERNEL_TYPE::StackVariables stack;
468 
469  kernelComponent.setup( k, stack );
470  for( integer q=0; q<KERNEL_TYPE::numQuadraturePointsPerElem; ++q )
471  {
472  kernelComponent.quadraturePointKernel( k, q, stack );
473  }
474  maxResidual.max( kernelComponent.complete( k, stack ) );
475  } );
476  return maxResidual.get();
477 }
478 
479 
480 } // namespace poromechanicsKernels
481 
482 } // namespace geos
483 
484 #endif // GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSKERNELS_SINGLEPHASEPOROMECHANICS_IMPL_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
#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.
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 computeFluidIncrement(localIndex const k, localIndex const q, real64 const &porosity, real64 const &porosity_n, real64 const &dPorosity_dVolStrain, real64 const &dPorosity_dPressure, real64 const &dPorosity_dTemperature, StackVariables &stack) const
Helper function to compute the fluid mass increment and its derivatives wrt primary variables.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
GEOS_HOST_DEVICE void computeBodyForce(localIndex const k, localIndex const q, real64 const &porosity, real64 const &dPorosity_dVolStrain, real64 const &dPorosity_dPressure, real64 const &dPorosity_dTemperature, real64 const &dSolidDensity_dPressure, StackVariables &stack) const
Helper function to compute the body force term (\rho g) and its derivatives wrt primary variables.
GEOS_HOST_DEVICE void assembleElementBasedFlowTerms(real64 const (&dNdX)[numNodesPerElem][3], real64 const &detJxW, StackVariables &stack) const
Assemble the local mass balance residual and derivatives using fluid mass/energy increment.
SinglePhasePoromechanics(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.
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
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...
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
std::string string
String type.
Definition: DataTypes.hpp:90
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Kernel variables (dof numbers, jacobian and residual) located on the stack.
real64 dFluidMassIncrement_dVolStrainIncrement
Derivative of mass accumulation wrt volumetric strain increment.
real64 dLocalResidualMass_dPressure[1][1]
Derivative of mass balance residual wrt pressure.
real64 dLocalResidualMass_dDisplacement[1][Base::StackVariables::numDispDofPerElem]
Derivative of mass balance residual wrt displacement.
real64 dFluidMassIncrement_dPressure
Derivative of mass accumulation wrt pressure.