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>
66template <
class TypeTag>
74 using MaterialState = MaterialStateTPSA<Evaluation>;
75 using Toolbox = MathToolbox<Evaluation>;
77 enum { conti0EqIdx = Indices::conti0EqIdx };
78 enum { contiRotEqIdx = Indices::contiRotEqIdx };
79 enum { contiSolidPresEqIdx = Indices::contiSolidPresEqIdx };
80 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
94 template <
class LhsEval>
96 const MaterialState& materialState,
97 const Problem& problem,
98 const unsigned globalIndex)
104 const Scalar sModulus = problem.shearModulus(globalIndex);
105 for (
unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
106 volTerm[contiRotEqIdx + dirIdx] +=
107 Toolbox::template decay<LhsEval>(materialState.rotation(dirIdx))
112 const Scalar lame = problem.lame(globalIndex);
113 volTerm[contiSolidPresEqIdx] +=
114 Toolbox::template decay<LhsEval>(materialState.solidPressure())
132 const MaterialState& materialStateIn,
133 const MaterialState& materialStateEx,
135 const unsigned globalIndexIn,
136 const unsigned globalIndexEx)
142 const Scalar weightAvgIn = problem.weightAverage(globalIndexIn, globalIndexEx);
143 const Scalar weightAvgEx = problem.weightAverage(globalIndexEx, globalIndexIn);
144 const Scalar weightProd = problem.weightProduct(globalIndexIn, globalIndexEx);
145 const Scalar normDist = problem.normalDistance(globalIndexIn, globalIndexEx);
146 const auto& faceNormal = problem.cellFaceNormal(globalIndexIn, globalIndexEx);
149 const Scalar sModulusIn = problem.shearModulus(globalIndexIn);
150 const Scalar sModulusEx = problem.shearModulus(globalIndexEx);
151 const Scalar eff_sModulus = weightAvgIn * sModulusIn + weightAvgEx * sModulusEx;
154 const Scalar distRatio = 0.5 * weightProd / normDist;
157 const Evaluation& solidPIn = materialStateIn.solidPressure();
158 const Scalar solidPEx = decay<Scalar>(materialStateEx.solidPressure());
163 faceTerm[contiSolidPresEqIdx] +=
164 distRatio * eff_sModulus * (solidPIn - solidPEx);
171 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
174 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
176 unsigned dirIdxNeg = modNeg(dirIdx - 1);
177 unsigned dirIdxPos = modNeg(dirIdx + 1);
180 const Scalar faceNormalDir = faceNormal[dirIdx];
181 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
182 const Scalar faceNormalPos = faceNormal[dirIdxPos];
184 const Evaluation& dispIn = materialStateIn.displacement(dirIdx);
185 const Scalar dispEx = decay<Scalar>(materialStateEx.displacement(dirIdx));
187 const Evaluation& rotInNeg = materialStateIn.rotation(dirIdxNeg);
188 const Evaluation& rotInPos = materialStateIn.rotation(dirIdxPos);
189 const Scalar rotExNeg = decay<Scalar>(materialStateEx.rotation(dirIdxNeg));
190 const Scalar rotExPos = decay<Scalar>(materialStateEx.rotation(dirIdxPos));
192 faceTerm[conti0EqIdx + dirIdx] +=
193 2.0 * (eff_sModulus / normDist) * (dispIn - dispEx)
194 - weightAvgIn * (faceNormalNeg * rotInPos - faceNormalPos * rotInNeg)
195 - weightAvgEx * (faceNormalNeg * rotExPos - faceNormalPos * rotExNeg)
196 - faceNormalDir * (weightAvgIn * solidPIn + weightAvgEx * solidPEx);
199 const Evaluation& dispInNeg = materialStateIn.displacement(dirIdxNeg);
200 const Evaluation& dispInPos = materialStateIn.displacement(dirIdxPos);
201 const Scalar dispExNeg = decay<Scalar>(materialStateEx.displacement(dirIdxNeg));
202 const Scalar dispExPos = decay<Scalar>(materialStateEx.displacement(dirIdxPos));
204 faceTerm[contiRotEqIdx + dirIdx] +=
205 - weightAvgEx * (faceNormalNeg * dispInPos - faceNormalPos * dispInNeg)
206 - weightAvgIn * (faceNormalNeg * dispExPos - faceNormalPos * dispExNeg);
209 faceTerm[contiSolidPresEqIdx] +=
210 - faceNormalDir * (weightAvgEx * dispIn + weightAvgIn * dispEx);
223 template <
class BoundaryConditionData>
225 const MaterialState& materialState,
226 const BoundaryConditionData& bdyInfo,
228 unsigned globalIndex)
231 switch (bdyInfo.type) {
240 case BCMECHType::FIXED:
241 throw std::runtime_error(
"BCTYPE FIXED has not been implemented in TPSA");
242 case BCMECHType::FREE:
250 throw std::logic_error(
"Unknown boundary condition type " +
252 " in computeBoundaryFlux()." );
268 template <
class BoundaryConditionData>
270 const MaterialState& materialState,
271 const BoundaryConditionData& bdyInfo,
273 unsigned globalIndex)
283 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
284 const Scalar weightAvg = problem.weightAverageBoundary(globalIndex, bfIdx);
285 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
286 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
289 const Scalar eff_sModulus = problem.shearModulus(globalIndex);
292 const Evaluation& solidP = materialState.solidPressure();
299 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
302 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
304 unsigned dirIdxNeg = modNeg(dirIdx - 1);
305 unsigned dirIdxPos = modNeg(dirIdx + 1);
308 const Scalar faceNormalDir = faceNormal[dirIdx];
309 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
310 const Scalar faceNormalPos = faceNormal[dirIdxPos];
312 const Evaluation& disp = materialState.displacement(dirIdx);
314 const Evaluation& rotNeg = materialState.rotation(dirIdxNeg);
315 const Evaluation& rotPos = materialState.rotation(dirIdxPos);
317 bndryTerm[conti0EqIdx + dirIdx] +=
318 2.0 * (eff_sModulus / normDist) * disp
319 - weightAvg * (faceNormalNeg * rotPos - faceNormalPos * rotNeg)
320 - faceNormalDir * weightAvg * solidP;
336 template <
class BoundaryConditionData>
338 const MaterialState& materialState,
339 const BoundaryConditionData& bdyInfo,
341 unsigned globalIndex)
347 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
348 const Scalar weightAvg = 1.0;
349 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
350 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
352 const Scalar sModulus = problem.shearModulus(globalIndex);
357 const Evaluation& solidP = materialState.solidPressure();
358 bndryTerm[contiSolidPresEqIdx] +=
359 0.5 * (normDist / sModulus) * solidP;
366 auto modNeg = [](
int i) {
return ((i % 3) + 3) % 3; };
369 Evaluation dotProd = 0;
370 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
371 dotProd += faceNormal[dirIdx] * materialState.rotation(dirIdx);
375 for (
int dirIdx = 0; dirIdx < 3; ++dirIdx) {
377 unsigned dirIdxNeg = modNeg(dirIdx - 1);
378 unsigned dirIdxPos = modNeg(dirIdx + 1);
381 const Scalar faceNormalDir = faceNormal[dirIdx];
382 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
383 const Scalar faceNormalPos = faceNormal[dirIdxPos];
385 const Evaluation& dispNeg = materialState.displacement(dirIdxNeg);
386 const Evaluation& dispPos = materialState.displacement(dirIdxPos);
388 const Evaluation& rot = materialState.rotation(dirIdx);
390 bndryTerm[contiRotEqIdx + dirIdx] +=
391 - weightAvg * (faceNormalNeg * dispPos - faceNormalPos * dispNeg)
392 + 0.5 * (normDist / sModulus) * (dotProd * faceNormalDir - rot);
395 const Evaluation& disp = materialState.displacement(dirIdx);
397 bndryTerm[contiSolidPresEqIdx] +=
398 - faceNormalDir * weightAvg * disp;
412 unsigned globalSpaceIdex,
420 problem.tpsaSource(sourceTerm, globalSpaceIdex, timeIdx);
Calculation of (linear) elasticity model terms for the residual.
Definition: elasticitylocalresidualtpsa.hpp:68
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:131
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:269
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:337
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:95
static void computeSourceTerm(Dune::FieldVector< Evaluation, numEq > &sourceTerm, Problem &problem, unsigned globalSpaceIdex, unsigned timeIdx)
Calculate source term in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:410
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:224
@ 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)