25#ifndef ELASTICITY_LOCAL_RESIDUAL_TPSA_HPP
26#define ELASTICITY_LOCAL_RESIDUAL_TPSA_HPP
28#include <dune/common/fvector.hh>
30#include <opm/input/eclipse/Schedule/BCProp.hpp>
32#include <opm/material/common/MathToolbox.hpp>
33#include <opm/material/materialstates/MaterialStateTPSA.hpp>
67template <
class TypeTag>
75 using MaterialState = MaterialStateTPSA<Evaluation>;
76 using Toolbox = MathToolbox<Evaluation>;
78 enum { conti0EqIdx = Indices::conti0EqIdx };
79 enum { contiRotEqIdx = Indices::contiRotEqIdx };
80 enum { contiSolidPresEqIdx = Indices::contiSolidPresEqIdx };
81 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
95 template <
class LhsEval>
97 const MaterialState& materialState,
98 const Problem& problem,
99 const unsigned globalIndex)
105 const Scalar sModulus = problem.shearModulus(globalIndex);
106 for (
unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
107 volTerm[contiRotEqIdx + dirIdx] +=
108 Toolbox::template decay<LhsEval>(materialState.rotation(dirIdx))
113 const Scalar lame = problem.lame(globalIndex);
114 volTerm[contiSolidPresEqIdx] +=
115 Toolbox::template decay<LhsEval>(materialState.solidPressure())
133 const MaterialState& materialStateIn,
134 const MaterialState& materialStateEx,
136 const unsigned globalIndexIn,
137 const unsigned globalIndexEx)
143 const Scalar weightAvgIn = problem.weightAverage(globalIndexIn, globalIndexEx);
144 const Scalar weightAvgEx = problem.weightAverage(globalIndexEx, globalIndexIn);
145 const Scalar weightProd = problem.weightProduct(globalIndexIn, globalIndexEx);
146 const Scalar normDist = problem.normalDistance(globalIndexIn, globalIndexEx);
147 const auto& faceNormal = problem.cellFaceNormal(globalIndexIn, globalIndexEx);
150 const Scalar sModulusIn = problem.shearModulus(globalIndexIn);
151 const Scalar sModulusEx = problem.shearModulus(globalIndexEx);
152 const Scalar eff_sModulus = weightAvgIn * sModulusIn + weightAvgEx * sModulusEx;
155 const Scalar distRatio = 0.5 * weightProd / normDist;
158 const Evaluation& solidPIn = materialStateIn.solidPressure();
159 const Scalar solidPEx = decay<Scalar>(materialStateEx.solidPressure());
164 faceTerm[contiSolidPresEqIdx] +=
165 distRatio * eff_sModulus * (solidPIn - solidPEx);
172 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
175 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
177 unsigned dirIdxNeg = modNeg(dirIdx - 1);
178 unsigned dirIdxPos = modNeg(dirIdx + 1);
181 const Scalar faceNormalDir = faceNormal[dirIdx];
182 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
183 const Scalar faceNormalPos = faceNormal[dirIdxPos];
185 const Evaluation& dispIn = materialStateIn.displacement(dirIdx);
186 const Scalar dispEx = decay<Scalar>(materialStateEx.displacement(dirIdx));
188 const Evaluation& rotInNeg = materialStateIn.rotation(dirIdxNeg);
189 const Evaluation& rotInPos = materialStateIn.rotation(dirIdxPos);
190 const Scalar rotExNeg = decay<Scalar>(materialStateEx.rotation(dirIdxNeg));
191 const Scalar rotExPos = decay<Scalar>(materialStateEx.rotation(dirIdxPos));
193 faceTerm[conti0EqIdx + dirIdx] +=
194 2.0 * (eff_sModulus / normDist) * (dispIn - dispEx)
195 - weightAvgIn * (faceNormalNeg * rotInPos - faceNormalPos * rotInNeg)
196 - weightAvgEx * (faceNormalNeg * rotExPos - faceNormalPos * rotExNeg)
197 - faceNormalDir * (weightAvgIn * solidPIn + weightAvgEx * solidPEx);
200 const Evaluation& dispInNeg = materialStateIn.displacement(dirIdxNeg);
201 const Evaluation& dispInPos = materialStateIn.displacement(dirIdxPos);
202 const Scalar dispExNeg = decay<Scalar>(materialStateEx.displacement(dirIdxNeg));
203 const Scalar dispExPos = decay<Scalar>(materialStateEx.displacement(dirIdxPos));
205 faceTerm[contiRotEqIdx + dirIdx] +=
206 - weightAvgEx * (faceNormalNeg * dispInPos - faceNormalPos * dispInNeg)
207 - weightAvgIn * (faceNormalNeg * dispExPos - faceNormalPos * dispExNeg);
210 faceTerm[contiSolidPresEqIdx] +=
211 - faceNormalDir * (weightAvgEx * dispIn + weightAvgIn * dispEx);
224 template <
class BoundaryConditionData>
226 const MaterialState& materialState,
227 const BoundaryConditionData& bdyInfo,
229 unsigned globalIndex)
232 switch (bdyInfo.type) {
241 case BCMECHType::FIXED:
242 throw std::runtime_error(
"BCTYPE FIXED has not been implemented in TPSA");
243 case BCMECHType::FREE:
251 throw std::logic_error(
"Unknown boundary condition type " +
253 " in computeBoundaryFlux()." );
269 template <
class BoundaryConditionData>
271 const MaterialState& materialState,
272 const BoundaryConditionData& bdyInfo,
274 unsigned globalIndex)
284 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
285 const Scalar weightAvg = problem.weightAverageBoundary(globalIndex, bfIdx);
286 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
287 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
290 const Scalar eff_sModulus = problem.shearModulus(globalIndex);
293 const Evaluation& solidP = materialState.solidPressure();
300 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
303 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
305 unsigned dirIdxNeg = modNeg(dirIdx - 1);
306 unsigned dirIdxPos = modNeg(dirIdx + 1);
309 const Scalar faceNormalDir = faceNormal[dirIdx];
310 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
311 const Scalar faceNormalPos = faceNormal[dirIdxPos];
313 const Evaluation& disp = materialState.displacement(dirIdx);
315 const Evaluation& rotNeg = materialState.rotation(dirIdxNeg);
316 const Evaluation& rotPos = materialState.rotation(dirIdxPos);
318 bndryTerm[conti0EqIdx + dirIdx] +=
319 2.0 * (eff_sModulus / normDist) * disp
320 - weightAvg * (faceNormalNeg * rotPos - faceNormalPos * rotNeg)
321 - faceNormalDir * weightAvg * solidP;
337 template <
class BoundaryConditionData>
339 const MaterialState& materialState,
340 const BoundaryConditionData& bdyInfo,
342 unsigned globalIndex)
348 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
349 const Scalar weightAvg = 1.0;
350 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
351 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
353 const Scalar sModulus = problem.shearModulus(globalIndex);
358 const Evaluation& solidP = materialState.solidPressure();
359 bndryTerm[contiSolidPresEqIdx] +=
360 0.5 * (normDist / sModulus) * solidP;
367 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
370 Evaluation dotProd = 0;
371 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
372 dotProd += faceNormal[dirIdx] * materialState.rotation(dirIdx);
376 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
378 unsigned dirIdxNeg = modNeg(dirIdx - 1);
379 unsigned dirIdxPos = modNeg(dirIdx + 1);
382 const Scalar faceNormalDir = faceNormal[dirIdx];
383 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
384 const Scalar faceNormalPos = faceNormal[dirIdxPos];
386 const Evaluation& dispNeg = materialState.displacement(dirIdxNeg);
387 const Evaluation& dispPos = materialState.displacement(dirIdxPos);
389 const Evaluation& rot = materialState.rotation(dirIdx);
391 bndryTerm[contiRotEqIdx + dirIdx] +=
392 - weightAvg * (faceNormalNeg * dispPos - faceNormalPos * dispNeg)
393 + 0.5 * (normDist / sModulus) * (dotProd * faceNormalDir - rot);
396 const Evaluation& disp = materialState.displacement(dirIdx);
398 bndryTerm[contiSolidPresEqIdx] +=
399 - faceNormalDir * weightAvg * disp;
413 unsigned globalSpaceIdex,
421 problem.tpsaSource(sourceTerm, globalSpaceIdex, timeIdx);
Defines a type tags and some fundamental properties all models.
Calculation of (linear) elasticity model terms for the residual.
Definition: elasticitylocalresidualtpsa.hpp:69
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.
Definition: elasticitylocalresidualtpsa.hpp:132
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.
Definition: elasticitylocalresidualtpsa.hpp:270
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.
Definition: elasticitylocalresidualtpsa.hpp:338
static void computeVolumeTerm(Dune::FieldVector< LhsEval, numEq > &volTerm, const MaterialState &materialState, const Problem &problem, const unsigned globalIndex)
Calculate volume terms in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:96
static void computeSourceTerm(Dune::FieldVector< Evaluation, numEq > &sourceTerm, Problem &problem, unsigned globalSpaceIdex, unsigned timeIdx)
Calculate source term in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:411
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.
Definition: elasticitylocalresidualtpsa.hpp:225
Declare the properties used by the infrastructure code of the finite volume discretizations.
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)