25#ifndef FLOW_PROBLEM_TPSA_HPP
26#define FLOW_PROBLEM_TPSA_HPP
28#include <dune/common/fvector.hh>
30#include <dune/grid/common/rangegenerators.hh>
32#include <opm/common/OpmLog/OpmLog.hpp>
34#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/materialstates/MaterialStateTPSA.hpp>
39#include <opm/models/io/vtktpsamodule.hpp>
40#include <opm/models/tpsa/tpsabaseproperties.hpp>
41#include <opm/models/tpsa/tpsamodel.hpp>
44#include <opm/simulators/flow/FacePropertiesTPSA.hpp>
53#include <fmt/format.h>
61template <
class TypeTag>
79 enum { dimWorld = GridView::dimensionworld };
83 enum { numPhases = FluidSystem::numPhases };
85 enum { contiRotEqIdx = Indices::contiRotEqIdx };
86 enum { contiSolidPresEqIdx = Indices::contiSolidPresEqIdx };
87 enum { solidPres0Idx = Indices::solidPres0Idx };
90 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
92 using InitialMaterialState = MaterialStateTPSA<Scalar>;
93 using Toolbox = MathToolbox<Evaluation>;
104 : ParentType(simulator)
105 , faceProps_(simulator.vanguard().eclState(),
106 simulator.vanguard().gridView(),
107 simulator.vanguard().cartesianIndexMapper(),
108 simulator.vanguard().grid(),
109 simulator.vanguard().cellCentroids())
110 , geoMechModel_(simulator)
112 if constexpr(enableMech) {
117 const auto& mechSolver = simulator.vanguard().eclState().runspec().mechSolver();
118 if (!mechSolver.tpsa()) {
119 std::string msg =
"Simulator with Tpsa-geomechanics enabled compile time, but deck does not contain "
122 throw std::runtime_error(msg);
127 const auto& mechSolver = simulator.vanguard().eclState().runspec().mechSolver();
128 if (mechSolver.tpsa()) {
129 std::string msg =
"TPSA keyword in deck, but Tpsa-geomechanics disabled compile-time!";
131 throw std::runtime_error(msg);
145 GeomechModel::registerParameters();
163 faceProps_.finishInit();
169 geoMechModel_.finishInit();
183 auto& uCur = geoMechModel_.solution(0);
187 ElementContext elemCtx(this->simulator());
188 for (
const auto& elem : elements(this->gridView())) {
190 if (elem.partitionType() != Dune::InteriorEntity) {
195 elemCtx.updateStencil(elem);
196 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
197 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
198 auto& priVars = uCur[globalIdx];
199 priVars.assignNaive(initialMaterialState_[globalIdx]);
200 priVars.checkDefined();
205 geoMechModel_.syncOverlap();
208 for (
unsigned timeIdx = 1; timeIdx < historySize; ++timeIdx) {
209 geoMechModel_.solution(timeIdx) = geoMechModel_.solution(0);
213 geoMechModel_.updateMaterialState(0);
222 Scalar avgSmodulus = 0.0;
223 const auto& gridView = this->gridView();
224 ElementContext elemCtx(this->simulator());
225 for(
const auto& elem: elements(gridView, Dune::Partitions::interior)) {
226 elemCtx.updatePrimaryStencil(elem);
227 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
230 std::size_t numDof = this->model().numGridDof();
231 const auto& comm = this->simulator().vanguard().grid().comm();
232 avgSmodulus = comm.sum(avgSmodulus);
233 Scalar numTotalDof = comm.sum(numDof);
234 avgSmodulus /= numTotalDof;
235 avgSmodulus = std::sqrt(avgSmodulus);
237 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
238 if (eqIdx < contiRotEqIdx) {
239 geoMechModel_.setEqWeight(eqIdx, 1 / avgSmodulus);
242 geoMechModel_.setEqWeight(eqIdx, avgSmodulus);
257 if (this->nonTrivialBoundaryConditions()) {
258 geoMechModel_.linearizer().updateBoundaryConditionData();
274 std::pair<BCMECHType, Dune::FieldVector<Evaluation, 3>>
278 if (!this->nonTrivialBoundaryConditions()) {
279 return { BCMECHType::NONE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
283 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
284 const auto& schedule = this->simulator().vanguard().schedule();
285 if (this->bcindex_(dir)[globalSpaceIdx] == 0 || schedule[this->episodeIndex()].bcprop.size() == 0) {
286 return { BCMECHType::NONE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
290 const auto& bc = schedule[this->episodeIndex()].bcprop[this->bcindex_(dir)[globalSpaceIdx]];
291 if (bc.bcmechtype == BCMECHType::FREE) {
292 return { BCMECHType::FREE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
295 return { bc.bcmechtype, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
308 void tpsaSource(Dune::FieldVector<Evaluation, numEq>& sourceTerm,
309 unsigned globalSpaceIdx,
315 const auto biot = this->
biotCoeff(globalSpaceIdx);
316 const auto lameParam = this->
lame(globalSpaceIdx);
318 const auto& iq = this->model().intensiveQuantities(globalSpaceIdx, timeIdx);
319 const auto& fs = iq.fluidState();
320 const auto pres = decay<Scalar>(fs.pressure(this->refPressurePhaseIdx_()));
321 const auto initPres = this->initialFluidState(globalSpaceIdx).pressure(this->refPressurePhaseIdx_());
323 auto sourceFromFlow = -biot / lameParam * (pres - initPres);
324 sourceTerm[contiSolidPresEqIdx] += sourceFromFlow;
339 assert (timeIdx <= historySize);
340 const auto solidPres = (timeIdx == 0) ?
341 decay<Scalar>( geoMechModel_.materialState(elementIdx, timeIdx).solidPressure()) :
342 geoMechModel_.solution(timeIdx)[elementIdx][solidPres0Idx];
343 const auto biot = this->
biotCoeff(elementIdx);
344 const auto lameParam = this->
lame(elementIdx);
346 return biot / lameParam * solidPres;
361 return faceProps_.weightAverage(globalElemIdxIn, globalElemIdxOut);
373 return faceProps_.weightAverageBoundary(globalElemIdxIn, boundaryFaceIdx);
383 Scalar
weightProduct(
unsigned globalElemIdxIn,
unsigned globalElemIdxOut)
const
385 return faceProps_.weightProduct(globalElemIdxIn, globalElemIdxOut);
397 return faceProps_.normalDistance(globalElemIdxIn, globalElemIdxOut);
409 return faceProps_.normalDistanceBoundary(globalElemIdxIn, boundaryFaceIdx);
421 return faceProps_.cellFaceNormal(globalElemIdxIn, globalElemIdxOut);
433 return faceProps_.cellFaceNormalBoundary(globalElemIdxIn, boundaryFaceIdx);
444 return faceProps_.shearModulus(globalElemIdx);
454 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
455 return mechSolver.laggedScheme();
465 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
466 return mechSolver.fixedStressScheme();
476 return geoMechModel_;
486 return geoMechModel_;
496 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
497 return std::make_pair(mechSolver.fixedStressMinIter(), mechSolver.fixedStressMaxIter());
514 std::size_t numDof = this->model().numGridDof();
515 initialMaterialState_.resize(numDof);
516 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
517 auto& dofMaterialState = initialMaterialState_[dofIdx];
518 for (
unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
519 dofMaterialState.setDisplacement(dirIdx, 0.0);
520 dofMaterialState.setRotation(dirIdx, 0.0);
522 dofMaterialState.setSolidPressure(0.0);
527 FaceProperties faceProps_;
528 GeomechModel geoMechModel_;
530 std::vector<Scalar> biotcoeff_;
531 std::vector<InitialMaterialState> initialMaterialState_;
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition CollectDataOnIORank.hpp:49
Cell face properties needed in TPSA equation calculations.
Definition FacePropertiesTPSA.hpp:48
FlowProblemBlackoil(Simulator &simulator)
Definition FlowProblemBlackoil.hpp:187
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition FlowProblemBlackoil.hpp:303
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition FlowProblemBlackoil.hpp:614
void beginTimeStep() override
Called by the simulator before each time integration.
Definition FlowProblemBlackoil.hpp:294
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblemBlackoil.hpp:173
Scalar shearModulus(unsigned globalElemIdx) const
Direct access to shear modulus in an element.
Definition FlowProblemTPSA.hpp:442
Scalar weightAverageBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to normal distance at the boundary.
Definition FlowProblemTPSA.hpp:371
static void registerParameters()
Register runtime parameters.
Definition FlowProblemTPSA.hpp:139
const DimVector & cellFaceNormalBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to face normal at the boundary.
Definition FlowProblemTPSA.hpp:431
std::pair< BCMECHType, Dune::FieldVector< Evaluation, 3 > > mechBoundaryCondition(const unsigned int globalSpaceIdx, const int directionId)
Organize mechanics boundary conditions.
Definition FlowProblemTPSA.hpp:275
void initialSolutionApplied() override
Set initial solution for the problem.
Definition FlowProblemTPSA.hpp:177
std::pair< int, int > fixedStressParameters() const
Get fixed-stress iteration parameters.
Definition FlowProblemTPSA.hpp:494
void beginTimeStep() override
Called by the simulator before each time integration.
Definition FlowProblemTPSA.hpp:250
DimVector cellFaceNormal(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
Direct access to face normal between two elements.
Definition FlowProblemTPSA.hpp:419
bool fixedStressScheme() const
Flow-TPSA fixed-stress coupling scheme activated?
Definition FlowProblemTPSA.hpp:463
Scalar normalDistance(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
Direct access to normal distance between two elements.
Definition FlowProblemTPSA.hpp:395
void computeAndSetEqWeights_()
Compute weights to rescale the TPSA equations.
Definition FlowProblemTPSA.hpp:219
const GeomechModel & geoMechModel() const
Get TPSA model.
Definition FlowProblemTPSA.hpp:474
Scalar weightAverage(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
Direct access to average (half-)weight at interface between two elements.
Definition FlowProblemTPSA.hpp:359
void tpsaSource(Dune::FieldVector< Evaluation, numEq > &sourceTerm, unsigned globalSpaceIdx, unsigned timeIdx)
Set mechanics source term, in particular coupling terms.
Definition FlowProblemTPSA.hpp:308
GeomechModel & geoMechModel()
Get TPSA model.
Definition FlowProblemTPSA.hpp:484
FlowProblemTPSA(Simulator &simulator)
Constructor.
Definition FlowProblemTPSA.hpp:103
Scalar rockMechPoroChange(unsigned elementIdx, unsigned timeIdx) const
Pore volume change due to geomechanics.
Definition FlowProblemTPSA.hpp:336
void finishInit()
Initialize the problem.
Definition FlowProblemTPSA.hpp:154
bool laggedScheme() const
Flow-TPSA lagged coupling scheme activated?
Definition FlowProblemTPSA.hpp:452
Scalar normalDistanceBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to normal distance at the boundary.
Definition FlowProblemTPSA.hpp:407
Scalar weightProduct(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
Direct access to product of weights at interface between two elements.
Definition FlowProblemTPSA.hpp:383
void readInitalConditionsTPSA_()
Read initial conditions and generate material state for TPSA model.
Definition FlowProblemTPSA.hpp:507
Scalar lame(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:735
Scalar biotCoeff(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:745
VTK output module for TPSA quantities.
Definition vtktpsamodule.hpp:47
static void registerParameters()
Register runtime parameters.
Definition vtktpsamodule.hpp:81
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
This file provides the infrastructure to retrieve run-time parameters.