GEOS
ImplicitRateAndStateKernels.hpp
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) 2018-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 
16 #ifndef GEOS_PHYSICSSOLVERS_INDUCEDSEISMICITY_KERNELS_IMPLICITRATEANDSTATEKERNELS_HPP_
17 #define GEOS_PHYSICSSOLVERS_INDUCEDSEISMICITY_KERNELS_IMPLICITRATEANDSTATEKERNELS_HPP_
18 
19 #include "RateAndStateKernelsBase.hpp"
20 #include "denseLinearAlgebra/denseLASolvers.hpp"
21 
22 namespace geos
23 {
24 
25 namespace rateAndStateKernels
26 {
27 
35 template< typename FRICTION_LAW_TYPE >
37 {
38 public:
39 
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() )
52  {}
53 
59  {
60 public:
61 
62  StackVariables( ) = default;
63 
64  real64 jacobian[2][2]{};
65 
66  real64 rhs[2]{};
67 
68  };
69 
71  void setup( localIndex const k,
72  real64 const dt,
73  StackVariables & stack ) const
74  {
75  real64 const normalTraction = m_normalTraction[k];
76  real64 const shearTractionMagnitude = LvArray::tensorOps::l2Norm< 2 >( m_shearTraction[k] );
77 
78  // Eq 1: Scalar force balance for slipRate and shear traction magnitude
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] ) };
83 
84  // Eq 2: slip law
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] ) };
88 
89  // Assemble Jacobian matrix
90  stack.jacobian[0][0] = dFriction[0]; // derivative of Eq 1 w.r.t. stateVariable
91  stack.jacobian[0][1] = dFriction[1]; // derivative of Eq 1 w.r.t. slipRate
92  stack.jacobian[1][0] = dStateEvolutionLaw[0]; // derivative of Eq 2 w.r.t. stateVariable
93  stack.jacobian[1][1] = dStateEvolutionLaw[1]; // derivative of Eq 2 w.r.t. slipRate
94  }
95 
97  void solve( localIndex const k,
98  StackVariables & stack ) const
99  {
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 );
104 
105  // Slip rate is bracketed between [0, shear traction magnitude / shear impedance]
106  // Check that the update did not end outside of the bracket.
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];
110 
111  real64 scalingFactor = 1.0;
112  if( solution[1] > upperBound )
113  {
114  scalingFactor = 0.5 * upperBound / solution[1];
115  }
116  else if( solution[1] < lowerBound )
117  {
118  scalingFactor = 0.5 * lowerBound / solution[1];
119  }
120 
121  LvArray::tensorOps::scale< 2 >( solution, scalingFactor );
122 
123  // Update variables
124  m_stateVariable[k] += solution[0];
125  m_slipRate[k] += solution[1];
126  }
127 
129  void projectSlipRate( localIndex const k ) const
130  {
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 );
133  }
134 
136  void udpateVariables( localIndex const k ) const
137  {
138  projectSlipRate( k );
139  m_stateVariable_n[k] = m_stateVariable[k];
140  m_slipRate_n[k] = m_slipRate[k];
141  }
142 
144  camp::tuple< int, real64 > checkConvergence( StackVariables const & stack,
145  real64 const tol ) const
146  {
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 };
150  return result;
151  }
152 
154  void resetState( localIndex const k ) const
155  {
156  m_stateVariable[k] = m_stateVariable_n[k];
157  m_slipRate[k] = m_slipRate_n[k];
158  }
159 
160  template< typename POLICY >
161  static real64 solveRateAndStateEquation( SurfaceElementSubRegion & subRegion,
162  ImplicitFixedStressRateAndStateKernel & kernel,
163  real64 dt,
164  integer const maxNewtonIter,
165  real64 const newtonTol )
166  {
167  bool converged = false;
168  for( integer attempt = 0; attempt < 5; attempt++ )
169  {
170  if( attempt > 0 )
171  {
172  forAll< POLICY >( subRegion.size(), [=] GEOS_HOST_DEVICE ( localIndex const k )
173  {
174  kernel.resetState( k );
175  } );
176  }
177  GEOS_LOG_RANK_0( GEOS_FMT( " Attempt {} ", attempt ) );
178  converged = newtonSolve< POLICY >( subRegion, kernel, dt, maxNewtonIter, newtonTol );
179  if( converged )
180  {
181  forAll< POLICY >( subRegion.size(), [=] GEOS_HOST_DEVICE ( localIndex const k )
182  {
183  kernel.udpateVariables( k );
184  } );
185  return dt;
186  }
187  else
188  {
189  GEOS_LOG_RANK_0( GEOS_FMT( " Attempt {} failed. Halving dt and retrying.", attempt ) );
190  dt *= 0.5;
191  }
192  }
193  if( !converged )
194  {
195  GEOS_ERROR( "Maximum number of attempts reached without convergence." );
196  }
197  return dt;
198  }
199 
200 private:
201 
202  arrayView1d< real64 > const m_slipRate;
203 
204  arrayView1d< real64 > const m_stateVariable;
205 
206  arrayView1d< real64 > const m_stateVariable_n;
207 
208  arrayView1d< real64 > const m_slipRate_n;
209 
210  arrayView1d< real64 > const m_normalTraction;
211 
212  arrayView2d< real64 > const m_shearTraction;
213 
214  arrayView2d< real64 > const m_slipVelocity;
215 
216  real64 const m_shearImpedance;
217 
218  typename FRICTION_LAW_TYPE::KernelWrapper m_frictionLaw;
219 
220 };
221 
222 } /* namespace rateAndStateKernels */
223 
224 } /* namespace geos */
225 
226 #endif /* GEOS_PHYSICSSOLVERS_INDUCEDSEISMICITY_KERNELS_IMPLICITRATEANDSTATEKERNELS_HPP_ */
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#define GEOS_LOG_RANK_0(msg)
Log a message on screen on rank 0.
Definition: Logger.hpp:101
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.
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