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>
41#include <opm/material/fluidsystems/blackoilpvt/ConstantRsDeadOilPvt.hpp>
47#include <opm/output/eclipse/EclipseIO.hpp>
49#include <opm/input/eclipse/Units/Units.hpp>
85template <
class TypeTag>
141 enum { enableDissolvedGas = Indices::compositionSwitchIdx >= 0 };
142 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
143 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
144 enum { enableGeochemistry = getPropValue<TypeTag, Properties::EnableGeochemistry>() };
145 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
147 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
148 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
149 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
152 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
153 using FoamModule = BlackOilFoamModule<TypeTag, enableFoam>;
154 using PolymerModule = BlackOilPolymerModule<TypeTag, enablePolymer>;
155 using SolventModule = BlackOilSolventModule<TypeTag, enableSolvent>;
157 using EclWriterType = EclWriter<TypeTag, OutputBlackOilModule<TypeTag> >;
158 using IndexTraits =
typename FluidSystem::IndexTraitsType;
160 using HybridNewton = BlackOilHybridNewton<TypeTag>;
161 using ModuleParams = BlackoilModuleParams<ConvectiveMixingModuleParam<Scalar>>;
164 using DamarisWriterType = DamarisWriter<TypeTag>;
180 DamarisWriterType::registerParameters();
193 simulator.vanguard().schedule(),
194 simulator.vanguard().actionState(),
195 simulator.vanguard().summaryState(),
197 simulator.vanguard().grid().comm())
203 const auto& vanguard = simulator.vanguard();
205 if constexpr (enableBrine) {
207 brineParams.template initFromState<enableBrine,
208 enableSaltPrecipitation>(vanguard.eclState());
209 BrineModule::setParams(std::move(brineParams));
212 if constexpr (enableDiffusion) {
213 DiffusionModule::initFromState(vanguard.eclState());
216 if constexpr (enableDispersion) {
217 DispersionModule::initFromState(vanguard.eclState());
220 if constexpr (enableExtbo) {
222 extboParams.template initFromState<enableExtbo>(vanguard.eclState());
223 ExtboModule::setParams(std::move(extboParams));
226 if constexpr (enableFoam) {
228 foamParams.template initFromState<enableFoam>(vanguard.eclState());
229 FoamModule::setParams(std::move(foamParams));
232 if constexpr (enableBioeffects) {
234 bioeffectsParams.template initFromState<enableBioeffects, enableMICP>(vanguard.eclState());
235 BioeffectsModule::setParams(std::move(bioeffectsParams));
238 if constexpr (enablePolymer) {
240 polymerParams.template initFromState<enablePolymer, enablePolymerMolarWeight>(vanguard.eclState());
241 PolymerModule::setParams(std::move(polymerParams));
244 if constexpr (enableSolvent) {
246 solventParams.template initFromState<enableSolvent>(vanguard.eclState(), vanguard.schedule());
247 SolventModule::setParams(std::move(solventParams));
251 eclWriter_ = std::make_unique<EclWriterType>(simulator);
255 if constexpr (!enableGeochemistry) {
256 if (vanguard.eclState().runspec().geochem().enabled()) {
257 throw std::runtime_error(
"GEOCHEM keyword in the deck but geochemistry module "
258 "was disabled at compile time!");
263 if constexpr (!enableMech) {
264 const auto& rspec = vanguard.eclState().runspec();
265 if (rspec.mech() && rspec.mechSolver().tpsa()) {
266 throw std::runtime_error(
"TPSA solver enabled in the deck, but geomechanics "
267 "module was disabled at compile time!");
273 damarisWriter_ = std::make_unique<DamarisWriterType>(simulator);
274 enableDamarisOutput_ = Parameters::Get<Parameters::EnableDamarisOutput>();
285 auto& simulator = this->simulator();
287 const int episodeIdx = simulator.episodeIndex();
288 const auto& schedule = simulator.vanguard().schedule();
293 .evalUDQAssignments(episodeIdx, simulator.vanguard().udqState());
295 if (episodeIdx >= 0) {
296 const auto& oilVap = schedule[episodeIdx].oilvap();
297 if (oilVap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS) {
298 FluidSystem::setVapPars(oilVap.vap1(), oilVap.vap2());
301 FluidSystem::setVapPars(0.0, 0.0);
304 if constexpr (enableConvectiveMixing) {
305 ConvectiveMixingModule::beginEpisode(simulator.vanguard().eclState(), schedule, episodeIdx,
306 this->moduleParams_.convectiveMixingModuleParam);
328 FlowProblemType::finishInit();
330 auto& simulator = this->simulator();
332 auto finishTransmissibilities = [updated =
false,
this]()
mutable
334 if (updated) {
return; }
336 this->
transmissibilities_.finishInit([&vg = this->simulator().vanguard()](
const unsigned int it) {
337 return vg.gridIdxToEquilGridIdx(it);
348 if (simulator.vanguard().grid().comm().size() > 1) {
349 if (simulator.vanguard().grid().comm().rank() == 0)
350 eclWriter_->setTransmissibilities(&simulator.vanguard().globalTransmissibility());
352 finishTransmissibilities();
353 eclWriter_->setTransmissibilities(&simulator.problem().eclTransmissibilities());
356 std::function<
unsigned int(
unsigned int)> equilGridToGrid = [&simulator](
unsigned int i) {
357 return simulator.vanguard().gridEquilIdxToGridIdx(i);
360 this->
eclWriter_->extractOutputTransAndNNC(equilGridToGrid);
362 simulator.vanguard().releaseGlobalTransmissibilities();
364 const auto& eclState = simulator.vanguard().eclState();
365 const auto& schedule = simulator.vanguard().schedule();
368 simulator.setStartTime(schedule.getStartTime());
369 simulator.setEndTime(schedule.simTime(schedule.size() - 1));
375 simulator.setEpisodeIndex(-1);
376 simulator.setEpisodeLength(0.0);
381 this->gravity_ = 0.0;
382 if (Parameters::Get<Parameters::EnableGravity>() &&
383 eclState.getInitConfig().hasGravity())
386 this->gravity_[dim - 1] = unit::gravity;
392 const auto& tuning = schedule[0].tuning();
399 bool isThermal = eclState.getSimulationConfig().isThermal();
400 bool isTemp = eclState.getSimulationConfig().isTemp();
401 bool conserveInnerEnergy = isTemp || (isThermal && Parameters::Get<Parameters::ConserveInnerEnergyThermal>());
402 FluidSystem::setEnergyEqualEnthalpy(conserveInnerEnergy);
404 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
405 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
410 [&simulator](
const unsigned idx)
412 std::array<int,dim> coords;
413 simulator.vanguard().cartesianCoordinate(idx, coords);
414 std::ranges::transform(coords, coords.begin(),
415 [](const auto c) { return c + 1; });
427 finishTransmissibilities();
429 const auto& initconfig = eclState.getInitConfig();
431 if (initconfig.restartRequested()) {
442 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>()) {
443 const auto& vanguard = this->simulator().vanguard();
444 const auto& gridView = vanguard.gridView();
445 const int numElements = gridView.size(0);
446 this->
polymer_.maxAdsorption.resize(numElements, 0.0);
455 this->
drift_.resize(this->model().numGridDof());
462 if (!initconfig.restartRequested() && !eclState.getIOConfig().initOnly()) {
463 simulator.startNextEpisode(schedule.seconds(1));
464 simulator.setEpisodeIndex(0);
465 simulator.setTimeStepIndex(0);
468 if (Parameters::Get<Parameters::CheckSatfuncConsistency>() &&
478 this->simulator().vanguard().grid().comm().barrier();
480 throw std::domain_error {
481 "Saturation function end-points do not "
482 "meet requisite consistency conditions"
491 eclState.runspec().tabdims().getNumPVTTables());
493 if (this->enableVtkOutput_() && eclState.getIOConfig().initOnly()) {
494 simulator.setTimeStepSize(0.0);
495 simulator.model().applyInitialSolution();
499 if (!eclState.getIOConfig().initOnly()) {
500 if (!this->
enableTuning_ && eclState.getSimulationConfig().anyTUNING()) {
501 OpmLog::info(
"\nThe deck has TUNING in the SCHEDULE section, but "
502 "it is ignored due\nto the flag --enable-tuning=false. "
503 "Set this flag to true to activate it.\n"
504 "Manually tuning the simulator with the TUNING keyword may "
505 "increase run time.\nIt is recommended using the simulator's "
506 "default tuning (--enable-tuning=false).");
516 FlowProblemType::endTimeStep();
517 this->endStepApplyAction();
524 this->eclWriter().mutableOutputModule().invalidateLocalData();
527 const auto& grid = this->simulator().vanguard().gridView().grid();
529 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
530 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
531 if (!isCpGrid || (grid.maxLevel() == 0)) {
532 this->eclWriter_->evalSummaryState(!this->episodeWillBeOver());
536 OPM_TIMEBLOCK(applyActions);
538 const int episodeIdx = this->episodeIndex();
539 auto& simulator = this->simulator();
543 this->simulator().vanguard().schedule().clearEvents(episodeIdx);
547 .applyActions(episodeIdx, simulator.time() + simulator.timeStepSize(),
548 [
this](
const bool global)
550 using TransUpdateQuantities = typename
551 Vanguard::TransmissibilityType::TransUpdateQuantities;
553 this->transmissibilities_
554 .update(global, TransUpdateQuantities::All,
555 [&vg = this->simulator().vanguard()]
556 (const unsigned int i)
558 return vg.gridIdxToEquilGridIdx(i);
569 OPM_TIMEBLOCK(endEpisode);
582 .evalUDQAssignments(this->episodeIndex(), this->simulator().vanguard().udqState());
584 FlowProblemType::endEpisode();
589 if (this->enableEclOutput_) {
590 this->eclWriter_->writeReports(timer);
601 FlowProblemType::writeOutput(verbose);
603 const auto isSubStep = !this->episodeWillBeOver();
605 auto localCellData = data::Solution {};
610 if (this->enableDamarisOutput_ && (this->damarisWriter_ !=
nullptr)) {
611 this->damarisWriter_->writeOutput(localCellData, isSubStep);
615 if (this->enableEclOutput_ && (this->eclWriter_ !=
nullptr)) {
616 this->eclWriter_->writeOutput(std::move(localCellData), isSubStep,
617 this->simulator().vanguard().schedule()
618 .exitStatus().has_value());
624 OPM_TIMEBLOCK(finalizeOutput);
636 FlowProblemType::initialSolutionApplied();
641 this->thresholdPressures_.finishInit();
644 const auto& grid = this->simulator().vanguard().gridView().grid();
646 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
647 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
649 if (!isCpGrid || (grid.maxLevel() == 0)) {
650 if (this->simulator().episodeIndex() == 0) {
651 eclWriter_->writeInitialFIPReport();
657 unsigned globalDofIdx,
658 unsigned timeIdx)
const override
660 this->aquiferModel_.addToSource(rate, globalDofIdx, timeIdx);
663 const auto& source = this->simulator().vanguard().schedule()[this->episodeIndex()].source();
664 std::array<int,3> ijk;
665 this->simulator().vanguard().cartesianCoordinate(globalDofIdx, ijk);
667 if (source.hasSource(ijk)) {
668 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
669 static std::array<SourceComponent, 3> sc_map = {SourceComponent::WATER, SourceComponent::OIL, SourceComponent::GAS};
670 static std::array<int, 3> phidx_map = {FluidSystem::waterPhaseIdx, FluidSystem::oilPhaseIdx, FluidSystem::gasPhaseIdx};
671 static std::array<int, 3> cidx_map = {waterCompIdx, oilCompIdx, gasCompIdx};
673 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
674 const auto phaseIdx = phidx_map[i];
675 const auto sourceComp = sc_map[i];
676 const auto compIdx = cidx_map[i];
677 if (!FluidSystem::phaseIsActive(phaseIdx)) {
680 Scalar mass_rate = source.rate(ijk, sourceComp) / this->model().dofTotalVolume(globalDofIdx);
681 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
682 mass_rate /= FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
684 rate[FluidSystem::canonicalToActiveCompIdx(compIdx)] += mass_rate;
687 if constexpr (enableSolvent) {
688 Scalar mass_rate = source.rate(ijk, SourceComponent::SOLVENT) / this->model().dofTotalVolume(globalDofIdx);
689 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
690 const auto& solventPvt = SolventModule::solventPvt();
691 mass_rate /= solventPvt.referenceDensity(pvtRegionIdx);
693 rate[Indices::contiSolventEqIdx] += mass_rate;
695 if constexpr (enablePolymer) {
696 rate[Indices::polymerConcentrationIdx] += source.rate(ijk, SourceComponent::POLYMER) / this->model().dofTotalVolume(globalDofIdx);
698 if constexpr (enableMICP) {
699 rate[Indices::microbialConcentrationIdx] += source.rate(ijk, SourceComponent::MICR) / this->model().dofTotalVolume(globalDofIdx);
700 rate[Indices::oxygenConcentrationIdx] += source.rate(ijk, SourceComponent::OXYG) / this->model().dofTotalVolume(globalDofIdx);
701 rate[Indices::ureaConcentrationIdx] += source.rate(ijk, SourceComponent::UREA) / (this->model().dofTotalVolume(globalDofIdx));
703 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
704 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
705 const auto phaseIdx = phidx_map[i];
706 if (!FluidSystem::phaseIsActive(phaseIdx)) {
709 const auto sourceComp = sc_map[i];
710 const auto source_hrate = source.hrate(ijk, sourceComp);
712 rate[Indices::contiEnergyEqIdx] += source_hrate.value() / this->model().dofTotalVolume(globalDofIdx);
714 const auto& intQuants = this->simulator().model().intensiveQuantities(globalDofIdx, 0);
715 auto fs = intQuants.fluidState();
717 const auto source_temp = source.temperature(ijk, sourceComp);
719 Scalar temperature = source_temp.value();
720 fs.setTemperature(temperature);
722 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, pvtRegionIdx);
723 Scalar mass_rate = source.rate(ijk, sourceComp)/ this->model().dofTotalVolume(globalDofIdx);
724 Scalar energy_rate = getValue(h)*mass_rate;
725 rate[Indices::contiEnergyEqIdx] += energy_rate;
733 if (this->enableDriftCompensation_) {
734 const auto& simulator = this->simulator();
735 const auto& model = this->model();
740 Scalar maxCompensation = model.newtonMethod().tolerance()/10;
741 Scalar poro = this->porosity(globalDofIdx, timeIdx);
742 Scalar dt = simulator.timeStepSize();
743 EqVector dofDriftRate = this->drift_[globalDofIdx];
744 dofDriftRate /= dt*model.dofTotalVolume(globalDofIdx);
747 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
748 Scalar cnv = std::abs(dofDriftRate[eqIdx])*dt*model.eqWeight(globalDofIdx, eqIdx)/poro;
749 if (cnv > maxCompensation) {
750 dofDriftRate[eqIdx] *= maxCompensation/cnv;
754 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx)
755 rate[eqIdx] -= dofDriftRate[eqIdx];
762 template <
class LhsEval,
class Callback>
765 OPM_TIMEBLOCK_LOCAL(permFactTransMultiplier, Subsystem::PvtProps);
766 if constexpr (enableSaltPrecipitation) {
767 const auto& fs = intQuants.fluidState();
768 unsigned tableIdx = this->simulator().problem().satnumRegionIndex(elementIdx);
769 LhsEval porosityFactor = obtain(1. - fs.saltSaturation());
770 porosityFactor = min(porosityFactor, 1.0);
771 const auto& permfactTable = BrineModule::permfactTable(tableIdx);
772 return permfactTable.eval(porosityFactor,
true);
774 else if constexpr (enableBioeffects) {
775 return obtain(intQuants.permFactor());
784 {
return initialFluidStates_[globalDofIdx]; }
787 {
return initialFluidStates_; }
790 {
return initialFluidStates_; }
793 {
return eclWriter_->eclIO(); }
796 {
return eclWriter_->setSubStepReport(report); }
799 {
return eclWriter_->setSimulationReport(report); }
803 OPM_TIMEBLOCK_LOCAL(boundaryFluidState, Subsystem::Assembly);
804 const auto& bcprop = this->simulator().vanguard().schedule()[this->episodeIndex()].bcprop;
805 if (bcprop.size() > 0) {
806 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
810 if (this->bcindex_(dir)[globalDofIdx] == 0)
811 return initialFluidStates_[globalDofIdx];
813 const auto& bc = bcprop[this->bcindex_(dir)[globalDofIdx]];
814 if (bc.bctype == BCType::DIRICHLET )
816 InitialFluidState fluidState;
817 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
818 fluidState.setPvtRegionIndex(pvtRegionIdx);
820 switch (bc.component) {
821 case BCComponent::OIL:
822 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx))
823 throw std::logic_error(
"oil is not active and you're trying to add oil BC");
825 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0);
827 case BCComponent::GAS:
828 if (!FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
829 throw std::logic_error(
"gas is not active and you're trying to add gas BC");
831 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 1.0);
833 case BCComponent::WATER:
834 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
835 throw std::logic_error(
"water is not active and you're trying to add water BC");
837 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0);
839 case BCComponent::SOLVENT:
840 case BCComponent::POLYMER:
841 case BCComponent::MICR:
842 case BCComponent::OXYG:
843 case BCComponent::UREA:
845 throw std::logic_error(
"you need to specify a valid component (OIL, WATER or GAS) when DIRICHLET type is set in BC");
847 fluidState.setTotalSaturation(1.0);
848 double pressure = initialFluidStates_[globalDofIdx].pressure(this->refPressurePhaseIdx_());
849 const auto pressure_input = bc.pressure;
850 if (pressure_input) {
851 pressure = *pressure_input;
854 std::array<Scalar, numPhases> pc = {0};
855 const auto& matParams = this->materialLawParams(globalDofIdx);
856 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
857 Valgrind::CheckDefined(pressure);
858 Valgrind::CheckDefined(pc);
859 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
860 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
861 if (Indices::oilEnabled)
862 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
863 else if (Indices::gasEnabled)
864 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
865 else if (Indices::waterEnabled)
867 fluidState.setPressure(phaseIdx, pressure);
869 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
870 double temperature = initialFluidStates_[globalDofIdx].temperature(0);
871 const auto temperature_input = bc.temperature;
872 if(temperature_input)
873 temperature = *temperature_input;
874 fluidState.setTemperature(temperature);
877 if constexpr (enableDissolvedGas) {
878 if (FluidSystem::enableDissolvedGas()) {
879 fluidState.setRs(0.0);
880 fluidState.setRv(0.0);
883 if constexpr (enableDisgasInWater) {
884 if (FluidSystem::enableDissolvedGasInWater()) {
885 fluidState.setRsw(0.0);
888 if constexpr (enableVapwat) {
889 if (FluidSystem::enableVaporizedWater()) {
890 fluidState.setRvw(0.0);
894 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
895 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
897 const auto& b = FluidSystem::inverseFormationVolumeFactor(fluidState, phaseIdx, pvtRegionIdx);
898 fluidState.setInvB(phaseIdx, b);
900 const auto& rho = FluidSystem::density(fluidState, phaseIdx, pvtRegionIdx);
901 fluidState.setDensity(phaseIdx, rho);
902 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
903 const auto& h = FluidSystem::enthalpy(fluidState, phaseIdx, pvtRegionIdx);
904 fluidState.setEnthalpy(phaseIdx, h);
907 fluidState.checkDefined();
911 return initialFluidStates_[globalDofIdx];
916 {
return *eclWriter_; }
919 {
return *eclWriter_; }
927 return this->mixControls_.maxGasDissolutionFactor(timeIdx, globalDofIdx,
928 this->episodeIndex(),
929 this->pvtRegionIndex(globalDofIdx));
938 return this->mixControls_.maxOilVaporizationFactor(timeIdx, globalDofIdx,
939 this->episodeIndex(),
940 this->pvtRegionIndex(globalDofIdx));
953 int episodeIdx = this->episodeIndex();
954 return !this->mixControls_.drsdtActive(episodeIdx) &&
955 !this->mixControls_.drvdtActive(episodeIdx) &&
956 this->rockCompPoroMultWc_.empty() &&
957 this->rockCompPoroMult_.empty();
966 template <
class Context>
969 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
971 values.setPvtRegionIndex(pvtRegionIndex(context, spaceIdx, timeIdx));
972 values.assignNaive(initialFluidStates_[globalDofIdx]);
974 if constexpr (enableSolvent) {
975 SolventModule::assignPrimaryVars(values,
976 this->solventSaturation_[globalDofIdx],
977 this->solventRsw_[globalDofIdx]);
980 if constexpr (enablePolymer) {
981 values[Indices::polymerConcentrationIdx] = this->polymer_.concentration[globalDofIdx];
984 if constexpr (enablePolymerMolarWeight) {
985 values[Indices::polymerMoleWeightIdx]= this->polymer_.moleWeight[globalDofIdx];
988 if constexpr (enableBrine) {
989 if (enableSaltPrecipitation && values.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
990 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltSaturation();
993 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltConcentration();
997 if constexpr (enableBioeffects) {
998 values[Indices::microbialConcentrationIdx] = this->bioeffects_.microbialConcentration[globalDofIdx];
999 values[Indices::biofilmVolumeFractionIdx]= this->bioeffects_.biofilmVolumeFraction[globalDofIdx];
1000 if constexpr (enableMICP) {
1001 values[Indices::oxygenConcentrationIdx]= this->bioeffects_.oxygenConcentration[globalDofIdx];
1002 values[Indices::ureaConcentrationIdx]= this->bioeffects_.ureaConcentration[globalDofIdx];
1003 values[Indices::calciteVolumeFractionIdx]= this->bioeffects_.calciteVolumeFraction[globalDofIdx];
1007 values.checkDefined();
1013 return this->mixControls_.drsdtcon(elemIdx, episodeIdx,
1014 this->pvtRegionIndex(elemIdx));
1019 return this->mixControls_.drsdtConvective(episodeIdx, this->pvtRegionIndex(elemIdx));
1027 template <
class Context>
1029 const Context& context,
1031 unsigned timeIdx)
const
1033 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
1034 if (!context.intersection(spaceIdx).boundary())
1037 if constexpr (energyModuleType != EnergyModules::FullyImplicitThermal || !enableThermalFluxBoundaries)
1045 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1046 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1047 values.setThermalFlow(context, spaceIdx, timeIdx, this->initialFluidStates_[globalDofIdx] );
1050 if (this->nonTrivialBoundaryConditions()) {
1051 unsigned indexInInside = context.intersection(spaceIdx).indexInInside();
1052 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1053 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1054 unsigned pvtRegionIdx = pvtRegionIndex(context, spaceIdx, timeIdx);
1055 const auto [type, massrate] = this->boundaryCondition(globalDofIdx, indexInInside);
1056 if (type == BCType::THERMAL)
1057 values.setThermalFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1058 else if (type == BCType::FREE || type == BCType::DIRICHLET)
1059 values.setFreeFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1060 else if (type == BCType::RATE)
1061 values.setMassRate(massrate, pvtRegionIdx);
1071 auto& simulator = this->simulator();
1072 const auto& eclState = simulator.vanguard().eclState();
1074 std::size_t numElems = this->model().numGridDof();
1075 this->initialFluidStates_.resize(numElems);
1076 if constexpr (enableSolvent) {
1077 this->solventSaturation_.resize(numElems, 0.0);
1078 this->solventRsw_.resize(numElems, 0.0);
1081 if constexpr (enablePolymer)
1082 this->polymer_.concentration.resize(numElems, 0.0);
1084 if constexpr (enablePolymerMolarWeight) {
1085 const std::string msg {
"Support of the RESTART for polymer molecular weight "
1086 "is not implemented yet. The polymer weight value will be "
1087 "zero when RESTART begins"};
1088 OpmLog::warning(
"NO_POLYMW_RESTART", msg);
1089 this->polymer_.moleWeight.resize(numElems, 0.0);
1092 if constexpr (enableBioeffects) {
1093 this->bioeffects_.resize(numElems);
1097 this->mixControls_.init(numElems, restart_step, eclState.runspec().tabdims().getNumPVTTables());
1099 if constexpr (enableBioeffects) {
1100 this->bioeffects_ = this->eclWriter_->outputModule().getBioeffects().getSolution();
1103 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1104 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1105 elemFluidState.setPvtRegionIndex(pvtRegionIndex(elemIdx));
1106 this->eclWriter_->outputModule().initHysteresisParams(simulator, elemIdx);
1107 this->eclWriter_->outputModule().assignToFluidState(elemFluidState, elemIdx);
1116 auto ssol = enableSolvent
1117 ? this->eclWriter_->outputModule().getSolventSaturation(elemIdx)
1120 this->processRestartSaturations_(elemFluidState, ssol);
1122 if constexpr (enableSolvent) {
1123 this->solventSaturation_[elemIdx] = ssol;
1124 this->solventRsw_[elemIdx] = this->eclWriter_->outputModule().getSolventRsw(elemIdx);
1129 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1130 bool needTemperature = (eclState.runspec().co2Storage() || eclState.runspec().h2Storage());
1131 if (needTemperature) {
1132 const auto& fp = simulator.vanguard().eclState().fieldProps();
1133 elemFluidState.setTemperature(fp.get_double(
"TEMPI")[elemIdx]);
1137 this->mixControls_.updateLastValues(elemIdx, elemFluidState.Rs(), elemFluidState.Rv());
1139 if constexpr (enablePolymer)
1140 this->polymer_.concentration[elemIdx] = this->eclWriter_->outputModule().getPolymerConcentration(elemIdx);
1144 const int episodeIdx = this->episodeIndex();
1145 this->mixControls_.updateMaxValues(episodeIdx, simulator.timeStepSize());
1150 auto& sol = this->model().solution(0);
1151 const auto& gridView = this->gridView();
1153 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1154 elemCtx.updatePrimaryStencil(elem);
1155 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1156 this->initial(sol[elemIdx], elemCtx, 0, 0);
1164 this->model().syncOverlap();
1167 this->updateReferencePorosity_();
1168 this->mixControls_.init(this->model().numGridDof(),
1169 this->episodeIndex(),
1170 eclState.runspec().tabdims().getNumPVTTables());
1178 {
return thresholdPressures_.thresholdPressure(elem1Idx, elem2Idx); }
1181 {
return thresholdPressures_; }
1184 {
return thresholdPressures_; }
1188 return moduleParams_;
1191 template<
class Serializer>
1195 serializer(mixControls_);
1196 serializer(*eclWriter_);
1202 this->updateExplicitQuantities_(first_step_after_restart);
1204 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>())
1205 updateMaxPolymerAdsorption_();
1207 mixControls_.updateExplicitQuantities(episodeIdx, timeStepSize);
1213 this->updateProperty_(
"FlowProblemBlackoil::updateMaxPolymerAdsorption_() failed:",
1216 this->updateMaxPolymerAdsorption_(compressedDofIdx,iq);
1222 const Scalar pa = scalarValue(iq.polymerAdsorption());
1223 auto& mpa = this->polymer_.maxAdsorption;
1224 if (mpa[compressedDofIdx] < pa) {
1225 mpa[compressedDofIdx] = pa;
1234 std::vector<Scalar> sumInvB(numPhases, 0.0);
1235 const auto& gridView = this->gridView();
1237 for(
const auto& elem: elements(gridView, Dune::Partitions::interior)) {
1238 elemCtx.updatePrimaryStencil(elem);
1239 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1240 const auto& dofFluidState = this->initialFluidStates_[elemIdx];
1241 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1242 if (!FluidSystem::phaseIsActive(phaseIdx))
1245 sumInvB[phaseIdx] += dofFluidState.invB(phaseIdx);
1249 std::size_t numDof = this->model().numGridDof();
1250 const auto& comm = this->simulator().vanguard().grid().comm();
1251 comm.sum(sumInvB.data(),sumInvB.size());
1252 Scalar numTotalDof = comm.sum(numDof);
1254 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1255 if (!FluidSystem::phaseIsActive(phaseIdx))
1258 Scalar avgB = numTotalDof / sumInvB[phaseIdx];
1259 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
1260 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
1261 this->model().setEqWeight(activeSolventCompIdx, avgB);
1268 OPM_TIMEBLOCK(updateCompositionChangeLimits);
1271 int episodeIdx = this->episodeIndex();
1272 std::array<bool,3> active{this->mixControls_.drsdtConvective(episodeIdx),
1273 this->mixControls_.drsdtActive(episodeIdx),
1274 this->mixControls_.drvdtActive(episodeIdx)};
1275 if (!active[0] && !active[1] && !active[2]) {
1279 this->updateProperty_(
"FlowProblemBlackoil::updateCompositionChangeLimits_()) failed:",
1280 [
this,episodeIdx,active](
unsigned compressedDofIdx,
1283 const DimMatrix& perm = this->intrinsicPermeability(compressedDofIdx);
1284 const Scalar distZ = active[0] ? this->simulator().vanguard().cellThickness(compressedDofIdx) : 0.0;
1285 const int pvtRegionIdx = this->pvtRegionIndex(compressedDofIdx);
1286 this->mixControls_.update(compressedDofIdx,
1289 this->gravity_[
dim - 1],
1302 if(this->simulator().vanguard().grid().maxLevel() > 0) {
1303 throw std::invalid_argument(
"Refined grids are not yet supported for restart ");
1307 auto& simulator = this->simulator();
1308 const auto& schedule = simulator.vanguard().schedule();
1309 const auto& eclState = simulator.vanguard().eclState();
1310 const auto& initconfig = eclState.getInitConfig();
1311 const int restart_step = initconfig.getRestartStep();
1313 simulator.setTime(schedule.seconds(restart_step));
1315 simulator.startNextEpisode(simulator.startTime() + simulator.time(),
1316 schedule.stepLength(restart_step));
1317 simulator.setEpisodeIndex(restart_step);
1319 this->eclWriter_->beginRestart();
1321 Scalar dt = std::min(this->eclWriter_->restartTimeStepSize(), simulator.episodeLength());
1322 simulator.setTimeStepSize(dt);
1324 this->readSolutionFromOutputModule(restart_step,
false);
1326 this->eclWriter_->endRestart();
1331 const auto& simulator = this->simulator();
1336 std::size_t numElems = this->model().numGridDof();
1337 this->initialFluidStates_.resize(numElems);
1338 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1339 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1346 const auto& simulator = this->simulator();
1347 const auto& vanguard = simulator.vanguard();
1348 const auto& eclState = vanguard.eclState();
1349 const auto& fp = eclState.fieldProps();
1350 bool has_swat = fp.has_double(
"SWAT");
1351 bool has_sgas = fp.has_double(
"SGAS");
1352 bool has_rs = fp.has_double(
"RS");
1353 bool has_rsw = fp.has_double(
"RSW");
1354 bool has_rv = fp.has_double(
"RV");
1355 bool has_rvw = fp.has_double(
"RVW");
1356 bool has_pressure = fp.has_double(
"PRESSURE");
1357 bool has_salt = fp.has_double(
"SALT");
1358 bool has_saltp = fp.has_double(
"SALTP");
1361 if (Indices::numPhases > 1) {
1362 if (FluidSystem::phaseIsActive(waterPhaseIdx) && !has_swat)
1363 throw std::runtime_error(
"The ECL input file requires the presence of the SWAT keyword if "
1364 "the water phase is active");
1365 if (FluidSystem::phaseIsActive(gasPhaseIdx) && !has_sgas && FluidSystem::phaseIsActive(oilPhaseIdx))
1366 throw std::runtime_error(
"The ECL input file requires the presence of the SGAS keyword if "
1367 "the gas phase is active");
1370 throw std::runtime_error(
"The ECL input file requires the presence of the PRESSURE "
1371 "keyword if the model is initialized explicitly");
1372 if (FluidSystem::enableDissolvedGas() && !has_rs)
1373 throw std::runtime_error(
"The ECL input file requires the RS keyword to be present if"
1374 " dissolved gas is enabled and the model is initialized explicitly");
1375 if (FluidSystem::enableDissolvedGasInWater() && !has_rsw)
1376 OpmLog::warning(
"The model is initialized explicitly and the RSW keyword is not present in the"
1377 " ECL input file. The RSW values are set equal to 0");
1378 if (FluidSystem::enableVaporizedOil() && !has_rv)
1379 throw std::runtime_error(
"The ECL input file requires the RV keyword to be present if"
1380 " vaporized oil is enabled and the model is initialized explicitly");
1381 if (FluidSystem::enableVaporizedWater() && !has_rvw)
1382 throw std::runtime_error(
"The ECL input file requires the RVW keyword to be present if"
1383 " vaporized water is enabled and the model is initialized explicitly");
1384 if (enableBrine && !has_salt)
1385 throw std::runtime_error(
"The ECL input file requires the SALT keyword to be present if"
1386 " brine is enabled and the model is initialized explicitly");
1387 if (enableSaltPrecipitation && !has_saltp)
1388 throw std::runtime_error(
"The ECL input file requires the SALTP keyword to be present if"
1389 " salt precipitation is enabled and the model is initialized explicitly");
1391 std::size_t numDof = this->model().numGridDof();
1393 initialFluidStates_.resize(numDof);
1395 std::vector<double> waterSaturationData;
1396 std::vector<double> gasSaturationData;
1397 std::vector<double> pressureData;
1398 std::vector<double> rsData;
1399 std::vector<double> rswData;
1400 std::vector<double> rvData;
1401 std::vector<double> rvwData;
1402 std::vector<double> tempiData;
1403 std::vector<double> saltData;
1404 std::vector<double> saltpData;
1406 if (FluidSystem::phaseIsActive(waterPhaseIdx) && Indices::numPhases > 1)
1407 waterSaturationData = fp.get_double(
"SWAT");
1409 waterSaturationData.resize(numDof);
1411 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx))
1412 gasSaturationData = fp.get_double(
"SGAS");
1414 gasSaturationData.resize(numDof);
1416 pressureData = fp.get_double(
"PRESSURE");
1417 if (FluidSystem::enableDissolvedGas())
1418 rsData = fp.get_double(
"RS");
1420 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1421 rswData = fp.get_double(
"RSW");
1423 if (FluidSystem::enableVaporizedOil())
1424 rvData = fp.get_double(
"RV");
1426 if (FluidSystem::enableVaporizedWater())
1427 rvwData = fp.get_double(
"RVW");
1430 tempiData = fp.get_double(
"TEMPI");
1433 if constexpr (enableBrine)
1434 saltData = fp.get_double(
"SALT");
1437 if constexpr (enableSaltPrecipitation)
1438 saltpData = fp.get_double(
"SALTP");
1441 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1442 auto& dofFluidState = initialFluidStates_[dofIdx];
1444 dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
1449 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1450 Scalar temperatureLoc = tempiData[dofIdx];
1451 if (!std::isfinite(temperatureLoc) || temperatureLoc <= 0)
1452 temperatureLoc = FluidSystem::surfaceTemperature;
1453 dofFluidState.setTemperature(temperatureLoc);
1459 if constexpr (enableBrine)
1460 dofFluidState.setSaltConcentration(saltData[dofIdx]);
1465 if constexpr (enableSaltPrecipitation)
1466 dofFluidState.setSaltSaturation(saltpData[dofIdx]);
1471 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
1472 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
1473 waterSaturationData[dofIdx]);
1475 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
1476 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
1477 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1479 - waterSaturationData[dofIdx]);
1482 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1483 gasSaturationData[dofIdx]);
1485 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1486 const Scalar soil = 1.0 - waterSaturationData[dofIdx] - gasSaturationData[dofIdx];
1487 if (soil < smallSaturationTolerance_) {
1488 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
1491 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, soil);
1498 Scalar pressure = pressureData[dofIdx];
1502 std::array<Scalar, numPhases> pc = {0};
1503 const auto& matParams = this->materialLawParams(dofIdx);
1504 MaterialLaw::capillaryPressures(pc, matParams, dofFluidState);
1505 Valgrind::CheckDefined(pressure);
1506 Valgrind::CheckDefined(pc);
1507 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1508 if (!FluidSystem::phaseIsActive(phaseIdx))
1511 if (Indices::oilEnabled)
1512 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
1513 else if (Indices::gasEnabled)
1514 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
1515 else if (Indices::waterEnabled)
1517 dofFluidState.setPressure(phaseIdx, pressure);
1520 if constexpr (enableDissolvedGas) {
1521 if (FluidSystem::enableDissolvedGas())
1522 dofFluidState.setRs(rsData[dofIdx]);
1523 else if (Indices::gasEnabled && Indices::oilEnabled)
1524 dofFluidState.setRs(0.0);
1525 if (FluidSystem::enableVaporizedOil())
1526 dofFluidState.setRv(rvData[dofIdx]);
1527 else if (Indices::gasEnabled && Indices::oilEnabled)
1528 dofFluidState.setRv(0.0);
1531 if constexpr (enableDisgasInWater) {
1532 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1533 dofFluidState.setRsw(rswData[dofIdx]);
1536 if constexpr (enableVapwat) {
1537 if (FluidSystem::enableVaporizedWater())
1538 dofFluidState.setRvw(rvwData[dofIdx]);
1544 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1545 if (!FluidSystem::phaseIsActive(phaseIdx))
1548 const auto& b = FluidSystem::inverseFormationVolumeFactor(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1549 dofFluidState.setInvB(phaseIdx, b);
1551 const auto& rho = FluidSystem::density(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1552 dofFluidState.setDensity(phaseIdx, rho);
1563 Scalar sumSaturation = 0.0;
1564 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1565 if (FluidSystem::phaseIsActive(phaseIdx)) {
1566 if (elemFluidState.saturation(phaseIdx) < smallSaturationTolerance_)
1567 elemFluidState.setSaturation(phaseIdx, 0.0);
1569 sumSaturation += elemFluidState.saturation(phaseIdx);
1573 if constexpr (enableSolvent) {
1574 if (solventSaturation < smallSaturationTolerance_)
1575 solventSaturation = 0.0;
1577 sumSaturation += solventSaturation;
1580 assert(sumSaturation > 0.0);
1582 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1583 if (FluidSystem::phaseIsActive(phaseIdx)) {
1584 const Scalar saturation = elemFluidState.saturation(phaseIdx) / sumSaturation;
1585 elemFluidState.setSaturation(phaseIdx, saturation);
1588 if constexpr (enableSolvent) {
1589 solventSaturation = solventSaturation / sumSaturation;
1595 FlowProblemType::readInitialCondition_();
1597 if constexpr (enableSolvent || enablePolymer || enablePolymerMolarWeight || enableBioeffects)
1598 this->readBlackoilExtentionsInitialConditions_(this->model().numGridDof(),
1601 enablePolymerMolarWeight,
1609 if constexpr (!enableSolvent)
1610 throw std::logic_error(
"solvent is disabled and you're trying to add solvent to BC");
1612 rate[Indices::solventSaturationIdx] = bc.rate;
1617 if constexpr (!enablePolymer)
1618 throw std::logic_error(
"polymer is disabled and you're trying to add polymer to BC");
1620 rate[Indices::polymerConcentrationIdx] = bc.rate;
1625 if constexpr (!enableMICP)
1626 throw std::logic_error(
"MICP is disabled and you're trying to add microbes to BC");
1628 rate[Indices::microbialConcentrationIdx] = bc.rate;
1633 if constexpr (!enableMICP)
1634 throw std::logic_error(
"MICP is disabled and you're trying to add oxygen to BC");
1636 rate[Indices::oxygenConcentrationIdx] = bc.rate;
1641 if constexpr (!enableMICP)
1642 throw std::logic_error(
"MICP is disabled and you're trying to add urea to BC");
1644 rate[Indices::ureaConcentrationIdx] = bc.rate;
1646 rate[Indices::ureaConcentrationIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
1651 OPM_TIMEBLOCK(updateExplicitQuantities);
1652 const bool invalidateFromMaxWaterSat = this->updateMaxWaterSaturation_();
1653 const bool invalidateFromMinPressure = this->updateMinPressure_();
1656 const bool invalidateFromHyst = this->updateHysteresis_();
1657 const bool invalidateFromMaxOilSat = this->updateMaxOilSaturation_();
1660 const bool invalidateDRDT = !first_step_after_restart && this->updateCompositionChangeLimits_();
1663 const bool invalidateIntensiveQuantities
1664 = invalidateFromMaxWaterSat || invalidateFromMinPressure || invalidateFromHyst || invalidateFromMaxOilSat || invalidateDRDT;
1665 if (invalidateIntensiveQuantities) {
1666 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
1667 this->model().invalidateAndUpdateIntensiveQuantities(0);
1670 this->updateRockCompTransMultVal_();
1675 if (
const auto nph = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
1676 + FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
1677 + FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1686 const auto numSamplePoints =
static_cast<std::size_t
>
1687 (Parameters::Get<Parameters::NumSatfuncConsistencySamplePoints>());
1689 auto sfuncConsistencyChecks =
1691 numSamplePoints, this->simulator().vanguard().eclState(),
1692 [&cmap = this->simulator().vanguard().cartesianIndexMapper()](
const int elemIdx)
1693 {
return cmap.cartesianIndex(elemIdx); }
1696 const auto ioRank = 0;
1697 const auto isIoRank = this->simulator().vanguard()
1698 .grid().comm().rank() == ioRank;
1704 .
run(this->simulator().vanguard().grid().levelGridView(0),
1705 [&vg = this->simulator().vanguard(),
1706 &emap = this->simulator().model().elementMapper()]
1708 {
return vg.gridIdxToEquilGridIdx(emap.index(elem)); });
1710 using ViolationLevel =
typename Satfunc::PhaseChecks::
1711 SatfuncConsistencyCheckManager<Scalar>::ViolationLevel;
1713 auto reportFailures = [&sfuncConsistencyChecks]
1714 (
const ViolationLevel level)
1716 sfuncConsistencyChecks.reportFailures
1717 (level, [](std::string_view record)
1718 { OpmLog::info(std::string { record }); });
1721 if (sfuncConsistencyChecks.anyFailedStandardChecks()) {
1723 OpmLog::warning(
"Saturation Function "
1724 "End-point Consistency Problems");
1726 reportFailures(ViolationLevel::Standard);
1730 if (sfuncConsistencyChecks.anyFailedCriticalChecks()) {
1732 OpmLog::error(
"Saturation Function "
1733 "End-point Consistency Failures");
1735 reportFailures(ViolationLevel::Critical);
1755 const Scalar smallSaturationTolerance_ = 1.e-6;
1757 bool enableDamarisOutput_ = false ;
1758 std::unique_ptr<DamarisWriterType> damarisWriter_;
1779 bool episodeWillBeOver()
const override
1781 const auto currTime = this->simulator().time()
1782 + this->simulator().timeStepSize();
1784 const auto nextReportStep =
1785 this->simulator().vanguard().schedule()
1786 .seconds(this->simulator().episodeIndex() + 1);
1788 const auto isSubStep = (nextReportStep - currTime)
1789 > (2 * std::numeric_limits<float>::epsilon()) * nextReportStep;
Classes required for dynamic convective mixing.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Class handling Action support in simulator.
Definition: ActionHandler.hpp:52
Provides the auxiliary methods required for consideration of the diffusion equation.
Provides the auxiliary methods required for consideration of the dispersion equation.
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
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:119
static void registerParameters()
Definition: EclWriter.hpp:144
Computes the initial condition based on the EQUIL keyword from ECL.
Definition: EquilInitializer.hpp:58
const ScalarFluidState & initialFluidState(unsigned elemIdx) const
Return the initial thermodynamic state which should be used as the initial condition.
Definition: EquilInitializer.hpp:199
BlackOilFluidState< Scalar, FluidSystem, energyModuleType !=EnergyModules::NoTemperature, energyModuleType==EnergyModules::FullyImplicitThermal, enableDissolution, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, enableSolvent, Indices::numPhases > ScalarFluidState
Definition: EquilInitializer.hpp:99
PolymerSolutionContainer< Scalar > polymer_
Definition: FlowGenericProblem.hpp:353
Scalar initialTimeStepSize_
Definition: FlowGenericProblem.hpp:365
bool enableDriftCompensation_
Definition: FlowGenericProblem.hpp:371
bool enableDriftCompensationTemp_
Definition: FlowGenericProblem.hpp:372
bool enableTuning_
Definition: FlowGenericProblem.hpp:364
void readRockParameters_(const std::vector< Scalar > &cellCenterDepths, std::function< std::array< int, 3 >(const unsigned)> ijkIndex)
Definition: FlowGenericProblem_impl.hpp:137
std::vector< Scalar > maxOilSaturation_
Definition: FlowGenericProblem.hpp:354
Scalar maxTimeStepAfterWellEvent_
Definition: FlowGenericProblem.hpp:366
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblemBlackoil.hpp:87
HybridNewton hybridNewton_
Definition: FlowProblemBlackoil.hpp:1766
void updateExplicitQuantities_(int episodeIdx, int timeStepSize, const bool first_step_after_restart) override
Definition: FlowProblemBlackoil.hpp:1200
bool updateMaxPolymerAdsorption_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblemBlackoil.hpp:1220
void handlePolymerBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1615
void writeOutput(const bool verbose) override
Write the requested quantities of the current solution into the output files.
Definition: FlowProblemBlackoil.hpp:599
void readInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1593
void handleOxygBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1631
void readEquilInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1329
void handleSolventBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1607
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:925
const std::vector< InitialFluidState > & initialFluidStates() const
Definition: FlowProblemBlackoil.hpp:789
void processRestartSaturations_(InitialFluidState &elemFluidState, Scalar &solventSaturation)
Definition: FlowProblemBlackoil.hpp:1559
std::vector< InitialFluidState > & initialFluidStates()
Definition: FlowProblemBlackoil.hpp:786
FlowProblemBlackoil(Simulator &simulator)
Definition: FlowProblemBlackoil.hpp:188
bool enableEclOutput_
Definition: FlowProblemBlackoil.hpp:1752
Scalar drsdtcon(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:1011
void endStepApplyAction()
Definition: FlowProblemBlackoil.hpp:520
bool drsdtconIsActive(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:1017
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:936
std::vector< InitialFluidState > initialFluidStates_
Definition: FlowProblemBlackoil.hpp:1750
void endTimeStep() override
Called by the simulator after each time integration.
Definition: FlowProblemBlackoil.hpp:514
void updateMaxPolymerAdsorption_()
Definition: FlowProblemBlackoil.hpp:1210
const InitialFluidState & initialFluidState(unsigned globalDofIdx) const
Definition: FlowProblemBlackoil.hpp:783
void endEpisode() override
Called by the simulator after the end of an episode.
Definition: FlowProblemBlackoil.hpp:567
void setSubStepReport(const SimulatorReportSingle &report)
Definition: FlowProblemBlackoil.hpp:795
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: FlowProblemBlackoil.hpp:967
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: FlowProblemBlackoil.hpp:323
void handleMicrBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1623
const FlowThresholdPressure< TypeTag > & thresholdPressure() const
Definition: FlowProblemBlackoil.hpp:1180
void finalizeOutput()
Definition: FlowProblemBlackoil.hpp:622
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: FlowProblemBlackoil.hpp:1028
InitialFluidState boundaryFluidState(unsigned globalDofIdx, const int directionId) const
Definition: FlowProblemBlackoil.hpp:801
std::unique_ptr< EclWriterType > eclWriter_
Definition: FlowProblemBlackoil.hpp:1753
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblemBlackoil.hpp:634
const ModuleParams & moduleParams() const
Definition: FlowProblemBlackoil.hpp:1186
void readEclRestartSolution_()
Definition: FlowProblemBlackoil.hpp:1299
FlowThresholdPressure< TypeTag > & thresholdPressure()
Definition: FlowProblemBlackoil.hpp:1183
FlowThresholdPressure< TypeTag > thresholdPressures_
Definition: FlowProblemBlackoil.hpp:1748
void readExplicitInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1344
void beginEpisode() override
Called by the simulator before an episode begins.
Definition: FlowProblemBlackoil.hpp:281
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:951
LhsEval permFactTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Calculate the transmissibility multiplier due to porosity reduction.
Definition: FlowProblemBlackoil.hpp:763
void serializeOp(Serializer &serializer)
Definition: FlowProblemBlackoil.hpp:1192
MixingRateControls< FluidSystem > mixControls_
Definition: FlowProblemBlackoil.hpp:1760
void writeReports(const SimulatorTimer &timer)
Definition: FlowProblemBlackoil.hpp:587
ModuleParams moduleParams_
Definition: FlowProblemBlackoil.hpp:1764
const EclWriterType & eclWriter() const
Definition: FlowProblemBlackoil.hpp:915
void setSimulationReport(const SimulatorReport &report)
Definition: FlowProblemBlackoil.hpp:798
void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const override
Definition: FlowProblemBlackoil.hpp:656
void computeAndSetEqWeights_()
Definition: FlowProblemBlackoil.hpp:1232
void beginTimeStep() override
Called by the simulator before each time integration.
Definition: FlowProblemBlackoil.hpp:314
void updateExplicitQuantities_(const bool first_step_after_restart)
Definition: FlowProblemBlackoil.hpp:1649
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblemBlackoil.hpp:174
ActionHandler< Scalar, IndexTraits > actionHandler_
Definition: FlowProblemBlackoil.hpp:1762
void readSolutionFromOutputModule(const int restart_step, bool fip_init)
Read simulator solution state from the outputmodule (used with restart)
Definition: FlowProblemBlackoil.hpp:1069
const EclipseIO & eclIO() const
Definition: FlowProblemBlackoil.hpp:792
Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
Definition: FlowProblemBlackoil.hpp:1177
void handleUreaBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1639
EclWriterType & eclWriter()
Definition: FlowProblemBlackoil.hpp:918
bool satfuncConsistencyRequirementsMet() const
Definition: FlowProblemBlackoil.hpp:1673
bool updateCompositionChangeLimits_()
Definition: FlowProblemBlackoil.hpp:1266
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:95
static constexpr bool enableFoam
Definition: FlowProblem.hpp:126
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:500
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:883
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:683
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:108
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:102
@ numPhases
Definition: FlowProblem.hpp:117
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:107
@ dimWorld
Definition: FlowProblem.hpp:113
@ numComponents
Definition: FlowProblem.hpp:118
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:152
@ enableExperiments
Definition: FlowProblem.hpp:133
GlobalEqVector drift_
Definition: FlowProblem.hpp:1814
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:149
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:166
int episodeIndex() const
Definition: FlowProblem.hpp:299
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:109
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:106
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:150
@ oilCompIdx
Definition: FlowProblem.hpp:145
static constexpr bool enableDiffusion
Definition: FlowProblem.hpp:123
@ numEq
Definition: FlowProblem.hpp:116
TracerModel tracerModel_
Definition: FlowProblem.hpp:1820
@ waterCompIdx
Definition: FlowProblem.hpp:146
WellModel wellModel_
Definition: FlowProblem.hpp:1816
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:307
static constexpr bool enablePolymerMolarWeight
Definition: FlowProblem.hpp:128
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:366
static constexpr bool enableSolvent
Definition: FlowProblem.hpp:129
@ gasPhaseIdx
Definition: FlowProblem.hpp:138
@ waterPhaseIdx
Definition: FlowProblem.hpp:140
static constexpr bool enablePolymer
Definition: FlowProblem.hpp:127
void readThermalParameters_()
Definition: FlowProblem.hpp:1459
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:161
@ dim
Definition: FlowProblem.hpp:112
TemperatureModel temperatureModel_
Definition: FlowProblem.hpp:1821
static constexpr bool enableExtbo
Definition: FlowProblem.hpp:125
static constexpr bool enableConvectiveMixing
Definition: FlowProblem.hpp:122
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:186
void updatePffDofData_()
Definition: FlowProblem.hpp:1607
static constexpr bool enableDispersion
Definition: FlowProblem.hpp:124
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:148
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:135
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1638
@ oilPhaseIdx
Definition: FlowProblem.hpp:139
@ enableMICP
Definition: FlowProblem.hpp:134
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1809
@ gasCompIdx
Definition: FlowProblem.hpp:144
static constexpr EnergyModules energyModuleType
Definition: FlowProblem.hpp:131
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:105
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:158
static constexpr bool enableBioeffects
Definition: FlowProblem.hpp:120
void readMaterialParameters_()
Definition: FlowProblem.hpp:1419
static constexpr bool enableBrine
Definition: FlowProblem.hpp:121
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:136
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:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Struct holding the parameters for the BlackOilBioeffectsModule class.
Definition: blackoilbioeffectsparams.hpp:42
Struct holding the parameters for the BlackoilBrineModule class.
Definition: blackoilbrineparams.hpp:42
Struct holding the parameters for the BlackoilExtboModule class.
Definition: blackoilextboparams.hpp:47
Struct holding the parameters for the BlackoilFoamModule class.
Definition: blackoilfoamparams.hpp:44
Struct holding the parameters for the BlackOilPolymerModule class.
Definition: blackoilpolymerparams.hpp:43
Struct holding the parameters for the BlackOilSolventModule class.
Definition: blackoilsolventparams.hpp:47
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34