31#ifndef OPM_FLOW_PROBLEM_BLACK_HPP
32#define OPM_FLOW_PROBLEM_BLACK_HPP
34#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
35#include <opm/material/fluidsystems/blackoilpvt/DryGasPvt.hpp>
36#include <opm/material/fluidsystems/blackoilpvt/WetGasPvt.hpp>
37#include <opm/material/fluidsystems/blackoilpvt/LiveOilPvt.hpp>
38#include <opm/material/fluidsystems/blackoilpvt/DeadOilPvt.hpp>
39#include <opm/material/fluidsystems/blackoilpvt/ConstantCompressibilityOilPvt.hpp>
40#include <opm/material/fluidsystems/blackoilpvt/ConstantCompressibilityWaterPvt.hpp>
44#include <opm/output/eclipse/EclipseIO.hpp>
46#include <opm/input/eclipse/Units/Units.hpp>
82template <
class TypeTag>
138 enum { enableDissolvedGas = Indices::compositionSwitchIdx >= 0 };
139 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
140 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
142 using SolventModule = BlackOilSolventModule<TypeTag>;
143 using PolymerModule = BlackOilPolymerModule<TypeTag>;
144 using FoamModule = BlackOilFoamModule<TypeTag>;
145 using BrineModule = BlackOilBrineModule<TypeTag>;
146 using ExtboModule = BlackOilExtboModule<TypeTag>;
147 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
148 using DispersionModule = BlackOilDispersionModule<TypeTag, enableDispersion>;
149 using DiffusionModule = BlackOilDiffusionModule<TypeTag, enableDiffusion>;
150 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
151 using ModuleParams =
typename BlackOilLocalResidualTPFA<TypeTag>::ModuleParams;
152 using HybridNewton = BlackOilHybridNewton<TypeTag>;
155 using EclWriterType = EclWriter<TypeTag, OutputBlackOilModule<TypeTag> >;
156 using IndexTraits =
typename FluidSystem::IndexTraitsType;
158 using DamarisWriterType = DamarisWriter<TypeTag>;
175 DamarisWriterType::registerParameters();
188 simulator.vanguard().schedule(),
189 simulator.vanguard().actionState(),
190 simulator.vanguard().summaryState(),
192 simulator.vanguard().grid().comm())
198 const auto& vanguard = simulator.vanguard();
201 brineParams.template initFromState<enableBrine,
202 enableSaltPrecipitation>(vanguard.eclState());
205 DiffusionModule::initFromState(vanguard.eclState());
206 DispersionModule::initFromState(vanguard.eclState());
209 extboParams.template initFromState<enableExtbo>(vanguard.eclState());
213 foamParams.template initFromState<enableFoam>(vanguard.eclState());
217 bioeffectsParams.template initFromState<enableBioeffects, enableMICP>(vanguard.eclState());
221 polymerParams.template initFromState<enablePolymer, enablePolymerMolarWeight>(vanguard.eclState());
225 solventParams.template initFromState<enableSolvent>(vanguard.eclState(), vanguard.schedule());
229 eclWriter_ = std::make_unique<EclWriterType>(simulator);
234 damarisWriter_ = std::make_unique<DamarisWriterType>(simulator);
235 enableDamarisOutput_ = Parameters::Get<Parameters::EnableDamarisOutput>();
246 auto& simulator = this->simulator();
248 const int episodeIdx = simulator.episodeIndex();
249 const auto& schedule = simulator.vanguard().schedule();
254 .evalUDQAssignments(episodeIdx, simulator.vanguard().udqState());
256 if (episodeIdx >= 0) {
257 const auto& oilVap = schedule[episodeIdx].oilvap();
258 if (oilVap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS) {
259 FluidSystem::setVapPars(oilVap.vap1(), oilVap.vap2());
262 FluidSystem::setVapPars(0.0, 0.0);
266 ConvectiveMixingModule::beginEpisode(simulator.vanguard().eclState(), schedule, episodeIdx,
267 this->moduleParams_.convectiveMixingModuleParam);
287 FlowProblemType::finishInit();
289 auto& simulator = this->simulator();
291 auto finishTransmissibilities = [updated =
false,
this]()
mutable
293 if (updated) {
return; }
295 this->
transmissibilities_.finishInit([&vg = this->simulator().vanguard()](
const unsigned int it) {
296 return vg.gridIdxToEquilGridIdx(it);
307 if (simulator.vanguard().grid().comm().size() > 1) {
308 if (simulator.vanguard().grid().comm().rank() == 0)
309 eclWriter_->setTransmissibilities(&simulator.vanguard().globalTransmissibility());
311 finishTransmissibilities();
312 eclWriter_->setTransmissibilities(&simulator.problem().eclTransmissibilities());
315 std::function<
unsigned int(
unsigned int)> equilGridToGrid = [&simulator](
unsigned int i) {
316 return simulator.vanguard().gridEquilIdxToGridIdx(i);
319 this->
eclWriter_->extractOutputTransAndNNC(equilGridToGrid);
321 simulator.vanguard().releaseGlobalTransmissibilities();
323 const auto& eclState = simulator.vanguard().eclState();
324 const auto& schedule = simulator.vanguard().schedule();
327 simulator.setStartTime(schedule.getStartTime());
328 simulator.setEndTime(schedule.simTime(schedule.size() - 1));
334 simulator.setEpisodeIndex(-1);
335 simulator.setEpisodeLength(0.0);
340 this->gravity_ = 0.0;
341 if (Parameters::Get<Parameters::EnableGravity>() &&
342 eclState.getInitConfig().hasGravity())
345 this->gravity_[dim - 1] = unit::gravity;
351 const auto& tuning = schedule[0].tuning();
358 bool isThermal = eclState.getSimulationConfig().isThermal();
359 bool isTemp = eclState.getSimulationConfig().isTemp();
360 bool conserveInnerEnergy = isTemp || (isThermal && Parameters::Get<Parameters::ConserveInnerEnergyThermal>());
361 FluidSystem::setEnergyEqualEnthalpy(conserveInnerEnergy);
363 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
364 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
369 [&simulator](
const unsigned idx)
371 std::array<int,dim> coords;
372 simulator.vanguard().cartesianCoordinate(idx, coords);
373 std::transform(coords.begin(), coords.end(), coords.begin(),
374 [](const auto c) { return c + 1; });
386 finishTransmissibilities();
388 const auto& initconfig = eclState.getInitConfig();
390 if (initconfig.restartRequested()) {
401 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>()) {
402 const auto& vanguard = this->simulator().vanguard();
403 const auto& gridView = vanguard.gridView();
404 const int numElements = gridView.size(0);
405 this->
polymer_.maxAdsorption.resize(numElements, 0.0);
414 this->
drift_.resize(this->model().numGridDof());
421 if (!initconfig.restartRequested() && !eclState.getIOConfig().initOnly()) {
422 simulator.startNextEpisode(schedule.seconds(1));
423 simulator.setEpisodeIndex(0);
424 simulator.setTimeStepIndex(0);
427 if (Parameters::Get<Parameters::CheckSatfuncConsistency>() &&
437 this->simulator().vanguard().grid().comm().barrier();
439 throw std::domain_error {
440 "Saturation function end-points do not "
441 "meet requisite consistency conditions"
450 eclState.runspec().tabdims().getNumPVTTables());
452 if (this->enableVtkOutput_() && eclState.getIOConfig().initOnly()) {
453 simulator.setTimeStepSize(0.0);
454 simulator.model().applyInitialSolution();
458 if (!eclState.getIOConfig().initOnly()) {
459 if (!this->
enableTuning_ && eclState.getSimulationConfig().anyTUNING()) {
460 OpmLog::info(
"\nThe deck has TUNING in the SCHEDULE section, but "
461 "it is ignored due\nto the flag --enable-tuning=false. "
462 "Set this flag to true to activate it.\n"
463 "Manually tuning the simulator with the TUNING keyword may "
464 "increase run time.\nIt is recommended using the simulator's "
465 "default tuning (--enable-tuning=false).");
475 FlowProblemType::endTimeStep();
476 this->endStepApplyAction();
483 this->eclWriter().mutableOutputModule().invalidateLocalData();
486 const auto& grid = this->simulator().vanguard().gridView().grid();
488 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
489 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
490 if (!isCpGrid || (grid.maxLevel() == 0)) {
491 this->eclWriter_->evalSummaryState(!this->episodeWillBeOver());
495 OPM_TIMEBLOCK(applyActions);
497 const int episodeIdx = this->episodeIndex();
498 auto& simulator = this->simulator();
502 this->simulator().vanguard().schedule().clearEvents(episodeIdx);
506 .applyActions(episodeIdx, simulator.time() + simulator.timeStepSize(),
507 [
this](
const bool global)
509 using TransUpdateQuantities = typename
510 Vanguard::TransmissibilityType::TransUpdateQuantities;
512 this->transmissibilities_
513 .update(global, TransUpdateQuantities::All,
514 [&vg = this->simulator().vanguard()]
515 (const unsigned int i)
517 return vg.gridIdxToEquilGridIdx(i);
528 OPM_TIMEBLOCK(endEpisode);
541 .evalUDQAssignments(this->episodeIndex(), this->simulator().vanguard().udqState());
543 FlowProblemType::endEpisode();
548 if (this->enableEclOutput_) {
549 this->eclWriter_->writeReports(timer);
560 FlowProblemType::writeOutput(verbose);
562 const auto isSubStep = !this->episodeWillBeOver();
564 auto localCellData = data::Solution {};
569 if (this->enableDamarisOutput_ && (this->damarisWriter_ !=
nullptr)) {
570 this->damarisWriter_->writeOutput(localCellData, isSubStep);
574 if (this->enableEclOutput_ && (this->eclWriter_ !=
nullptr)) {
575 this->eclWriter_->writeOutput(std::move(localCellData), isSubStep);
581 OPM_TIMEBLOCK(finalizeOutput);
593 FlowProblemType::initialSolutionApplied();
598 this->thresholdPressures_.finishInit();
601 const auto& grid = this->simulator().vanguard().gridView().grid();
603 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
604 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
606 if (!isCpGrid || (grid.maxLevel() == 0)) {
607 if (this->simulator().episodeIndex() == 0) {
608 eclWriter_->writeInitialFIPReport();
614 unsigned globalDofIdx,
615 unsigned timeIdx)
const override
617 this->aquiferModel_.addToSource(rate, globalDofIdx, timeIdx);
620 const auto& source = this->simulator().vanguard().schedule()[this->episodeIndex()].source();
621 std::array<int,3> ijk;
622 this->simulator().vanguard().cartesianCoordinate(globalDofIdx, ijk);
624 if (source.hasSource(ijk)) {
625 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
626 static std::array<SourceComponent, 3> sc_map = {SourceComponent::WATER, SourceComponent::OIL, SourceComponent::GAS};
627 static std::array<int, 3> phidx_map = {FluidSystem::waterPhaseIdx, FluidSystem::oilPhaseIdx, FluidSystem::gasPhaseIdx};
628 static std::array<int, 3> cidx_map = {waterCompIdx, oilCompIdx, gasCompIdx};
630 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
631 const auto phaseIdx = phidx_map[i];
632 const auto sourceComp = sc_map[i];
633 const auto compIdx = cidx_map[i];
634 if (!FluidSystem::phaseIsActive(phaseIdx)) {
637 Scalar mass_rate = source.rate(ijk, sourceComp) / this->model().dofTotalVolume(globalDofIdx);
638 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
639 mass_rate /= FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
641 rate[FluidSystem::canonicalToActiveCompIdx(compIdx)] += mass_rate;
644 if constexpr (enableSolvent) {
645 Scalar mass_rate = source.rate(ijk, SourceComponent::SOLVENT) / this->model().dofTotalVolume(globalDofIdx);
646 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
647 const auto& solventPvt = SolventModule::solventPvt();
648 mass_rate /= solventPvt.referenceDensity(pvtRegionIdx);
650 rate[Indices::contiSolventEqIdx] += mass_rate;
652 if constexpr (enablePolymer) {
653 rate[Indices::polymerConcentrationIdx] += source.rate(ijk, SourceComponent::POLYMER) / this->model().dofTotalVolume(globalDofIdx);
655 if constexpr (enableMICP) {
656 rate[Indices::microbialConcentrationIdx] += source.rate(ijk, SourceComponent::MICR) / this->model().dofTotalVolume(globalDofIdx);
657 rate[Indices::oxygenConcentrationIdx] += source.rate(ijk, SourceComponent::OXYG) / this->model().dofTotalVolume(globalDofIdx);
658 rate[Indices::ureaConcentrationIdx] += source.rate(ijk, SourceComponent::UREA) / (this->model().dofTotalVolume(globalDofIdx));
660 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
661 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
662 const auto phaseIdx = phidx_map[i];
663 if (!FluidSystem::phaseIsActive(phaseIdx)) {
666 const auto sourceComp = sc_map[i];
667 const auto source_hrate = source.hrate(ijk, sourceComp);
669 rate[Indices::contiEnergyEqIdx] += source_hrate.value() / this->model().dofTotalVolume(globalDofIdx);
671 const auto& intQuants = this->simulator().model().intensiveQuantities(globalDofIdx, 0);
672 auto fs = intQuants.fluidState();
674 const auto source_temp = source.temperature(ijk, sourceComp);
676 Scalar temperature = source_temp.value();
677 fs.setTemperature(temperature);
679 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, pvtRegionIdx);
680 Scalar mass_rate = source.rate(ijk, sourceComp)/ this->model().dofTotalVolume(globalDofIdx);
681 Scalar energy_rate = getValue(h)*mass_rate;
682 rate[Indices::contiEnergyEqIdx] += energy_rate;
690 if (this->enableDriftCompensation_) {
691 const auto& simulator = this->simulator();
692 const auto& model = this->model();
697 Scalar maxCompensation = model.newtonMethod().tolerance()/10;
698 Scalar poro = this->porosity(globalDofIdx, timeIdx);
699 Scalar dt = simulator.timeStepSize();
700 EqVector dofDriftRate = this->drift_[globalDofIdx];
701 dofDriftRate /= dt*model.dofTotalVolume(globalDofIdx);
704 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
705 Scalar cnv = std::abs(dofDriftRate[eqIdx])*dt*model.eqWeight(globalDofIdx, eqIdx)/poro;
706 if (cnv > maxCompensation) {
707 dofDriftRate[eqIdx] *= maxCompensation/cnv;
711 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx)
712 rate[eqIdx] -= dofDriftRate[eqIdx];
719 template <
class LhsEval,
class Callback>
722 OPM_TIMEBLOCK_LOCAL(permFactTransMultiplier, Subsystem::PvtProps);
723 if constexpr (enableSaltPrecipitation) {
724 const auto& fs = intQuants.fluidState();
725 unsigned tableIdx = this->simulator().problem().satnumRegionIndex(elementIdx);
726 LhsEval porosityFactor = obtain(1. - fs.saltSaturation());
727 porosityFactor = min(porosityFactor, 1.0);
728 const auto& permfactTable = BrineModule::permfactTable(tableIdx);
729 return permfactTable.eval(porosityFactor,
true);
731 else if constexpr (enableBioeffects) {
732 return obtain(intQuants.permFactor());
741 {
return initialFluidStates_[globalDofIdx]; }
744 {
return initialFluidStates_; }
747 {
return initialFluidStates_; }
750 {
return eclWriter_->eclIO(); }
753 {
return eclWriter_->setSubStepReport(report); }
756 {
return eclWriter_->setSimulationReport(report); }
760 OPM_TIMEBLOCK_LOCAL(boundaryFluidState, Subsystem::Assembly);
761 const auto& bcprop = this->simulator().vanguard().schedule()[this->episodeIndex()].bcprop;
762 if (bcprop.size() > 0) {
763 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
767 if (this->bcindex_(dir)[globalDofIdx] == 0)
768 return initialFluidStates_[globalDofIdx];
770 const auto& bc = bcprop[this->bcindex_(dir)[globalDofIdx]];
771 if (bc.bctype == BCType::DIRICHLET )
773 InitialFluidState fluidState;
774 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
775 fluidState.setPvtRegionIndex(pvtRegionIdx);
777 switch (bc.component) {
778 case BCComponent::OIL:
779 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx))
780 throw std::logic_error(
"oil is not active and you're trying to add oil BC");
782 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0);
784 case BCComponent::GAS:
785 if (!FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
786 throw std::logic_error(
"gas is not active and you're trying to add gas BC");
788 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 1.0);
790 case BCComponent::WATER:
791 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
792 throw std::logic_error(
"water is not active and you're trying to add water BC");
794 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0);
796 case BCComponent::SOLVENT:
797 case BCComponent::POLYMER:
798 case BCComponent::MICR:
799 case BCComponent::OXYG:
800 case BCComponent::UREA:
802 throw std::logic_error(
"you need to specify a valid component (OIL, WATER or GAS) when DIRICHLET type is set in BC");
804 fluidState.setTotalSaturation(1.0);
805 double pressure = initialFluidStates_[globalDofIdx].pressure(this->refPressurePhaseIdx_());
806 const auto pressure_input = bc.pressure;
807 if (pressure_input) {
808 pressure = *pressure_input;
811 std::array<Scalar, numPhases> pc = {0};
812 const auto& matParams = this->materialLawParams(globalDofIdx);
813 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
814 Valgrind::CheckDefined(pressure);
815 Valgrind::CheckDefined(pc);
816 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
817 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
818 if (Indices::oilEnabled)
819 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
820 else if (Indices::gasEnabled)
821 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
822 else if (Indices::waterEnabled)
824 fluidState.setPressure(phaseIdx, pressure);
826 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
827 double temperature = initialFluidStates_[globalDofIdx].temperature(0);
828 const auto temperature_input = bc.temperature;
829 if(temperature_input)
830 temperature = *temperature_input;
831 fluidState.setTemperature(temperature);
834 if constexpr (enableDissolvedGas) {
835 if (FluidSystem::enableDissolvedGas()) {
836 fluidState.setRs(0.0);
837 fluidState.setRv(0.0);
840 if constexpr (enableDisgasInWater) {
841 if (FluidSystem::enableDissolvedGasInWater()) {
842 fluidState.setRsw(0.0);
845 if constexpr (enableVapwat) {
846 if (FluidSystem::enableVaporizedWater()) {
847 fluidState.setRvw(0.0);
851 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
852 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
854 const auto& b = FluidSystem::inverseFormationVolumeFactor(fluidState, phaseIdx, pvtRegionIdx);
855 fluidState.setInvB(phaseIdx, b);
857 const auto& rho = FluidSystem::density(fluidState, phaseIdx, pvtRegionIdx);
858 fluidState.setDensity(phaseIdx, rho);
859 if constexpr (energyModuleType == EnergyModules::SequentialImplicitThermal || energyModuleType == EnergyModules::FullyImplicitThermal) {
860 const auto& h = FluidSystem::enthalpy(fluidState, phaseIdx, pvtRegionIdx);
861 fluidState.setEnthalpy(phaseIdx, h);
864 fluidState.checkDefined();
868 return initialFluidStates_[globalDofIdx];
873 {
return *eclWriter_; }
876 {
return *eclWriter_; }
884 return this->mixControls_.maxGasDissolutionFactor(timeIdx, globalDofIdx,
885 this->episodeIndex(),
886 this->pvtRegionIndex(globalDofIdx));
895 return this->mixControls_.maxOilVaporizationFactor(timeIdx, globalDofIdx,
896 this->episodeIndex(),
897 this->pvtRegionIndex(globalDofIdx));
910 int episodeIdx = this->episodeIndex();
911 return !this->mixControls_.drsdtActive(episodeIdx) &&
912 !this->mixControls_.drvdtActive(episodeIdx) &&
913 this->rockCompPoroMultWc_.empty() &&
914 this->rockCompPoroMult_.empty();
923 template <
class Context>
926 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
928 values.setPvtRegionIndex(pvtRegionIndex(context, spaceIdx, timeIdx));
929 values.assignNaive(initialFluidStates_[globalDofIdx]);
931 SolventModule::assignPrimaryVars(values,
932 enableSolvent ? this->solventSaturation_[globalDofIdx] : 0.0,
933 enableSolvent ? this->solventRsw_[globalDofIdx] : 0.0);
935 if constexpr (enablePolymer)
936 values[Indices::polymerConcentrationIdx] = this->polymer_.concentration[globalDofIdx];
938 if constexpr (enablePolymerMolarWeight)
939 values[Indices::polymerMoleWeightIdx]= this->polymer_.moleWeight[globalDofIdx];
941 if constexpr (enableBrine) {
942 if (enableSaltPrecipitation && values.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
943 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltSaturation();
946 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltConcentration();
950 if constexpr (enableBioeffects) {
951 values[Indices::microbialConcentrationIdx] = this->bioeffects_.microbialConcentration[globalDofIdx];
952 values[Indices::biofilmVolumeFractionIdx]= this->bioeffects_.biofilmVolumeFraction[globalDofIdx];
953 if constexpr (enableMICP) {
954 values[Indices::oxygenConcentrationIdx]= this->bioeffects_.oxygenConcentration[globalDofIdx];
955 values[Indices::ureaConcentrationIdx]= this->bioeffects_.ureaConcentration[globalDofIdx];
956 values[Indices::calciteVolumeFractionIdx]= this->bioeffects_.calciteVolumeFraction[globalDofIdx];
960 values.checkDefined();
966 return this->mixControls_.drsdtcon(elemIdx, episodeIdx,
967 this->pvtRegionIndex(elemIdx));
972 return this->mixControls_.drsdtConvective(episodeIdx, this->pvtRegionIndex(elemIdx));
980 template <
class Context>
982 const Context& context,
984 unsigned timeIdx)
const
986 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
987 if (!context.intersection(spaceIdx).boundary())
990 if constexpr (energyModuleType != EnergyModules::FullyImplicitThermal || !enableThermalFluxBoundaries)
998 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
999 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1000 values.setThermalFlow(context, spaceIdx, timeIdx, this->initialFluidStates_[globalDofIdx] );
1003 if (this->nonTrivialBoundaryConditions()) {
1004 unsigned indexInInside = context.intersection(spaceIdx).indexInInside();
1005 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1006 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1007 unsigned pvtRegionIdx = pvtRegionIndex(context, spaceIdx, timeIdx);
1008 const auto [type, massrate] = this->boundaryCondition(globalDofIdx, indexInInside);
1009 if (type == BCType::THERMAL)
1010 values.setThermalFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1011 else if (type == BCType::FREE || type == BCType::DIRICHLET)
1012 values.setFreeFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1013 else if (type == BCType::RATE)
1014 values.setMassRate(massrate, pvtRegionIdx);
1024 auto& simulator = this->simulator();
1025 const auto& eclState = simulator.vanguard().eclState();
1027 std::size_t numElems = this->model().numGridDof();
1028 this->initialFluidStates_.resize(numElems);
1029 if constexpr (enableSolvent) {
1030 this->solventSaturation_.resize(numElems, 0.0);
1031 this->solventRsw_.resize(numElems, 0.0);
1034 if constexpr (enablePolymer)
1035 this->polymer_.concentration.resize(numElems, 0.0);
1037 if constexpr (enablePolymerMolarWeight) {
1038 const std::string msg {
"Support of the RESTART for polymer molecular weight "
1039 "is not implemented yet. The polymer weight value will be "
1040 "zero when RESTART begins"};
1041 OpmLog::warning(
"NO_POLYMW_RESTART", msg);
1042 this->polymer_.moleWeight.resize(numElems, 0.0);
1045 if constexpr (enableBioeffects) {
1046 this->bioeffects_.resize(numElems);
1050 this->mixControls_.init(numElems, restart_step, eclState.runspec().tabdims().getNumPVTTables());
1052 if constexpr (enableBioeffects) {
1053 this->bioeffects_ = this->eclWriter_->outputModule().getBioeffects().getSolution();
1056 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1057 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1058 elemFluidState.setPvtRegionIndex(pvtRegionIndex(elemIdx));
1059 this->eclWriter_->outputModule().initHysteresisParams(simulator, elemIdx);
1060 this->eclWriter_->outputModule().assignToFluidState(elemFluidState, elemIdx);
1069 auto ssol = enableSolvent
1070 ? this->eclWriter_->outputModule().getSolventSaturation(elemIdx)
1073 this->processRestartSaturations_(elemFluidState, ssol);
1075 if constexpr (enableSolvent) {
1076 this->solventSaturation_[elemIdx] = ssol;
1077 this->solventRsw_[elemIdx] = this->eclWriter_->outputModule().getSolventRsw(elemIdx);
1082 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1083 bool needTemperature = (eclState.runspec().co2Storage() || eclState.runspec().h2Storage());
1084 if (needTemperature) {
1085 const auto& fp = simulator.vanguard().eclState().fieldProps();
1086 elemFluidState.setTemperature(fp.get_double(
"TEMPI")[elemIdx]);
1090 this->mixControls_.updateLastValues(elemIdx, elemFluidState.Rs(), elemFluidState.Rv());
1092 if constexpr (enablePolymer)
1093 this->polymer_.concentration[elemIdx] = this->eclWriter_->outputModule().getPolymerConcentration(elemIdx);
1097 const int episodeIdx = this->episodeIndex();
1098 this->mixControls_.updateMaxValues(episodeIdx, simulator.timeStepSize());
1103 auto& sol = this->model().solution(0);
1104 const auto& gridView = this->gridView();
1106 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1107 elemCtx.updatePrimaryStencil(elem);
1108 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1109 this->initial(sol[elemIdx], elemCtx, 0, 0);
1117 this->model().syncOverlap();
1120 this->updateReferencePorosity_();
1121 this->mixControls_.init(this->model().numGridDof(),
1122 this->episodeIndex(),
1123 eclState.runspec().tabdims().getNumPVTTables());
1131 {
return thresholdPressures_.thresholdPressure(elem1Idx, elem2Idx); }
1134 {
return thresholdPressures_; }
1137 {
return thresholdPressures_; }
1141 return moduleParams_;
1144 template<
class Serializer>
1148 serializer(mixControls_);
1149 serializer(*eclWriter_);
1155 this->updateExplicitQuantities_(first_step_after_restart);
1157 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>())
1158 updateMaxPolymerAdsorption_();
1160 mixControls_.updateExplicitQuantities(episodeIdx, timeStepSize);
1166 this->updateProperty_(
"FlowProblemBlackoil::updateMaxPolymerAdsorption_() failed:",
1169 this->updateMaxPolymerAdsorption_(compressedDofIdx,iq);
1175 const Scalar pa = scalarValue(iq.polymerAdsorption());
1176 auto& mpa = this->polymer_.maxAdsorption;
1177 if (mpa[compressedDofIdx] < pa) {
1178 mpa[compressedDofIdx] = pa;
1187 std::vector<Scalar> sumInvB(numPhases, 0.0);
1188 const auto& gridView = this->gridView();
1190 for(
const auto& elem: elements(gridView, Dune::Partitions::interior)) {
1191 elemCtx.updatePrimaryStencil(elem);
1192 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1193 const auto& dofFluidState = this->initialFluidStates_[elemIdx];
1194 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1195 if (!FluidSystem::phaseIsActive(phaseIdx))
1198 sumInvB[phaseIdx] += dofFluidState.invB(phaseIdx);
1202 std::size_t numDof = this->model().numGridDof();
1203 const auto& comm = this->simulator().vanguard().grid().comm();
1204 comm.sum(sumInvB.data(),sumInvB.size());
1205 Scalar numTotalDof = comm.sum(numDof);
1207 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1208 if (!FluidSystem::phaseIsActive(phaseIdx))
1211 Scalar avgB = numTotalDof / sumInvB[phaseIdx];
1212 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
1213 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
1214 this->model().setEqWeight(activeSolventCompIdx, avgB);
1221 OPM_TIMEBLOCK(updateCompositionChangeLimits);
1224 int episodeIdx = this->episodeIndex();
1225 std::array<bool,3> active{this->mixControls_.drsdtConvective(episodeIdx),
1226 this->mixControls_.drsdtActive(episodeIdx),
1227 this->mixControls_.drvdtActive(episodeIdx)};
1228 if (!active[0] && !active[1] && !active[2]) {
1232 this->updateProperty_(
"FlowProblemBlackoil::updateCompositionChangeLimits_()) failed:",
1233 [
this,episodeIdx,active](
unsigned compressedDofIdx,
1236 const DimMatrix& perm = this->intrinsicPermeability(compressedDofIdx);
1237 const Scalar distZ = active[0] ? this->simulator().vanguard().cellThickness(compressedDofIdx) : 0.0;
1238 const int pvtRegionIdx = this->pvtRegionIndex(compressedDofIdx);
1239 this->mixControls_.update(compressedDofIdx,
1242 this->gravity_[
dim - 1],
1255 if(this->simulator().vanguard().grid().maxLevel() > 0) {
1256 throw std::invalid_argument(
"Refined grids are not yet supported for restart ");
1260 auto& simulator = this->simulator();
1261 const auto& schedule = simulator.vanguard().schedule();
1262 const auto& eclState = simulator.vanguard().eclState();
1263 const auto& initconfig = eclState.getInitConfig();
1264 const int restart_step = initconfig.getRestartStep();
1266 simulator.setTime(schedule.seconds(restart_step));
1268 simulator.startNextEpisode(simulator.startTime() + simulator.time(),
1269 schedule.stepLength(restart_step));
1270 simulator.setEpisodeIndex(restart_step);
1272 this->eclWriter_->beginRestart();
1274 Scalar dt = std::min(this->eclWriter_->restartTimeStepSize(), simulator.episodeLength());
1275 simulator.setTimeStepSize(dt);
1277 this->readSolutionFromOutputModule(restart_step,
false);
1279 this->eclWriter_->endRestart();
1284 const auto& simulator = this->simulator();
1289 std::size_t numElems = this->model().numGridDof();
1290 this->initialFluidStates_.resize(numElems);
1291 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1292 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1299 const auto& simulator = this->simulator();
1300 const auto& vanguard = simulator.vanguard();
1301 const auto& eclState = vanguard.eclState();
1302 const auto& fp = eclState.fieldProps();
1303 bool has_swat = fp.has_double(
"SWAT");
1304 bool has_sgas = fp.has_double(
"SGAS");
1305 bool has_rs = fp.has_double(
"RS");
1306 bool has_rsw = fp.has_double(
"RSW");
1307 bool has_rv = fp.has_double(
"RV");
1308 bool has_rvw = fp.has_double(
"RVW");
1309 bool has_pressure = fp.has_double(
"PRESSURE");
1310 bool has_salt = fp.has_double(
"SALT");
1311 bool has_saltp = fp.has_double(
"SALTP");
1314 if (Indices::numPhases > 1) {
1315 if (FluidSystem::phaseIsActive(waterPhaseIdx) && !has_swat)
1316 throw std::runtime_error(
"The ECL input file requires the presence of the SWAT keyword if "
1317 "the water phase is active");
1318 if (FluidSystem::phaseIsActive(gasPhaseIdx) && !has_sgas && FluidSystem::phaseIsActive(oilPhaseIdx))
1319 throw std::runtime_error(
"The ECL input file requires the presence of the SGAS keyword if "
1320 "the gas phase is active");
1323 throw std::runtime_error(
"The ECL input file requires the presence of the PRESSURE "
1324 "keyword if the model is initialized explicitly");
1325 if (FluidSystem::enableDissolvedGas() && !has_rs)
1326 throw std::runtime_error(
"The ECL input file requires the RS keyword to be present if"
1327 " dissolved gas is enabled and the model is initialized explicitly");
1328 if (FluidSystem::enableDissolvedGasInWater() && !has_rsw)
1329 OpmLog::warning(
"The model is initialized explicitly and the RSW keyword is not present in the"
1330 " ECL input file. The RSW values are set equal to 0");
1331 if (FluidSystem::enableVaporizedOil() && !has_rv)
1332 throw std::runtime_error(
"The ECL input file requires the RV keyword to be present if"
1333 " vaporized oil is enabled and the model is initialized explicitly");
1334 if (FluidSystem::enableVaporizedWater() && !has_rvw)
1335 throw std::runtime_error(
"The ECL input file requires the RVW keyword to be present if"
1336 " vaporized water is enabled and the model is initialized explicitly");
1337 if (enableBrine && !has_salt)
1338 throw std::runtime_error(
"The ECL input file requires the SALT keyword to be present if"
1339 " brine is enabled and the model is initialized explicitly");
1340 if (enableSaltPrecipitation && !has_saltp)
1341 throw std::runtime_error(
"The ECL input file requires the SALTP keyword to be present if"
1342 " salt precipitation is enabled and the model is initialized explicitly");
1344 std::size_t numDof = this->model().numGridDof();
1346 initialFluidStates_.resize(numDof);
1348 std::vector<double> waterSaturationData;
1349 std::vector<double> gasSaturationData;
1350 std::vector<double> pressureData;
1351 std::vector<double> rsData;
1352 std::vector<double> rswData;
1353 std::vector<double> rvData;
1354 std::vector<double> rvwData;
1355 std::vector<double> tempiData;
1356 std::vector<double> saltData;
1357 std::vector<double> saltpData;
1359 if (FluidSystem::phaseIsActive(waterPhaseIdx) && Indices::numPhases > 1)
1360 waterSaturationData = fp.get_double(
"SWAT");
1362 waterSaturationData.resize(numDof);
1364 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx))
1365 gasSaturationData = fp.get_double(
"SGAS");
1367 gasSaturationData.resize(numDof);
1369 pressureData = fp.get_double(
"PRESSURE");
1370 if (FluidSystem::enableDissolvedGas())
1371 rsData = fp.get_double(
"RS");
1373 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1374 rswData = fp.get_double(
"RSW");
1376 if (FluidSystem::enableVaporizedOil())
1377 rvData = fp.get_double(
"RV");
1379 if (FluidSystem::enableVaporizedWater())
1380 rvwData = fp.get_double(
"RVW");
1383 tempiData = fp.get_double(
"TEMPI");
1386 if constexpr (enableBrine)
1387 saltData = fp.get_double(
"SALT");
1390 if constexpr (enableSaltPrecipitation)
1391 saltpData = fp.get_double(
"SALTP");
1394 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1395 auto& dofFluidState = initialFluidStates_[dofIdx];
1397 dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
1402 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1403 Scalar temperatureLoc = tempiData[dofIdx];
1404 if (!std::isfinite(temperatureLoc) || temperatureLoc <= 0)
1405 temperatureLoc = FluidSystem::surfaceTemperature;
1406 dofFluidState.setTemperature(temperatureLoc);
1412 if constexpr (enableBrine)
1413 dofFluidState.setSaltConcentration(saltData[dofIdx]);
1418 if constexpr (enableSaltPrecipitation)
1419 dofFluidState.setSaltSaturation(saltpData[dofIdx]);
1424 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
1425 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
1426 waterSaturationData[dofIdx]);
1428 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
1429 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
1430 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1432 - waterSaturationData[dofIdx]);
1435 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1436 gasSaturationData[dofIdx]);
1438 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1439 const Scalar soil = 1.0 - waterSaturationData[dofIdx] - gasSaturationData[dofIdx];
1440 if (soil < smallSaturationTolerance_) {
1441 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
1444 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, soil);
1451 Scalar pressure = pressureData[dofIdx];
1455 std::array<Scalar, numPhases> pc = {0};
1456 const auto& matParams = this->materialLawParams(dofIdx);
1457 MaterialLaw::capillaryPressures(pc, matParams, dofFluidState);
1458 Valgrind::CheckDefined(pressure);
1459 Valgrind::CheckDefined(pc);
1460 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1461 if (!FluidSystem::phaseIsActive(phaseIdx))
1464 if (Indices::oilEnabled)
1465 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
1466 else if (Indices::gasEnabled)
1467 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
1468 else if (Indices::waterEnabled)
1470 dofFluidState.setPressure(phaseIdx, pressure);
1473 if constexpr (enableDissolvedGas) {
1474 if (FluidSystem::enableDissolvedGas())
1475 dofFluidState.setRs(rsData[dofIdx]);
1476 else if (Indices::gasEnabled && Indices::oilEnabled)
1477 dofFluidState.setRs(0.0);
1478 if (FluidSystem::enableVaporizedOil())
1479 dofFluidState.setRv(rvData[dofIdx]);
1480 else if (Indices::gasEnabled && Indices::oilEnabled)
1481 dofFluidState.setRv(0.0);
1484 if constexpr (enableDisgasInWater) {
1485 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1486 dofFluidState.setRsw(rswData[dofIdx]);
1489 if constexpr (enableVapwat) {
1490 if (FluidSystem::enableVaporizedWater())
1491 dofFluidState.setRvw(rvwData[dofIdx]);
1497 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1498 if (!FluidSystem::phaseIsActive(phaseIdx))
1501 const auto& b = FluidSystem::inverseFormationVolumeFactor(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1502 dofFluidState.setInvB(phaseIdx, b);
1504 const auto& rho = FluidSystem::density(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1505 dofFluidState.setDensity(phaseIdx, rho);
1516 Scalar sumSaturation = 0.0;
1517 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1518 if (FluidSystem::phaseIsActive(phaseIdx)) {
1519 if (elemFluidState.saturation(phaseIdx) < smallSaturationTolerance_)
1520 elemFluidState.setSaturation(phaseIdx, 0.0);
1522 sumSaturation += elemFluidState.saturation(phaseIdx);
1526 if constexpr (enableSolvent) {
1527 if (solventSaturation < smallSaturationTolerance_)
1528 solventSaturation = 0.0;
1530 sumSaturation += solventSaturation;
1533 assert(sumSaturation > 0.0);
1535 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1536 if (FluidSystem::phaseIsActive(phaseIdx)) {
1537 const Scalar saturation = elemFluidState.saturation(phaseIdx) / sumSaturation;
1538 elemFluidState.setSaturation(phaseIdx, saturation);
1541 if constexpr (enableSolvent) {
1542 solventSaturation = solventSaturation / sumSaturation;
1548 FlowProblemType::readInitialCondition_();
1550 if constexpr (enableSolvent || enablePolymer || enablePolymerMolarWeight || enableBioeffects)
1551 this->readBlackoilExtentionsInitialConditions_(this->model().numGridDof(),
1554 enablePolymerMolarWeight,
1562 if constexpr (!enableSolvent)
1563 throw std::logic_error(
"solvent is disabled and you're trying to add solvent to BC");
1565 rate[Indices::solventSaturationIdx] = bc.rate;
1570 if constexpr (!enablePolymer)
1571 throw std::logic_error(
"polymer is disabled and you're trying to add polymer to BC");
1573 rate[Indices::polymerConcentrationIdx] = bc.rate;
1578 if constexpr (!enableMICP)
1579 throw std::logic_error(
"MICP is disabled and you're trying to add microbes to BC");
1581 rate[Indices::microbialConcentrationIdx] = bc.rate;
1586 if constexpr (!enableMICP)
1587 throw std::logic_error(
"MICP is disabled and you're trying to add oxygen to BC");
1589 rate[Indices::oxygenConcentrationIdx] = bc.rate;
1594 if constexpr (!enableMICP)
1595 throw std::logic_error(
"MICP is disabled and you're trying to add urea to BC");
1597 rate[Indices::ureaConcentrationIdx] = bc.rate;
1599 rate[Indices::ureaConcentrationIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
1604 OPM_TIMEBLOCK(updateExplicitQuantities);
1605 const bool invalidateFromMaxWaterSat = this->updateMaxWaterSaturation_();
1606 const bool invalidateFromMinPressure = this->updateMinPressure_();
1609 const bool invalidateFromHyst = this->updateHysteresis_();
1610 const bool invalidateFromMaxOilSat = this->updateMaxOilSaturation_();
1613 const bool invalidateDRDT = !first_step_after_restart && this->updateCompositionChangeLimits_();
1616 const bool invalidateIntensiveQuantities
1617 = invalidateFromMaxWaterSat || invalidateFromMinPressure || invalidateFromHyst || invalidateFromMaxOilSat || invalidateDRDT;
1618 if (invalidateIntensiveQuantities) {
1619 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
1620 this->model().invalidateAndUpdateIntensiveQuantities(0);
1623 this->updateRockCompTransMultVal_();
1628 if (
const auto nph = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
1629 + FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
1630 + FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1639 const auto numSamplePoints =
static_cast<std::size_t
>
1640 (Parameters::Get<Parameters::NumSatfuncConsistencySamplePoints>());
1642 auto sfuncConsistencyChecks =
1644 numSamplePoints, this->simulator().vanguard().eclState(),
1645 [&cmap = this->simulator().vanguard().cartesianIndexMapper()](
const int elemIdx)
1646 {
return cmap.cartesianIndex(elemIdx); }
1649 const auto ioRank = 0;
1650 const auto isIoRank = this->simulator().vanguard()
1651 .grid().comm().rank() == ioRank;
1657 .
run(this->simulator().vanguard().grid().levelGridView(0),
1658 [&vg = this->simulator().vanguard(),
1659 &emap = this->simulator().model().elementMapper()]
1661 {
return vg.gridIdxToEquilGridIdx(emap.index(elem)); });
1663 using ViolationLevel =
typename Satfunc::PhaseChecks::
1664 SatfuncConsistencyCheckManager<Scalar>::ViolationLevel;
1666 auto reportFailures = [&sfuncConsistencyChecks]
1667 (
const ViolationLevel level)
1669 sfuncConsistencyChecks.reportFailures
1670 (level, [](std::string_view record)
1671 { OpmLog::info(std::string { record }); });
1674 if (sfuncConsistencyChecks.anyFailedStandardChecks()) {
1676 OpmLog::warning(
"Saturation Function "
1677 "End-point Consistency Problems");
1679 reportFailures(ViolationLevel::Standard);
1683 if (sfuncConsistencyChecks.anyFailedCriticalChecks()) {
1685 OpmLog::error(
"Saturation Function "
1686 "End-point Consistency Failures");
1688 reportFailures(ViolationLevel::Critical);
1708 const Scalar smallSaturationTolerance_ = 1.e-6;
1710 bool enableDamarisOutput_ = false ;
1711 std::unique_ptr<DamarisWriterType> damarisWriter_;
1732 bool episodeWillBeOver()
const override
1734 const auto currTime = this->simulator().time()
1735 + this->simulator().timeStepSize();
1737 const auto nextReportStep =
1738 this->simulator().vanguard().schedule()
1739 .seconds(this->simulator().episodeIndex() + 1);
1741 const auto isSubStep = (nextReportStep - currTime)
1742 > (2 * std::numeric_limits<float>::epsilon()) * nextReportStep;
Class handling Action support in simulator.
Definition: ActionHandler.hpp:52
static void setParams(BlackOilBioeffectsParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilbioeffectsmodules.hh:131
static void setParams(BlackOilBrineParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilbrinemodules.hh:88
static void setParams(BlackOilExtboParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilextbomodules.hh:87
static void setParams(BlackOilFoamParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilfoammodules.hh:88
Hybrid Newton solver extension for the black-oil model.
Definition: HybridNewton.hpp:59
void tryApplyHybridNewton()
Attempt to apply the Hybrid Newton correction at the current timestep.
Definition: HybridNewton.hpp:99
static void setParams(BlackOilPolymerParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilpolymermodules.hh:95
static void setParams(BlackOilSolventParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilsolventmodules.hh:99
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
static void registerParameters()
Definition: EclWriter.hpp:138
Computes the initial condition based on the EQUIL keyword from ECL.
Definition: EquilInitializer.hpp:59
BlackOilFluidState< Scalar, FluidSystem, energyModuleType==EnergyModules::ConstantTemperature,(energyModuleType==EnergyModules::FullyImplicitThermal||energyModuleType==EnergyModules::SequentialImplicitThermal), enableDissolution, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, Indices::numPhases > ScalarFluidState
Definition: EquilInitializer.hpp:98
const ScalarFluidState & initialFluidState(unsigned elemIdx) const
Return the initial thermodynamic state which should be used as the initial condition.
Definition: EquilInitializer.hpp:198
PolymerSolutionContainer< Scalar > polymer_
Definition: FlowGenericProblem.hpp:335
Scalar initialTimeStepSize_
Definition: FlowGenericProblem.hpp:346
bool enableDriftCompensation_
Definition: FlowGenericProblem.hpp:352
bool enableTuning_
Definition: FlowGenericProblem.hpp:345
void readRockParameters_(const std::vector< Scalar > &cellCenterDepths, std::function< std::array< int, 3 >(const unsigned)> ijkIndex)
Definition: FlowGenericProblem_impl.hpp:136
std::vector< Scalar > maxOilSaturation_
Definition: FlowGenericProblem.hpp:336
Scalar maxTimeStepAfterWellEvent_
Definition: FlowGenericProblem.hpp:347
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblemBlackoil.hpp:84
HybridNewton hybridNewton_
Definition: FlowProblemBlackoil.hpp:1719
void updateExplicitQuantities_(int episodeIdx, int timeStepSize, const bool first_step_after_restart) override
Definition: FlowProblemBlackoil.hpp:1153
bool updateMaxPolymerAdsorption_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblemBlackoil.hpp:1173
void handlePolymerBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1568
void writeOutput(const bool verbose) override
Write the requested quantities of the current solution into the output files.
Definition: FlowProblemBlackoil.hpp:558
void readInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1546
void handleOxygBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1584
void readEquilInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1282
void handleSolventBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1560
Scalar maxGasDissolutionFactor(unsigned timeIdx, unsigned globalDofIdx) const
Returns the maximum value of the gas dissolution factor at the current time for a given degree of fre...
Definition: FlowProblemBlackoil.hpp:882
const std::vector< InitialFluidState > & initialFluidStates() const
Definition: FlowProblemBlackoil.hpp:746
void processRestartSaturations_(InitialFluidState &elemFluidState, Scalar &solventSaturation)
Definition: FlowProblemBlackoil.hpp:1512
std::vector< InitialFluidState > & initialFluidStates()
Definition: FlowProblemBlackoil.hpp:743
FlowProblemBlackoil(Simulator &simulator)
Definition: FlowProblemBlackoil.hpp:183
bool enableEclOutput_
Definition: FlowProblemBlackoil.hpp:1705
Scalar drsdtcon(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:964
void endStepApplyAction()
Definition: FlowProblemBlackoil.hpp:479
bool drsdtconIsActive(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:970
Scalar maxOilVaporizationFactor(unsigned timeIdx, unsigned globalDofIdx) const
Returns the maximum value of the oil vaporization factor at the current time for a given degree of fr...
Definition: FlowProblemBlackoil.hpp:893
std::vector< InitialFluidState > initialFluidStates_
Definition: FlowProblemBlackoil.hpp:1703
void endTimeStep() override
Called by the simulator after each time integration.
Definition: FlowProblemBlackoil.hpp:473
void updateMaxPolymerAdsorption_()
Definition: FlowProblemBlackoil.hpp:1163
const InitialFluidState & initialFluidState(unsigned globalDofIdx) const
Definition: FlowProblemBlackoil.hpp:740
void endEpisode() override
Called by the simulator after the end of an episode.
Definition: FlowProblemBlackoil.hpp:526
void setSubStepReport(const SimulatorReportSingle &report)
Definition: FlowProblemBlackoil.hpp:752
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: FlowProblemBlackoil.hpp:924
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: FlowProblemBlackoil.hpp:282
void handleMicrBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1576
const FlowThresholdPressure< TypeTag > & thresholdPressure() const
Definition: FlowProblemBlackoil.hpp:1133
void finalizeOutput()
Definition: FlowProblemBlackoil.hpp:579
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: FlowProblemBlackoil.hpp:981
InitialFluidState boundaryFluidState(unsigned globalDofIdx, const int directionId) const
Definition: FlowProblemBlackoil.hpp:758
std::unique_ptr< EclWriterType > eclWriter_
Definition: FlowProblemBlackoil.hpp:1706
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblemBlackoil.hpp:591
const ModuleParams & moduleParams() const
Definition: FlowProblemBlackoil.hpp:1139
void readEclRestartSolution_()
Definition: FlowProblemBlackoil.hpp:1252
FlowThresholdPressure< TypeTag > & thresholdPressure()
Definition: FlowProblemBlackoil.hpp:1136
FlowThresholdPressure< TypeTag > thresholdPressures_
Definition: FlowProblemBlackoil.hpp:1701
void readExplicitInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1297
void beginEpisode() override
Called by the simulator before an episode begins.
Definition: FlowProblemBlackoil.hpp:242
bool recycleFirstIterationStorage() const
Return if the storage term of the first iteration is identical to the storage term for the solution o...
Definition: FlowProblemBlackoil.hpp:908
LhsEval permFactTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Calculate the transmissibility multiplier due to porosity reduction.
Definition: FlowProblemBlackoil.hpp:720
void serializeOp(Serializer &serializer)
Definition: FlowProblemBlackoil.hpp:1145
MixingRateControls< FluidSystem > mixControls_
Definition: FlowProblemBlackoil.hpp:1713
void writeReports(const SimulatorTimer &timer)
Definition: FlowProblemBlackoil.hpp:546
ModuleParams moduleParams_
Definition: FlowProblemBlackoil.hpp:1717
const EclWriterType & eclWriter() const
Definition: FlowProblemBlackoil.hpp:872
void setSimulationReport(const SimulatorReport &report)
Definition: FlowProblemBlackoil.hpp:755
void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const override
Definition: FlowProblemBlackoil.hpp:613
void computeAndSetEqWeights_()
Definition: FlowProblemBlackoil.hpp:1185
void beginTimeStep() override
Called by the simulator before each time integration.
Definition: FlowProblemBlackoil.hpp:273
void updateExplicitQuantities_(const bool first_step_after_restart)
Definition: FlowProblemBlackoil.hpp:1602
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblemBlackoil.hpp:169
ActionHandler< Scalar, IndexTraits > actionHandler_
Definition: FlowProblemBlackoil.hpp:1715
void readSolutionFromOutputModule(const int restart_step, bool fip_init)
Read simulator solution state from the outputmodule (used with restart)
Definition: FlowProblemBlackoil.hpp:1022
const EclipseIO & eclIO() const
Definition: FlowProblemBlackoil.hpp:749
Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
Definition: FlowProblemBlackoil.hpp:1130
void handleUreaBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1592
EclWriterType & eclWriter()
Definition: FlowProblemBlackoil.hpp:875
bool satfuncConsistencyRequirementsMet() const
Definition: FlowProblemBlackoil.hpp:1626
bool updateCompositionChangeLimits_()
Definition: FlowProblemBlackoil.hpp:1219
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:95
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:133
@ numComponents
Definition: FlowProblem.hpp:118
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:486
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:836
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:669
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:108
@ gasPhaseIdx
Definition: FlowProblem.hpp:137
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:102
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:107
@ enablePolymer
Definition: FlowProblem.hpp:131
@ waterPhaseIdx
Definition: FlowProblem.hpp:139
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:151
@ numPhases
Definition: FlowProblem.hpp:117
GlobalEqVector drift_
Definition: FlowProblem.hpp:1722
@ enableExperiments
Definition: FlowProblem.hpp:127
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:148
@ enableFoam
Definition: FlowProblem.hpp:129
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:165
@ numEq
Definition: FlowProblem.hpp:116
int episodeIndex() const
Definition: FlowProblem.hpp:289
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:109
@ dimWorld
Definition: FlowProblem.hpp:113
@ gasCompIdx
Definition: FlowProblem.hpp:143
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:106
@ enableMICP
Definition: FlowProblem.hpp:130
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:149
@ enableBrine
Definition: FlowProblem.hpp:121
@ enableSolvent
Definition: FlowProblem.hpp:134
TracerModel tracerModel_
Definition: FlowProblem.hpp:1728
WellModel wellModel_
Definition: FlowProblem.hpp:1724
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:297
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:132
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:354
@ enableExtbo
Definition: FlowProblem.hpp:128
void readThermalParameters_()
Definition: FlowProblem.hpp:1401
@ waterCompIdx
Definition: FlowProblem.hpp:145
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:135
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:160
@ enableConvectiveMixing
Definition: FlowProblem.hpp:122
@ enableBioeffects
Definition: FlowProblem.hpp:120
TemperatureModel temperatureModel_
Definition: FlowProblem.hpp:1729
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:103
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:182
void updatePffDofData_()
Definition: FlowProblem.hpp:1515
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:147
@ oilCompIdx
Definition: FlowProblem.hpp:144
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1546
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1717
static constexpr EnergyModules energyModuleType
Definition: FlowProblem.hpp:125
@ enableDispersion
Definition: FlowProblem.hpp:124
@ oilPhaseIdx
Definition: FlowProblem.hpp:138
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:105
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:157
@ dim
Definition: FlowProblem.hpp:112
@ enableDiffusion
Definition: FlowProblem.hpp:123
void readMaterialParameters_()
Definition: FlowProblem.hpp:1367
This class calculates the threshold pressure for grid faces according to the Eclipse Reference Manual...
Definition: FlowThresholdPressure.hpp:59
Class handling mixing rate controls for a FlowProblemBlackoil.
Definition: MixingRateControls.hpp:46
Definition: SatfuncConsistencyCheckManager.hpp:58
SatfuncConsistencyCheckManager & collectFailuresTo(const int root)
Definition: SatfuncConsistencyCheckManager.hpp:99
void run(const GridView &gv, GetCellIndex &&getCellIndex)
Definition: SatfuncConsistencyCheckManager.hpp:128
Definition: SimulatorTimer.hpp:39
VTK output module for the tracer model's parameters.
Definition: VtkTracerModule.hpp:58
static void registerParameters()
Register all run-time parameters for the tracer VTK output module.
Definition: VtkTracerModule.hpp:84
@ NONE
Definition: DeferredLogger.hpp:46
static constexpr int dim
Definition: structuredgridvanguard.hh:68
Definition: blackoilbioeffectsmodules.hh:43
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
Struct holding the parameters for the BlackOilBioeffectsModule class.
Definition: blackoilbioeffectsparams.hpp:44
Struct holding the parameters for the BlackoilBrineModule class.
Definition: blackoilbrineparams.hpp:44
Struct holding the parameters for the BlackoilExtboModule class.
Definition: blackoilextboparams.hpp:49
Struct holding the parameters for the BlackoilFoamModule class.
Definition: blackoilfoamparams.hpp:46
Struct holding the parameters for the BlackOilPolymerModule class.
Definition: blackoilpolymerparams.hpp:45
Struct holding the parameters for the BlackOilSolventModule class.
Definition: blackoilsolventparams.hpp:49
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34