Calculation of (linear) elasticity model terms for the residual.
More...
#include <elasticitylocalresidualtpsa.hpp>
|
| template<class LhsEval > |
| static void | computeVolumeTerm (Dune::FieldVector< LhsEval, numEq > &volTerm, const MaterialState &materialState, const Problem &problem, const unsigned globalIndex) |
| | Calculate volume terms in TPSA formulation. More...
|
| |
| static void | computeFaceTerm (Dune::FieldVector< Evaluation, numEq > &faceTerm, const MaterialState &materialStateIn, const MaterialState &materialStateEx, Problem &problem, const unsigned globalIndexIn, const unsigned globalIndexEx) |
| | Calculate terms across cell faces in TPSA formulation. More...
|
| |
| template<class BoundaryConditionData > |
| static void | computeBoundaryTerm (Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex) |
| | Calculate boundary conditions in TPSA formulation given by BCCON/BCPROP. More...
|
| |
| template<class BoundaryConditionData > |
| static void | computeBoundaryTermFixed (Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex) |
| | Calculate fixed displacement boundary condition in TPSA formulation. More...
|
| |
| template<class BoundaryConditionData > |
| static void | computeBoundaryTermFree (Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex) |
| | Calculate free (or zero traction) boundary condition in TPSA formulation. More...
|
| |
| static void | computeSourceTerm (Dune::FieldVector< Evaluation, numEq > &sourceTerm, Problem &problem, unsigned globalSpaceIdex, unsigned timeIdx) |
| | Calculate source term in TPSA formulation. More...
|
| |
template<class TypeTag>
class Opm::ElasticityLocalResidual< TypeTag > Calculation of (linear) elasticity model terms for the residual.
The linearized Biot model is solved where it is assumed that solid mechanics are governed by Hooke's law and conservation of linear momentum:
where \f&\sigma\f& is Cauchy stress tensor, \f&u\f& is displacement, \f&\epsilon(\cdot)\f& is thesymmetric gradient, \f&\mu\f& and \f&\lambda\f& are Lame's first (aka shear modulus) and second parameters, \f&\alpha\f& is the Biot-Willis parameter, \f&(p_f-p_0)\f& is fluid pressure difference wrt hydrostatic, and \f&f_u\f& are body forces.
The equations are discretized using two-point stress approximation following Boon et al. (2025), Solving Biot poroelasticity by coupling OPM Flow with the two-point stress approximation finite volume method, arXiv:2510.23432v1.
The resulting equations contain a volume term where only single-cell variables are used; face terms where variables across cell faces are calculated; boundary terms, similar to face terms, but cell faces are at the boundary; and source terms where coupling and potential body forces are calculated.
◆ computeBoundaryTerm()
template<class TypeTag >
template<class BoundaryConditionData >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeBoundaryTerm |
( |
Dune::FieldVector< Evaluation, numEq > & |
bndryTerm, |
|
|
const MaterialState & |
materialState, |
|
|
const BoundaryConditionData & |
bdyInfo, |
|
|
Problem & |
problem, |
|
|
unsigned |
globalIndex |
|
) |
| |
|
inlinestatic |
◆ computeBoundaryTermFixed()
template<class TypeTag >
template<class BoundaryConditionData >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeBoundaryTermFixed |
( |
Dune::FieldVector< Evaluation, numEq > & |
bndryTerm, |
|
|
const MaterialState & |
materialState, |
|
|
const BoundaryConditionData & |
bdyInfo, |
|
|
Problem & |
problem, |
|
|
unsigned |
globalIndex |
|
) |
| |
|
inlinestatic |
Calculate fixed displacement boundary condition in TPSA formulation.
- Parameters
-
| bndryTerm | Boundary term vector |
| materialState | Material state container |
| bdyInfo | Boundary condition info container |
| problem | Flow problem |
| globalIndex | Cell index |
- Note
- BCMECHType::NONE is a implemented as a specialization of BCMECHType::FIXED with displacement equal to zero on the boundary face
Referenced by Opm::ElasticityLocalResidual< TypeTag >::computeBoundaryTerm().
◆ computeBoundaryTermFree()
template<class TypeTag >
template<class BoundaryConditionData >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeBoundaryTermFree |
( |
Dune::FieldVector< Evaluation, numEq > & |
bndryTerm, |
|
|
const MaterialState & |
materialState, |
|
|
const BoundaryConditionData & |
bdyInfo, |
|
|
Problem & |
problem, |
|
|
unsigned |
globalIndex |
|
) |
| |
|
inlinestatic |
Calculate free (or zero traction) boundary condition in TPSA formulation.
- Parameters
-
| bndryTerm | Boundary term vector |
| materialState | Material state container |
| bdyInfo | Boundary condition info container |
| problem | Flow problem |
| globalIndex | Cell index |
- Note
- Free, or zero traction, BC is equivalent of having a spring at infinity where we have assumed all (primary) variables and parameters (e.g., shear modulus) are zero.
Referenced by Opm::ElasticityLocalResidual< TypeTag >::computeBoundaryTerm().
◆ computeFaceTerm()
template<class TypeTag >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeFaceTerm |
( |
Dune::FieldVector< Evaluation, numEq > & |
faceTerm, |
|
|
const MaterialState & |
materialStateIn, |
|
|
const MaterialState & |
materialStateEx, |
|
|
Problem & |
problem, |
|
|
const unsigned |
globalIndexIn, |
|
|
const unsigned |
globalIndexEx |
|
) |
| |
|
inlinestatic |
Calculate terms across cell faces in TPSA formulation.
- Parameters
-
| faceTerm | Face term vector |
| materialStateIn | Material state container of inside cell |
| materialStateEx | Material state container of outside cel |
| problem | Flow problem |
| globalIndexIn | Inside cell index |
| globalIndexEx | Outside cell index |
Material state, problem input and global index here might/should be merged in "IntensiveQuantity" and "NeighborInfo" containers as in BlackOilLocalResidualTPFA
◆ computeSourceTerm()
template<class TypeTag >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeSourceTerm |
( |
Dune::FieldVector< Evaluation, numEq > & |
sourceTerm, |
|
|
Problem & |
problem, |
|
|
unsigned |
globalSpaceIdex, |
|
|
unsigned |
timeIdx |
|
) |
| |
|
inlinestatic |
Calculate source term in TPSA formulation.
- Parameters
-
| sourceTerm | Source term vector |
| problem | Flow problem |
| globalSpaceIdex | Cell index |
| timeIdx | Time index |
◆ computeVolumeTerm()
template<class TypeTag >
template<class LhsEval >
| static void Opm::ElasticityLocalResidual< TypeTag >::computeVolumeTerm |
( |
Dune::FieldVector< LhsEval, numEq > & |
volTerm, |
|
|
const MaterialState & |
materialState, |
|
|
const Problem & |
problem, |
|
|
const unsigned |
globalIndex |
|
) |
| |
|
inlinestatic |
Calculate volume terms in TPSA formulation.
- Parameters
-
| volTerm | Volume term vector |
| materialState | Material state container |
| problem | Flow problem |
| globalIndex | Cell index |
Material state, problem input and global index here might/should be merged in an "IntensiveQuantity" container as in BlackOilLocalResidualTPFA
The documentation for this class was generated from the following file:
|