16 #ifndef GEOS_PHYSICSSOLVERS_INDUCEDSEISMICITY_KERNELS_IMPLICITRATEANDSTATEKERNELS_HPP_
17 #define GEOS_PHYSICSSOLVERS_INDUCEDSEISMICITY_KERNELS_IMPLICITRATEANDSTATEKERNELS_HPP_
19 #include "RateAndStateKernelsBase.hpp"
20 #include "denseLinearAlgebra/denseLASolvers.hpp"
25 namespace rateAndStateKernels
35 template<
typename FRICTION_LAW_TYPE >
41 FRICTION_LAW_TYPE
const & frictionLaw,
42 real64 const shearImpedance ):
43 m_slipRate( subRegion.
getField< fields::rateAndState::slipRate >() ),
44 m_stateVariable( subRegion.
getField< fields::rateAndState::stateVariable >() ),
45 m_stateVariable_n( subRegion.
getField< fields::rateAndState::stateVariable_n >() ),
46 m_slipRate_n( subRegion.
getField< fields::rateAndState::slipRate_n >() ),
47 m_normalTraction( subRegion.
getField< fields::rateAndState::normalTraction >() ),
48 m_shearTraction( subRegion.
getField< fields::rateAndState::shearTraction >() ),
49 m_slipVelocity( subRegion.
getField< fields::rateAndState::slipVelocity >() ),
50 m_shearImpedance( shearImpedance ),
51 m_frictionLaw( frictionLaw.createKernelUpdates() )
75 real64 const normalTraction = m_normalTraction[k];
76 real64 const shearTractionMagnitude = LvArray::tensorOps::l2Norm< 2 >( m_shearTraction[k] );
79 stack.rhs[0] = shearTractionMagnitude - m_shearImpedance * m_slipRate[k]
80 - normalTraction * m_frictionLaw.frictionCoefficient( k, m_slipRate[k], m_stateVariable[k] );
81 real64 const dFriction[2] = { -normalTraction * m_frictionLaw.dFrictionCoefficient_dStateVariable( k, m_slipRate[k], m_stateVariable[k] ),
82 -m_shearImpedance - normalTraction * m_frictionLaw.dFrictionCoefficient_dSlipRate( k, m_slipRate[k], m_stateVariable[k] ) };
85 stack.rhs[1] = (m_stateVariable[k] - m_stateVariable_n[k]) / dt - m_frictionLaw.stateEvolution( k, m_slipRate[k], m_stateVariable[k] );
86 real64 const dStateEvolutionLaw[2] = { 1.0 / dt - m_frictionLaw.dStateEvolution_dStateVariable( k, m_slipRate[k], m_stateVariable[k] ),
87 -m_frictionLaw.dStateEvolution_dSlipRate( k, m_slipRate[k], m_stateVariable[k] ) };
90 stack.jacobian[0][0] = dFriction[0];
91 stack.jacobian[0][1] = dFriction[1];
92 stack.jacobian[1][0] = dStateEvolutionLaw[0];
93 stack.jacobian[1][1] = dStateEvolutionLaw[1];
101 real64 solution[2] = {0.0, 0.0};
102 LvArray::tensorOps::scale< 2 >( stack.rhs, -1.0 );
103 denseLinearAlgebra::solve< 2 >( stack.jacobian, stack.rhs, solution );
107 real64 const shearTractionMagnitude = LvArray::tensorOps::l2Norm< 2 >( m_shearTraction[k] );
108 real64 const upperBound = shearTractionMagnitude / m_shearImpedance - m_slipRate[k];
109 real64 const lowerBound = -m_slipRate[k];
111 real64 scalingFactor = 1.0;
112 if( solution[1] > upperBound )
114 scalingFactor = 0.5 * upperBound / solution[1];
116 else if( solution[1] < lowerBound )
118 scalingFactor = 0.5 * lowerBound / solution[1];
121 LvArray::tensorOps::scale< 2 >( solution, scalingFactor );
124 m_stateVariable[k] += solution[0];
125 m_slipRate[k] += solution[1];
129 void projectSlipRate(
localIndex const k )
const
131 real64 const frictionCoefficient = m_frictionLaw.frictionCoefficient( k, m_slipRate[k], m_stateVariable[k] );
132 projectSlipRateBase( k, frictionCoefficient, m_shearImpedance, m_normalTraction, m_shearTraction, m_slipRate, m_slipVelocity );
136 void udpateVariables(
localIndex const k )
const
138 projectSlipRate( k );
139 m_stateVariable_n[k] = m_stateVariable[k];
140 m_slipRate_n[k] = m_slipRate[k];
144 camp::tuple< int, real64 > checkConvergence(
StackVariables const & stack,
147 real64 const residualNorm = LvArray::tensorOps::l2Norm< 2 >( stack.rhs );
148 int const converged = residualNorm < tol ? 1 : 0;
149 camp::tuple< int, real64 > result { converged, residualNorm };
156 m_stateVariable[k] = m_stateVariable_n[k];
157 m_slipRate[k] = m_slipRate_n[k];
160 template<
typename POLICY >
161 static real64 solveRateAndStateEquation( SurfaceElementSubRegion & subRegion,
162 ImplicitFixedStressRateAndStateKernel & kernel,
167 bool converged =
false;
168 for(
integer attempt = 0; attempt < 5; attempt++ )
174 kernel.resetState( k );
178 converged = newtonSolve< POLICY >( subRegion, kernel, dt, maxNewtonIter, newtonTol );
183 kernel.udpateVariables( k );
189 GEOS_LOG_RANK_0( GEOS_FMT(
" Attempt {} failed. Halving dt and retrying.", attempt ) );
195 GEOS_ERROR(
"Maximum number of attempts reached without convergence." );
202 arrayView1d< real64 >
const m_slipRate;
204 arrayView1d< real64 >
const m_stateVariable;
206 arrayView1d< real64 >
const m_stateVariable_n;
208 arrayView1d< real64 >
const m_slipRate_n;
210 arrayView1d< real64 >
const m_normalTraction;
212 arrayView2d< real64 >
const m_shearTraction;
214 arrayView2d< real64 >
const m_slipVelocity;
216 real64 const m_shearImpedance;
218 typename FRICTION_LAW_TYPE::KernelWrapper m_frictionLaw;
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
#define GEOS_LOG_RANK_0(msg)
Log a message on screen on rank 0.
GEOS_DECLTYPE_AUTO_RETURN getField() const
Get a view to the field associated with a trait from this ObjectManagerBase.
GEOS_HOST_DEVICE void solve(localIndex const k, StackVariables &stack) const
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
int integer
Signed integer type.
Kernel variables (dof numbers, jacobian and residual) located on the stack.
Kernel variables located on the stack.