20 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSKERNELS_SINGLEPHASEPOROMECHANICS_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSKERNELS_SINGLEPHASEPOROMECHANICS_IMPL_HPP_
23 #include "constitutive/fluid/singlefluid/SingleFluidBase.hpp"
24 #include "finiteElement/BilinearFormUtilities.hpp"
25 #include "finiteElement/elementFormulations/FiniteElementOperators.hpp"
26 #include "finiteElement/LinearFormUtilities.hpp"
32 namespace poromechanicsKernels
35 template<
typename SUBREGION_TYPE,
36 typename CONSTITUTIVE_TYPE,
43 SUBREGION_TYPE
const & elementSubRegion,
44 FE_TYPE
const & finiteElementSpace,
45 CONSTITUTIVE_TYPE & inputConstitutiveType,
51 real64 const (&gravityVector)[3],
52 string const inputFlowDofKey,
53 integer const performStressInitialization,
54 string const fluidModelKey ):
61 inputConstitutiveType,
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 )
76 template<
typename SUBREGION_TYPE,
77 typename CONSTITUTIVE_TYPE,
88 real64 dPorosity_dVolStrain = 0.0;
89 real64 dPorosity_dPressure = 0.0;
90 real64 dPorosity_dTemperature = 0.0;
91 real64 dSolidDensity_dPressure = 0.0;
94 m_constitutiveUpdate.smallStrainUpdatePoromechanics( k, q,
98 stack.deltaTemperature,
99 stack.deltaTemperatureFromLastStep,
100 stack.strainIncrement,
102 stack.dTotalStress_dPressure,
103 stack.dTotalStress_dTemperature,
105 m_performStressInitialization,
108 dPorosity_dVolStrain,
110 dPorosity_dTemperature,
111 dSolidDensity_dPressure );
114 computeBodyForce( k, q,
116 dPorosity_dVolStrain,
118 dPorosity_dTemperature,
119 dSolidDensity_dPressure,
123 computeFluidIncrement( k, q,
126 dPorosity_dVolStrain,
128 dPorosity_dTemperature,
132 template<
typename SUBREGION_TYPE,
133 typename CONSTITUTIVE_TYPE,
141 real64 const & dPorosity_dVolStrain,
142 real64 const & dPorosity_dPressure,
143 real64 const & dPorosity_dTemperature,
144 real64 const & dSolidDensity_dPressure,
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 );
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 );
160 template<
typename SUBREGION_TYPE,
161 typename CONSTITUTIVE_TYPE,
169 real64 const & porosity_n,
170 real64 const & dPorosity_dVolStrain,
171 real64 const & dPorosity_dPressure,
172 real64 const & dPorosity_dTemperature,
177 stack.
fluidMassIncrement = porosity * m_fluidDensity( k, q ) - porosity_n * m_fluidDensity_n( k, q );
182 template<
typename SUBREGION_TYPE,
183 typename CONSTITUTIVE_TYPE,
189 real64 const ( &dNdX )[numNodesPerElem][3],
193 using namespace PDEUtilities;
195 constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
196 constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
197 constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
200 LinearFormUtilities::compute< displacementTestSpace,
201 DifferentialOperator::SymmetricGradient >
203 stack.localResidualMomentum,
208 LinearFormUtilities::compute< displacementTestSpace,
209 DifferentialOperator::Identity >
211 stack.localResidualMomentum,
217 BilinearFormUtilities::compute< displacementTestSpace,
218 displacementTrialSpace,
219 DifferentialOperator::SymmetricGradient,
220 DifferentialOperator::SymmetricGradient >
222 stack.dLocalResidualMomentum_dDisplacement,
228 BilinearFormUtilities::compute< displacementTestSpace,
229 displacementTrialSpace,
230 DifferentialOperator::Identity,
231 DifferentialOperator::Divergence >
233 stack.dLocalResidualMomentum_dDisplacement,
235 stack.dBodyForce_dVolStrainIncrement,
240 BilinearFormUtilities::compute< displacementTestSpace,
242 DifferentialOperator::SymmetricGradient,
243 DifferentialOperator::Identity >
245 stack.dLocalResidualMomentum_dPressure,
247 stack.dTotalStress_dPressure,
251 BilinearFormUtilities::compute< displacementTestSpace,
253 DifferentialOperator::Identity,
254 DifferentialOperator::Identity >
256 stack.dLocalResidualMomentum_dPressure,
258 stack.dBodyForce_dPressure,
263 template<
typename SUBREGION_TYPE,
264 typename CONSTITUTIVE_TYPE,
273 using namespace PDEUtilities;
275 constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
276 constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
277 constexpr FunctionSpace pressureTestSpace = pressureTrialSpace;
280 LinearFormUtilities::compute< pressureTestSpace,
281 DifferentialOperator::Identity >
289 BilinearFormUtilities::compute< pressureTestSpace,
290 displacementTrialSpace,
291 DifferentialOperator::Identity,
292 DifferentialOperator::Divergence >
301 BilinearFormUtilities::compute< pressureTestSpace,
303 DifferentialOperator::Identity,
304 DifferentialOperator::Identity >
313 template<
typename SUBREGION_TYPE,
314 typename CONSTITUTIVE_TYPE,
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 );
373 LvArray::tensorOps::fill< 6 >( stack.strainIncrement, 0.0 );
374 finiteElement::feOps::symmetricGradient( dNdX, stack.uhat_local, stack.strainIncrement );
379 smallStrainUpdate( k, q, stack );
383 assembleMomentumBalanceTerms( N, dNdX, detJxW, stack );
387 assembleElementBasedFlowTerms( dNdX, detJxW, stack );
390 template<
typename SUBREGION_TYPE,
391 typename CONSTITUTIVE_TYPE,
402 m_finiteElementSpace.getNumSupportPoints( stack.feStack );
403 integer numDisplacementDofs = numSupportPoints * numDofPerTestSupportPoint;
405 for(
int localNode = 0; localNode < numSupportPoints; ++localNode )
407 for(
int dim = 0; dim < numDofPerTestSupportPoint; ++dim )
410 localIndex const dof = LvArray::integerConversion< localIndex >( stack.localRowDofIndex[numDofPerTestSupportPoint*localNode + dim] - m_dofRankOffset );
413 if( dof < 0 || dof >= m_matrix.numRows() )
417 m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
418 stack.localRowDofIndex,
419 stack.dLocalResidualMomentum_dDisplacement[numDofPerTestSupportPoint * localNode + dim],
420 numDisplacementDofs );
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],
431 localIndex const dof = LvArray::integerConversion< localIndex >( stack.localPressureDofIndex - m_dofRankOffset );
434 if( 0 <= dof && dof < m_matrix.numRows() )
436 m_matrix.template addToRowBinarySearchUnsorted< serialAtomic >( dof,
437 stack.localRowDofIndex,
439 numDisplacementDofs );
440 m_matrix.template addToRow< serialAtomic >( dof,
441 &stack.localPressureDofIndex,
449 template<
typename SUBREGION_TYPE,
450 typename CONSTITUTIVE_TYPE,
452 template<
typename POLICY,
453 typename KERNEL_TYPE >
457 KERNEL_TYPE
const & kernelComponent )
462 RAJA::ReduceMax< ReducePolicy< POLICY >,
real64 > maxResidual( 0 );
464 forAll< POLICY >( numElems,
467 typename KERNEL_TYPE::StackVariables stack;
469 kernelComponent.setup( k, stack );
470 for(
integer q=0; q<KERNEL_TYPE::numQuadraturePointsPerElem; ++q )
472 kernelComponent.quadraturePointKernel( k, q, stack );
474 maxResidual.max( kernelComponent.complete( k, stack ) );
476 return maxResidual.get();
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
#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.
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.
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.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
std::string string
String type.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
int integer
Signed integer type.
Kernel variables (dof numbers, jacobian and residual) located on the stack.
Kernel variables allocated on the stack.
real64 dFluidMassIncrement_dVolStrainIncrement
Derivative of mass accumulation wrt volumetric strain increment.
real64 fluidMassIncrement
Mass accumulation.
real64 localResidualMass[1]
Mass balance residual.
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.