30#ifndef OPM_FLOW_PROBLEM_HPP
31#define OPM_FLOW_PROBLEM_HPP
33#include <dune/common/version.hh>
34#include <dune/common/fvector.hh>
35#include <dune/common/fmatrix.hh>
37#include <opm/common/utility/TimeService.hpp>
39#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
40#include <opm/input/eclipse/Schedule/Schedule.hpp>
41#include <opm/input/eclipse/Units/Units.hpp>
43#include <opm/material/common/ConditionalStorage.hpp>
44#include <opm/material/common/Valgrind.hpp>
45#include <opm/material/densead/Evaluation.hpp>
46#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
47#include <opm/material/thermal/EclThermalLawManager.hpp>
53#include <opm/output/eclipse/EclipseIO.hpp>
62#include <opm/simulators/flow/FlowUtils.hpp>
65#include <opm/simulators/timestepping/AdaptiveTimeStepping.hpp>
66#include <opm/simulators/timestepping/SimulatorReport.hpp>
68#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
69#include <opm/simulators/utils/ParallelSerialization.hpp>
70#include <opm/simulators/utils/satfunc/RelpermDiagnostics.hpp>
72#include <opm/utility/CopyablePtr.hpp>
90template <
class TypeTag>
93 GetPropType<TypeTag, Properties::FluidSystem>>
111 enum { dim = GridView::dimension };
112 enum { dimWorld = GridView::dimensionworld };
116 enum { numPhases = FluidSystem::numPhases };
117 enum { numComponents = FluidSystem::numComponents };
128 enum { enableMICP = Indices::enableMICP };
136 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
137 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
138 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
142 enum { gasCompIdx = FluidSystem::gasCompIdx };
143 enum { oilCompIdx = FluidSystem::oilCompIdx };
144 enum { waterCompIdx = FluidSystem::waterCompIdx };
149 using Element =
typename GridView::template Codim<0>::Entity;
153 using MaterialLawParams =
typename EclMaterialLawManager::MaterialLawParams;
154 using SolidEnergyLawParams =
typename EclThermalLawManager::SolidEnergyLawParams;
155 using ThermalConductionLawParams =
typename EclThermalLawManager::ThermalConductionLawParams;
164 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
168 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag>>;
183 ParentType::registerParameters();
198 const std::string&)>
addKey,
206 return detail::eclPositionalParameter(
addKey,
217 : ParentType(simulator)
218 ,
BaseType(simulator.vanguard().eclState(),
219 simulator.vanguard().schedule(),
220 simulator.vanguard().gridView())
221 , transmissibilities_(simulator.vanguard().eclState(),
222 simulator.vanguard().gridView(),
223 simulator.vanguard().cartesianIndexMapper(),
224 simulator.vanguard().grid(),
225 simulator.vanguard().cellCentroids(),
229 , wellModel_(simulator)
230 , aquiferModel_(simulator)
231 , pffDofData_(simulator.gridView(),
this->elementMapper())
232 , tracerModel_(simulator)
234 if (! Parameters::Get<Parameters::CheckSatfuncConsistency>()) {
240 simulator.vanguard().levelCartesianIndexMapper());
246 void prefetch(
const Element&
elem)
const
247 { this->pffDofData_.prefetch(
elem); }
260 template <
class Restarter>
267 wellModel_.deserialize(
res);
270 aquiferModel_.deserialize(
res);
279 template <
class Restarter>
282 wellModel_.serialize(
res);
284 aquiferModel_.serialize(
res);
287 int episodeIndex()
const
289 return std::max(this->simulator().episodeIndex(), 0);
299 auto& simulator = this->simulator();
301 auto& eclState = simulator.vanguard().eclState();
302 const auto& schedule = simulator.vanguard().schedule();
303 const auto& events = schedule[
episodeIdx].events();
305 if (
episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
313 const auto&
cc = simulator.vanguard().grid().comm();
314 eclState.apply_schedule_keywords(
miniDeck );
315 eclBroadcast(
cc, eclState.getTransMult() );
318 std::function<
unsigned int(
unsigned int)>
equilGridToGrid = [&simulator](
unsigned int i) {
319 return simulator.vanguard().gridEquilIdxToGridIdx(i);
323 using TransUpdateQuantities =
typename Vanguard::TransmissibilityType::TransUpdateQuantities;
324 transmissibilities_.update(
true, TransUpdateQuantities::All,
equilGridToGrid);
325 this->referencePorosity_[1] = this->referencePorosity_[0];
326 updateReferencePorosity_();
328 this->model().linearizer().updateDiscretizationParameters();
331 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
334 wellModel_.beginEpisode();
337 aquiferModel_.beginEpisode();
340 Scalar
dt = limitNextTimeStepSize_(simulator.episodeLength());
345 dt = std::min(
dt, this->initialTimeStepSize_);
346 simulator.setTimeStepSize(
dt);
356 const int timeStepSize = this->simulator().timeStepSize();
358 this->beginTimeStep_(enableExperiments,
360 this->simulator().timeStepIndex(),
361 this->simulator().startTime(),
362 this->simulator().time(),
364 this->simulator().endTime());
372 if (nonTrivialBoundaryConditions()) {
373 this->model().linearizer().updateBoundaryConditionData();
376 wellModel_.beginTimeStep();
377 aquiferModel_.beginTimeStep();
378 tracerModel_.beginTimeStep();
387 wellModel_.beginIteration();
388 aquiferModel_.beginIteration();
397 wellModel_.endIteration();
398 aquiferModel_.endIteration();
415 const int rank = this->simulator().gridView().comm().rank();
417 std::cout <<
"checking conservativeness of solution\n";
420 this->model().checkConservativeness(-1,
true);
422 std::cout <<
"solution is sufficiently conservative\n";
427 auto& simulator = this->simulator();
428 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
430 this->wellModel_.endTimeStep();
431 this->aquiferModel_.endTimeStep();
432 this->tracerModel_.endTimeStep();
435 this->model().linearizer().updateFlowsInfo();
437 if (this->enableDriftCompensation_) {
440 const auto& residual = this->model().linearizer().residual();
442 for (
unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
443 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
460 this->wellModel_.endEpisode();
461 this->aquiferModel_.endEpisode();
463 const auto& schedule = this->simulator().vanguard().schedule();
466 if (
episodeIdx + 1 >=
static_cast<int>(schedule.size()) - 1) {
467 this->simulator().setFinished(
true);
472 this->simulator().startNextEpisode(schedule.stepLength(
episodeIdx + 1));
483 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
487 ParentType::writeOutput(
verbose);
494 template <
class Context>
515 template <
class Context>
521 return pffDofData_.get(context.element(),
toDofLocalIdx).transmissibility;
535 template <
class Context>
541 return *pffDofData_.get(context.element(),
toDofLocalIdx).diffusivity;
573 template <
class Context>
577 unsigned elemIdx = elemCtx.globalSpaceIndex(0, 0);
578 return transmissibilities_.transmissibilityBoundary(elemIdx,
boundaryFaceIdx);
603 template <
class Context>
608 const auto&
face = context.stencil(
timeIdx).interiorFace(faceIdx);
610 return *pffDofData_.get(context.element(),
toDofLocalIdx).thermalHalfTransIn;
616 template <
class Context>
621 const auto&
face = context.stencil(
timeIdx).interiorFace(faceIdx);
623 return *pffDofData_.get(context.element(),
toDofLocalIdx).thermalHalfTransOut;
629 template <
class Context>
633 unsigned elemIdx = elemCtx.globalSpaceIndex(0, 0);
634 return transmissibilities_.thermalHalfTransBoundary(elemIdx,
boundaryFaceIdx);
641 {
return transmissibilities_; }
645 {
return tracerModel_; }
647 TracerModel& tracerModel()
648 {
return tracerModel_; }
658 template <
class Context>
671 template <
class Context>
686 return this->simulator().vanguard().cellCenterDepth(
globalSpaceIdx);
692 template <
class Context>
702 template <
class Context>
714 const auto&
rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
716 return asImp_().initialFluidState(
globalSpaceIdx).pressure(refPressurePhaseIdx_());
719 if (this->rockParams_.empty())
723 if (!this->rockTableIdx_.empty()) {
726 return this->rockParams_[
tableIdx].referencePressure;
733 template <
class Context>
743 return materialLawManager_->materialLawParams(globalDofIdx);
748 return materialLawManager_->materialLawParams(globalDofIdx,
facedir);
754 template <
class Context>
755 const SolidEnergyLawParams&
767 template <
class Context>
768 const ThermalConductionLawParams &
772 return thermalLawManager_->thermalConductionLawParams(
globalSpaceIdx);
782 {
return materialLawManager_; }
784 template <
class FluidState,
class ...Args>
786 std::array<Evaluation,numPhases> &mobility,
787 DirectionalMobilityPtr &
dirMob,
788 FluidState &fluidState,
791 using ContainerT = std::array<Evaluation, numPhases>;
798 Valgrind::CheckDefined(mobility);
800 if (materialLawManager_->hasDirectionalRelperms()
801 || materialLawManager_->hasDirectionalImbnum())
803 using Dir = FaceDir::DirEnum;
804 constexpr int ndim = 3;
805 dirMob = std::make_unique<DirectionalMobility<TypeTag>>();
806 Dir
facedirs[
ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
807 for (
int i = 0; i<
ndim; i++) {
819 {
return materialLawManager_; }
825 template <
class Context>
833 template <
class Context>
841 template <
class Context>
849 template <
class Context>
858 template <
class Context>
866 {
return this->simulator().vanguard().caseName(); }
871 template <
class Context>
877 return asImp_().initialFluidState(globalDofIdx).temperature(0);
881 Scalar
temperature(
unsigned globalDofIdx,
unsigned )
const
885 return asImp_().initialFluidState(globalDofIdx).temperature(0);
888 const SolidEnergyLawParams&
892 return this->thermalLawManager_->solidEnergyLawParams(
globalSpaceIdx);
894 const ThermalConductionLawParams &
898 return this->thermalLawManager_->thermalConductionLawParams(
globalSpaceIdx);
912 if (!this->vapparsActive(this->episodeIndex()))
915 return this->maxOilSaturation_[globalDofIdx];
929 if (!this->vapparsActive(this->episodeIndex()))
932 this->maxOilSaturation_[globalDofIdx] = value;
941 this->model().invalidateAndUpdateIntensiveQuantities(0);
947 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
948 this->model().invalidateAndUpdateIntensiveQuantities(1);
956 aquiferModel_.initialSolutionApplied();
961 this->model().invalidateAndUpdateIntensiveQuantities(0);
970 template <
class Context>
972 const Context& context,
976 const unsigned globalDofIdx = context.globalSpaceIndex(
spaceIdx,
timeIdx);
980 void source(RateVector& rate,
981 unsigned globalDofIdx,
988 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
993 rate[
eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
995 Valgrind::CheckDefined(rate[
eqIdx]);
1000 addToSourceDense(rate, globalDofIdx,
timeIdx);
1003 virtual void addToSourceDense(RateVector& rate,
1004 unsigned globalDofIdx,
1013 {
return wellModel_; }
1016 {
return wellModel_; }
1018 const AquiferModel& aquiferModel()
const
1019 {
return aquiferModel_; }
1021 AquiferModel& mutableAquiferModel()
1022 {
return aquiferModel_; }
1024 bool nonTrivialBoundaryConditions()
const
1025 {
return nonTrivialBoundaryConditions_; }
1037 if (this->nextTimeStepSize_ > 0.0)
1038 return this->nextTimeStepSize_;
1040 const auto& simulator = this->simulator();
1045 return this->initialTimeStepSize_;
1050 const auto& newtonMethod = this->model().newtonMethod();
1051 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1059 template <
class LhsEval>
1063 if (this->rockCompPoroMult_.empty() &&
this->rockCompPoroMultWc_.empty())
1067 if (!this->rockTableIdx_.empty())
1070 const auto& fs = intQuants.fluidState();
1072 const auto&
rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1073 if (!this->minRefPressure_.empty())
1079 if (!this->overburdenPressure_.empty())
1086 if (!this->rockCompPoroMult_.empty()) {
1091 assert(!this->rockCompPoroMultWc_.empty());
1103 template <
class LhsEval>
1106 const bool implicit = !this->explicitRockCompaction_;
1108 : this->simulator().problem().getRockCompTransMultVal(
elementIdx);
1115 template <
class LhsEval>
1120 const bool implicit = !this->explicitRockCompaction_;
1122 : this->simulator().problem().getRockCompTransMultVal(
elementIdx);
1131 if (!nonTrivialBoundaryConditions_) {
1132 return { BCType::NONE, RateVector(0.0) };
1134 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(
directionId);
1135 const auto& schedule = this->simulator().vanguard().schedule();
1137 return { BCType::NONE, RateVector(0.0) };
1139 if (schedule[this->episodeIndex()].
bcprop.size() == 0) {
1140 return { BCType::NONE, RateVector(0.0) };
1142 const auto&
bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[
globalSpaceIdx]];
1143 if (
bc.bctype!=BCType::RATE) {
1144 return {
bc.bctype, RateVector(0.0) };
1147 RateVector rate = 0.0;
1148 switch (
bc.component) {
1149 case BCComponent::OIL:
1150 rate[FluidSystem::canonicalToActiveCompIdx(oilCompIdx)] =
bc.rate;
1152 case BCComponent::GAS:
1153 rate[FluidSystem::canonicalToActiveCompIdx(gasCompIdx)] =
bc.rate;
1155 case BCComponent::WATER:
1156 rate[FluidSystem::canonicalToActiveCompIdx(waterCompIdx)] =
bc.rate;
1158 case BCComponent::SOLVENT:
1159 this->handleSolventBC(
bc, rate);
1161 case BCComponent::POLYMER:
1162 this->handlePolymerBC(
bc, rate);
1164 case BCComponent::MICR:
1165 this->handleMicrBC(
bc, rate);
1167 case BCComponent::OXYG:
1168 this->handleOxygBC(
bc, rate);
1170 case BCComponent::UREA:
1171 this->handleUreaBC(
bc, rate);
1173 case BCComponent::NONE:
1174 throw std::logic_error(
"you need to specify the component when RATE type is set in BC");
1178 return {
bc.bctype, rate};
1182 template<
class Serializer>
1194 Implementation& asImp_()
1195 {
return *
static_cast<Implementation *
>(
this); }
1197 const Implementation& asImp_()
const
1198 {
return *
static_cast<const Implementation *
>(
this); }
1201 template<
class UpdateFunc>
1202 void updateProperty_(
const std::string&
failureMsg,
1206 const auto& model = this->simulator().model();
1207 const auto& primaryVars = model.solution(0);
1208 const auto& vanguard = this->simulator().vanguard();
1209 std::size_t numGridDof = primaryVars.size();
1210 OPM_BEGIN_PARALLEL_TRY_CATCH();
1212#pragma omp parallel for
1214 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1215 const auto&
iq = *model.cachedIntensiveQuantities(dofIdx, 0);
1218 OPM_END_PARALLEL_TRY_CATCH(
failureMsg, vanguard.grid().comm());
1221 bool updateMaxOilSaturation_()
1228 this->updateProperty_(
"FlowProblem::updateMaxOilSaturation_() failed:",
1239 bool updateMaxOilSaturation_(
unsigned compressedDofIdx,
const IntensiveQuantities&
iq)
1242 const auto& fs =
iq.fluidState();
1243 const Scalar So =
decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1244 auto&
mos = this->maxOilSaturation_;
1253 bool updateMaxWaterSaturation_()
1257 if (this->maxWaterSaturation_.empty())
1260 this->maxWaterSaturation_[1] = this->maxWaterSaturation_[0];
1261 this->updateProperty_(
"FlowProblem::updateMaxWaterSaturation_() failed:",
1270 bool updateMaxWaterSaturation_(
unsigned compressedDofIdx,
const IntensiveQuantities&
iq)
1273 const auto& fs =
iq.fluidState();
1274 const Scalar Sw =
decay<Scalar>(fs.saturation(waterPhaseIdx));
1275 auto&
mow = this->maxWaterSaturation_;
1284 bool updateMinPressure_()
1288 if (this->minRefPressure_.empty())
1291 this->updateProperty_(
"FlowProblem::updateMinPressure_() failed:",
1301 const auto& fs =
iq.fluidState();
1316 std::function<std::vector<double>(
const FieldPropsManager&,
const std::string&)>
1317 fieldPropDoubleOnLeafAssigner_()
1319 const auto&
lookup = this->lookUpData_;
1330 template<
typename IntType>
1331 std::function<std::vector<IntType>(
const FieldPropsManager&,
const std::string&,
bool)>
1332 fieldPropIntTypeOnLeafAssigner_()
1334 const auto&
lookup = this->lookUpData_;
1341 void readMaterialParameters_()
1344 const auto& simulator = this->simulator();
1345 const auto& vanguard = simulator.vanguard();
1346 const auto& eclState = vanguard.eclState();
1349 OPM_BEGIN_PARALLEL_TRY_CATCH();
1350 this->updatePvtnum_();
1351 this->updateSatnum_();
1354 this->updateMiscnum_();
1356 this->updatePlmixnum_();
1358 OPM_END_PARALLEL_TRY_CATCH(
"Invalid region numbers: ", vanguard.gridView().comm());
1361 updateReferencePorosity_();
1362 this->referencePorosity_[1] = this->referencePorosity_[0];
1367 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1368 materialLawManager_->initFromState(eclState);
1369 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1371 this-> lookupIdxOnLevelZeroAssigner_());
1375 void readThermalParameters_()
1377 if constexpr (enableEnergy)
1379 const auto& simulator = this->simulator();
1380 const auto& vanguard = simulator.vanguard();
1381 const auto& eclState = vanguard.eclState();
1384 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1385 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1386 this-> fieldPropDoubleOnLeafAssigner_(),
1391 void updateReferencePorosity_()
1393 const auto& simulator = this->simulator();
1394 const auto& vanguard = simulator.vanguard();
1395 const auto& eclState = vanguard.eclState();
1397 std::size_t numDof = this->model().numGridDof();
1399 this->referencePorosity_[0].resize(numDof);
1401 const auto&
fp = eclState.fieldProps();
1402 const std::vector<double>
porvData =
this -> fieldPropDoubleOnLeafAssigner_()(
fp,
"PORV");
1403 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1404 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1410 Scalar dofVolume = simulator.model().dofTotalVolume(
sfcdofIdx);
1416 virtual void readInitialCondition_()
1419 const auto& simulator = this->simulator();
1420 const auto& vanguard = simulator.vanguard();
1421 const auto& eclState = vanguard.eclState();
1423 if (eclState.getInitConfig().hasEquil())
1424 readEquilInitialCondition_();
1426 readExplicitInitialCondition_();
1429 std::size_t
numElems = this->model().numGridDof();
1430 for (std::size_t elemIdx = 0; elemIdx <
numElems; ++elemIdx) {
1431 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1432 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1433 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1434 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1435 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1436 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1437 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1441 virtual void readEquilInitialCondition_() = 0;
1442 virtual void readExplicitInitialCondition_() = 0;
1445 bool updateHysteresis_()
1447 if (!materialLawManager_->enableHysteresis())
1452 this->updateProperty_(
"FlowProblem::updateHysteresis_() failed:",
1469 Scalar getRockCompTransMultVal(std::size_t dofIdx)
const
1471 if (this->rockCompTransMultVal_.empty())
1474 return this->rockCompTransMultVal_[dofIdx];
1484 Scalar transmissibility;
1488 void updatePffDofData_()
1492 const Stencil& stencil,
1496 const auto& elementMapper = this->model().elementMapper();
1503 if constexpr (enableEnergy) {
1507 if constexpr (enableDiffusion)
1509 if (enableDispersion)
1514 pffDofData_.update(
distFn);
1519 void readBoundaryConditions_()
1521 const auto& simulator = this->simulator();
1522 const auto& vanguard = simulator.vanguard();
1523 const auto&
bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1525 nonTrivialBoundaryConditions_ =
true;
1527 std::size_t
numCartDof = vanguard.cartesianSize();
1528 unsigned numElems = vanguard.gridView().size(0);
1531 for (
unsigned elemIdx = 0; elemIdx <
numElems; ++elemIdx)
1536 &vanguard](
const auto&
bcface,
1542 std::array<int, 3> tmp = {i,j,
k};
1551 std::vector<int>& data = bcindex_(
bcface.dir);
1552 const int index =
bcface.index;
1554 [&data,index](
int elemIdx)
1555 { data[elemIdx] = index; });
1562 Scalar limitNextTimeStepSize_(Scalar
dtNext)
const
1564 if constexpr (enableExperiments) {
1565 const auto& simulator = this->simulator();
1566 const auto& schedule = simulator.vanguard().schedule();
1570 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1572 if (this->enableTuning_) {
1574 maxTimeStepSize =
tuning.TSMAXZ;
1580 simulator.episodeStartTime() + simulator.episodeLength()
1581 - (simulator.startTime() + simulator.time());
1591 if (simulator.episodeStarts()) {
1594 const auto& events = simulator.vanguard().schedule()[
reportStepIdx].events();
1596 events.hasEvent(ScheduleEvents::NEW_WELL)
1597 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1598 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1599 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1601 dtNext = std::min(
dtNext, this->maxTimeStepAfterWellEvent_);
1608 int refPressurePhaseIdx_()
const {
1609 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1612 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1616 return waterPhaseIdx;
1620 void updateRockCompTransMultVal_()
1622 const auto& model = this->simulator().model();
1623 std::size_t numGridDof = this->model().numGridDof();
1624 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1626 const auto&
iq = *model.cachedIntensiveQuantities(
elementIdx, 0);
1637 template <
class LhsEval>
1641 if (this->rockCompTransMult_.empty() &&
this->rockCompTransMultWc_.empty())
1645 if (!this->rockTableIdx_.empty())
1648 const auto& fs = intQuants.fluidState();
1650 const auto&
rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1651 if (!this->minRefPressure_.empty())
1657 if (!this->overburdenPressure_.empty())
1664 if (!this->rockCompTransMult_.empty())
1668 assert(!this->rockCompTransMultWc_.empty());
1675 typename Vanguard::TransmissibilityType transmissibilities_;
1677 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1678 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1680 GlobalEqVector drift_;
1682 WellModel wellModel_;
1683 AquiferModel aquiferModel_;
1691 std::array<std::vector<T>,6> data;
1693 void resize(std::size_t size, T
defVal)
1695 for (
auto&
d : data)
1699 const std::vector<T>& operator()(FaceDir::DirEnum dir)
const
1701 if (dir == FaceDir::DirEnum::Unknown)
1702 throw std::runtime_error(
"Tried to access BC data for the 'Unknown' direction");
1704 int div =
static_cast<int>(dir);
1705 while ((
div /= 2) >= 1)
1711 std::vector<T>& operator()(FaceDir::DirEnum dir)
1713 return const_cast<std::vector<T>&
>(std::as_const(*
this)(dir));
1717 virtual void handleSolventBC(
const BCProp::BCFace&, RateVector&)
const = 0;
1719 virtual void handlePolymerBC(
const BCProp::BCFace&, RateVector&)
const = 0;
1721 virtual void handleMicrBC(
const BCProp::BCFace&, RateVector&)
const = 0;
1723 virtual void handleOxygBC(
const BCProp::BCFace&, RateVector&)
const = 0;
1725 virtual void handleUreaBC(
const BCProp::BCFace&, RateVector&)
const = 0;
1728 bool nonTrivialBoundaryConditions_ =
false;
1729 bool first_step_ =
true;
1735 return this->simulator().episodeWillBeOver();
Helper class for grid instantiation of ECL file-format using problems.
This is a "dummy" gradient calculator which does not do anything.
Collects necessary output values and pass it to opm-common's ECL output.
Computes the initial condition based on the EQUIL keyword from ECL.
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
A class which handles tracers as specified in by ECL.
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition FlowGenericProblem.hpp:61
Scalar maxPolymerAdsorption(unsigned elemIdx) const
Returns the max polymer adsorption value.
Definition FlowGenericProblem_impl.hpp:744
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition FlowGenericProblem_impl.hpp:703
static std::string briefDescription()
Returns a human readable description of the problem for the help message.
Definition FlowGenericProblem_impl.hpp:129
Scalar porosity(unsigned globalSpaceIdx, unsigned timeIdx) const
Direct indexed access to the porosity.
Definition FlowGenericProblem_impl.hpp:332
Scalar rockCompressibility(unsigned globalSpaceIdx) const
Direct access to rock compressibility.
Definition FlowGenericProblem_impl.hpp:317
unsigned miscnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant MISC region given a cell index.
Definition FlowGenericProblem_impl.hpp:723
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition FlowGenericProblem_impl.hpp:713
unsigned plmixnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant PLMIXNUM (for polymer module) region given a cell index.
Definition FlowGenericProblem_impl.hpp:733
bool shouldWriteOutput() const
Always returns true.
Definition FlowGenericProblem.hpp:277
static std::string helpPreamble(int, const char **argv)
Returns the string that is printed before the list of command line parameters in the help message.
Definition FlowGenericProblem_impl.hpp:114
bool shouldWriteRestartFile() const
Returns true if an eWoms restart file should be written to disk.
Definition FlowGenericProblem.hpp:286
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition FlowProblem.hpp:94
virtual bool episodeWillBeOver() const
Whether or not the current episode will end at the end of the current time step.
Definition FlowProblem.hpp:1733
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition FlowProblem.hpp:1012
Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
Direct access to the transmissibility between two elements.
Definition FlowProblem.hpp:527
static int handlePositionalParameter(std::function< void(const std::string &, const std::string &)> addKey, std::set< std::string > &seenParams, std::string &errorMsg, int, const char **argv, int paramIdx, int)
Handles positional command line parameters.
Definition FlowProblem.hpp:197
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition FlowProblem.hpp:479
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition FlowProblem.hpp:509
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:826
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:659
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition FlowProblem.hpp:384
unsigned satnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:834
Scalar thermalHalfTransmissibilityOut(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:617
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition FlowProblem.hpp:574
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition FlowProblem.hpp:594
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:693
Scalar maxOilSaturation(unsigned globalDofIdx) const
Returns an element's historic maximum oil phase saturation that was observed during the simulation.
Definition FlowProblem.hpp:910
std::string name() const
The problem name.
Definition FlowProblem.hpp:865
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition FlowProblem.hpp:1104
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition FlowProblem.hpp:394
FlowProblem(Simulator &simulator)
Definition FlowProblem.hpp:216
const Vanguard::TransmissibilityType & eclTransmissibilities() const
Return a reference to the object that handles the "raw" transmissibilities.
Definition FlowProblem.hpp:640
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition FlowProblem.hpp:1033
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition FlowProblem.hpp:1638
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:769
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition FlowProblem.hpp:516
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition FlowProblem.hpp:295
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition FlowProblem.hpp:352
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition FlowProblem.hpp:684
LhsEval wellTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Return the well transmissibility multiplier due to rock changes.
Definition FlowProblem.hpp:1116
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition FlowProblem.hpp:781
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:850
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:872
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition FlowProblem.hpp:630
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:842
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:495
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition FlowProblem.hpp:859
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:604
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition FlowProblem.hpp:927
const SolidEnergyLawParams & solidEnergyLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Return the parameters for the energy storage law of the rock.
Definition FlowProblem.hpp:756
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
give the dispersivity for a face i.e.
Definition FlowProblem.hpp:554
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblem.hpp:181
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition FlowProblem.hpp:456
void source(RateVector &rate, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition FlowProblem.hpp:971
virtual void endTimeStep()
Called by the simulator after each time integration.
Definition FlowProblem.hpp:404
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition FlowProblem.hpp:938
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:734
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition FlowProblem.hpp:561
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition FlowProblem.hpp:536
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition FlowProblem.hpp:712
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition FlowProblem.hpp:261
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition FlowProblem.hpp:818
Scalar dofCenterDepth(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the depth of an degree of freedom [m].
Definition FlowProblem.hpp:672
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition FlowProblem.hpp:1060
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition FlowProblem.hpp:584
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition FlowProblem.hpp:280
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:703
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
give the transmissibility for a face i.e.
Definition FlowProblem.hpp:547
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition pffgridvector.hh:50
This class is intend to be a relperm diagnostics, to detect wrong input of relperm table and endpoint...
Definition RelpermDiagnostics.hpp:51
void diagnosis(const EclipseState &eclState, const LevelCartesianIndexMapper &levelCartesianIndexMapper)
This function is used to diagnosis relperm in eclipse data file.
Definition RelpermDiagnostics.cpp:830
A class which handles tracers as specified in by ECL.
Definition TracerModel.hpp:75
This file contains definitions related to directional mobilities.
The base class for the element-centered finite-volume discretization scheme.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition propertysystem.hh:224
Definition FlowProblem.hpp:1690
Definition FlowProblem.hpp:1479