20 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_SINGLEPHASEPOROMECHANICSDAMAGE_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_SINGLEPHASEPOROMECHANICSDAMAGE_IMPL_HPP_
23 #include "constitutive/fluid/singlefluid/SingleFluidBase.hpp"
26 #include "physicsSolvers/multiphysics/poromechanicsKernels/SinglePhasePoromechanicsDamage.hpp"
31 namespace poromechanicsDamageKernels
34 template<
typename SUBREGION_TYPE,
35 typename CONSTITUTIVE_TYPE,
42 SUBREGION_TYPE
const & elementSubRegion,
43 FE_TYPE
const & finiteElementSpace,
44 CONSTITUTIVE_TYPE & inputConstitutiveType,
50 real64 const (&gravityVector)[3],
51 string const inputFlowDofKey,
52 integer const performStressInitialization,
53 string const fluidModelKey ):
60 inputConstitutiveType,
68 performStressInitialization,
70 m_fluidPressureGradient( elementSubRegion.template getReference<
array2d<
real64 > >(
"pressureGradient" ) )
73 template<
typename SUBREGION_TYPE,
74 typename CONSTITUTIVE_TYPE,
85 real64 dPorosity_dVolStrain = 0.0;
86 real64 dPorosity_dPressure = 0.0;
87 real64 dPorosity_dTemperature = 0.0;
88 real64 dSolidDensity_dPressure = 0.0;
90 real64 fluidPressureGradient[3]{};
94 fluidPressureGradient[i] = m_fluidPressureGradient( k, i );
98 m_constitutiveUpdate.smallStrainUpdatePoromechanics( k, q,
109 m_performStressInitialization,
112 dPorosity_dVolStrain,
114 dPorosity_dTemperature,
115 dSolidDensity_dPressure );
118 m_constitutiveUpdate.computeFractureFlowTerm( k, q,
119 fluidPressureGradient,
124 Base::computeBodyForce( k, q,
126 dPorosity_dVolStrain,
128 dPorosity_dTemperature,
129 dSolidDensity_dPressure,
133 Base::computeFluidIncrement( k, q,
136 dPorosity_dVolStrain,
138 dPorosity_dTemperature,
142 template<
typename SUBREGION_TYPE,
143 typename CONSTITUTIVE_TYPE,
149 real64 const ( &dNdX )[numNodesPerElem][3],
153 using namespace PDEUtilities;
155 constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
156 constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
157 constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
160 LinearFormUtilities::compute< displacementTestSpace,
161 DifferentialOperator::SymmetricGradient >
168 if( m_gravityAcceleration > 0.0 )
170 LinearFormUtilities::compute< displacementTestSpace,
171 DifferentialOperator::Identity >
179 LinearFormUtilities::compute< displacementTestSpace,
180 DifferentialOperator::Identity >
188 BilinearFormUtilities::compute< displacementTestSpace,
189 displacementTrialSpace,
190 DifferentialOperator::SymmetricGradient,
191 DifferentialOperator::SymmetricGradient >
199 if( m_gravityAcceleration > 0.0 )
201 BilinearFormUtilities::compute< displacementTestSpace,
202 displacementTrialSpace,
203 DifferentialOperator::Identity,
204 DifferentialOperator::Divergence >
214 BilinearFormUtilities::compute< displacementTestSpace,
216 DifferentialOperator::SymmetricGradient,
217 DifferentialOperator::Identity >
225 if( m_gravityAcceleration > 0.0 )
227 BilinearFormUtilities::compute< displacementTestSpace,
229 DifferentialOperator::Identity,
230 DifferentialOperator::Identity >
239 BilinearFormUtilities::compute< displacementTestSpace,
241 DifferentialOperator::Identity,
242 DifferentialOperator::Identity >
251 template<
typename SUBREGION_TYPE,
252 typename CONSTITUTIVE_TYPE,
259 StackVariables & stack )
const
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 );
270 LvArray::tensorOps::fill< 6 >( stack.strainIncrement, 0.0 );
271 FE_TYPE::symmetricGradient( dNdX, stack.uhat_local, stack.strainIncrement );
276 smallStrainUpdate( k, q, stack );
280 assembleMomentumBalanceTerms( N, dNdX, detJxW, stack );
284 Base::assembleElementBasedFlowTerms( dNdX, detJxW, stack );
287 template<
typename SUBREGION_TYPE,
288 typename CONSTITUTIVE_TYPE,
296 real64 const maxForce = Base::complete( k, stack );
301 template<
typename SUBREGION_TYPE,
302 typename CONSTITUTIVE_TYPE,
304 template<
typename POLICY,
305 typename KERNEL_TYPE >
308 KERNEL_TYPE
const & kernelComponent )
313 RAJA::ReduceMax< ReducePolicy< POLICY >,
real64 > maxResidual( 0 );
315 forAll< POLICY >( numElems,
318 typename KERNEL_TYPE::StackVariables stack;
320 kernelComponent.setup( k, stack );
321 for(
integer q=0; q<KERNEL_TYPE::numQuadraturePointsPerElem; ++q )
323 kernelComponent.quadraturePointKernel( k, q, stack );
325 maxResidual.max( kernelComponent.complete( k, stack ) );
327 return maxResidual.get();
#define GEOS_HOST_DEVICE
Marks a host-device function.
#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 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.
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
std::int32_t integer
Signed integer type.
real64 dFractureFlowTerm_dPressure[3]
Derivative of the fracture flow term wrt pressure.
real64 fractureFlowTerm[3]
Fracture flow term.
real64(& dLocalResidualMomentum_dDisplacement)[numDispDofPerElem][numDispDofPerElem]
Derivative of linear momentum balance residual wrt displacement.
real64 dBodyForce_dPressure[3]
Derivative of body force wrt pressure.
real64 bodyForce[3]
Body force.
real64 dTotalStress_dTemperature[6]
Derivative of total stress wrt temperature.
CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps stiffness
Derivative of total stress wrt displacement.
real64 strainIncrement[6]
Strain increment.
real64 temperature
Temperature.
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.
real64 totalStress[6]
Total stress.