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>
139 enum { enableDissolvedGas = Indices::compositionSwitchIdx >= 0 };
140 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
141 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
143 using SolventModule = BlackOilSolventModule<TypeTag>;
144 using PolymerModule = BlackOilPolymerModule<TypeTag>;
145 using FoamModule = BlackOilFoamModule<TypeTag>;
146 using BrineModule = BlackOilBrineModule<TypeTag>;
147 using ExtboModule = BlackOilExtboModule<TypeTag>;
148 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
149 using DispersionModule = BlackOilDispersionModule<TypeTag, enableDispersion>;
150 using DiffusionModule = BlackOilDiffusionModule<TypeTag, enableDiffusion>;
151 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
152 using ModuleParams =
typename BlackOilLocalResidualTPFA<TypeTag>::ModuleParams;
153 using HybridNewton = BlackOilHybridNewton<TypeTag>;
156 using EclWriterType = EclWriter<TypeTag, OutputBlackOilModule<TypeTag> >;
157 using IndexTraits =
typename FluidSystem::IndexTraitsType;
159 using DamarisWriterType = DamarisWriter<TypeTag>;
176 DamarisWriterType::registerParameters();
189 simulator.vanguard().schedule(),
190 simulator.vanguard().actionState(),
191 simulator.vanguard().summaryState(),
193 simulator.vanguard().grid().comm())
199 const auto& vanguard = simulator.vanguard();
202 brineParams.template initFromState<enableBrine,
203 enableSaltPrecipitation>(vanguard.eclState());
206 DiffusionModule::initFromState(vanguard.eclState());
207 DispersionModule::initFromState(vanguard.eclState());
210 extboParams.template initFromState<enableExtbo>(vanguard.eclState());
214 foamParams.template initFromState<enableFoam>(vanguard.eclState());
218 bioeffectsParams.template initFromState<enableBioeffects, enableMICP>(vanguard.eclState());
222 polymerParams.template initFromState<enablePolymer, enablePolymerMolarWeight>(vanguard.eclState());
226 solventParams.template initFromState<enableSolvent>(vanguard.eclState(), vanguard.schedule());
230 eclWriter_ = std::make_unique<EclWriterType>(simulator);
235 damarisWriter_ = std::make_unique<DamarisWriterType>(simulator);
236 enableDamarisOutput_ = Parameters::Get<Parameters::EnableDamarisOutput>();
247 auto& simulator = this->simulator();
249 const int episodeIdx = simulator.episodeIndex();
250 const auto& schedule = simulator.vanguard().schedule();
255 .evalUDQAssignments(episodeIdx, simulator.vanguard().udqState());
257 if (episodeIdx >= 0) {
258 const auto& oilVap = schedule[episodeIdx].oilvap();
259 if (oilVap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS) {
260 FluidSystem::setVapPars(oilVap.vap1(), oilVap.vap2());
263 FluidSystem::setVapPars(0.0, 0.0);
267 ConvectiveMixingModule::beginEpisode(simulator.vanguard().eclState(), schedule, episodeIdx,
268 this->moduleParams_.convectiveMixingModuleParam);
288 FlowProblemType::finishInit();
290 auto& simulator = this->simulator();
292 auto finishTransmissibilities = [updated =
false,
this]()
mutable
294 if (updated) {
return; }
296 this->
transmissibilities_.finishInit([&vg = this->simulator().vanguard()](
const unsigned int it) {
297 return vg.gridIdxToEquilGridIdx(it);
308 if (simulator.vanguard().grid().comm().size() > 1) {
309 if (simulator.vanguard().grid().comm().rank() == 0)
310 eclWriter_->setTransmissibilities(&simulator.vanguard().globalTransmissibility());
312 finishTransmissibilities();
313 eclWriter_->setTransmissibilities(&simulator.problem().eclTransmissibilities());
316 std::function<
unsigned int(
unsigned int)> equilGridToGrid = [&simulator](
unsigned int i) {
317 return simulator.vanguard().gridEquilIdxToGridIdx(i);
320 this->
eclWriter_->extractOutputTransAndNNC(equilGridToGrid);
322 simulator.vanguard().releaseGlobalTransmissibilities();
324 const auto& eclState = simulator.vanguard().eclState();
325 const auto& schedule = simulator.vanguard().schedule();
328 simulator.setStartTime(schedule.getStartTime());
329 simulator.setEndTime(schedule.simTime(schedule.size() - 1));
335 simulator.setEpisodeIndex(-1);
336 simulator.setEpisodeLength(0.0);
341 this->gravity_ = 0.0;
342 if (Parameters::Get<Parameters::EnableGravity>() &&
343 eclState.getInitConfig().hasGravity())
346 this->gravity_[dim - 1] = unit::gravity;
352 const auto& tuning = schedule[0].tuning();
357 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
358 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
363 [&simulator](
const unsigned idx)
365 std::array<int,dim> coords;
366 simulator.vanguard().cartesianCoordinate(idx, coords);
367 std::transform(coords.begin(), coords.end(), coords.begin(),
368 [](const auto c) { return c + 1; });
380 finishTransmissibilities();
382 const auto& initconfig = eclState.getInitConfig();
384 if (initconfig.restartRequested()) {
395 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>()) {
396 const auto& vanguard = this->simulator().vanguard();
397 const auto& gridView = vanguard.gridView();
398 const int numElements = gridView.size(0);
399 this->
polymer_.maxAdsorption.resize(numElements, 0.0);
408 this->
drift_.resize(this->model().numGridDof());
415 if (!initconfig.restartRequested() && !eclState.getIOConfig().initOnly()) {
416 simulator.startNextEpisode(schedule.seconds(1));
417 simulator.setEpisodeIndex(0);
418 simulator.setTimeStepIndex(0);
421 if (Parameters::Get<Parameters::CheckSatfuncConsistency>() &&
431 this->simulator().vanguard().grid().comm().barrier();
433 throw std::domain_error {
434 "Saturation function end-points do not "
435 "meet requisite consistency conditions"
444 eclState.runspec().tabdims().getNumPVTTables());
446 if (this->enableVtkOutput_() && eclState.getIOConfig().initOnly()) {
447 simulator.setTimeStepSize(0.0);
448 simulator.model().applyInitialSolution();
452 if (!eclState.getIOConfig().initOnly()) {
453 if (!this->
enableTuning_ && eclState.getSimulationConfig().anyTUNING()) {
454 OpmLog::info(
"\nThe deck has TUNING in the SCHEDULE section, but "
455 "it is ignored due\nto the flag --enable-tuning=false. "
456 "Set this flag to true to activate it.\n"
457 "Manually tuning the simulator with the TUNING keyword may "
458 "increase run time.\nIt is recommended using the simulator's "
459 "default tuning (--enable-tuning=false).");
469 FlowProblemType::endTimeStep();
470 this->endStepApplyAction();
477 this->eclWriter().mutableOutputModule().invalidateLocalData();
480 const auto& grid = this->simulator().vanguard().gridView().grid();
482 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
483 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
484 if (!isCpGrid || (grid.maxLevel() == 0)) {
485 this->eclWriter_->evalSummaryState(!this->episodeWillBeOver());
489 OPM_TIMEBLOCK(applyActions);
491 const int episodeIdx = this->episodeIndex();
492 auto& simulator = this->simulator();
496 this->simulator().vanguard().schedule().clearEvents(episodeIdx);
500 .applyActions(episodeIdx, simulator.time() + simulator.timeStepSize(),
501 [
this](
const bool global)
503 using TransUpdateQuantities = typename
504 Vanguard::TransmissibilityType::TransUpdateQuantities;
506 this->transmissibilities_
507 .update(global, TransUpdateQuantities::All,
508 [&vg = this->simulator().vanguard()]
509 (const unsigned int i)
511 return vg.gridIdxToEquilGridIdx(i);
522 OPM_TIMEBLOCK(endEpisode);
535 .evalUDQAssignments(this->episodeIndex(), this->simulator().vanguard().udqState());
537 FlowProblemType::endEpisode();
542 if (this->enableEclOutput_) {
543 this->eclWriter_->writeReports(timer);
554 FlowProblemType::writeOutput(verbose);
556 const auto isSubStep = !this->episodeWillBeOver();
558 auto localCellData = data::Solution {};
563 if (this->enableDamarisOutput_ && (this->damarisWriter_ !=
nullptr)) {
564 this->damarisWriter_->writeOutput(localCellData, isSubStep);
568 if (this->enableEclOutput_ && (this->eclWriter_ !=
nullptr)) {
569 this->eclWriter_->writeOutput(std::move(localCellData), isSubStep);
575 OPM_TIMEBLOCK(finalizeOutput);
587 FlowProblemType::initialSolutionApplied();
592 this->thresholdPressures_.finishInit();
595 const auto& grid = this->simulator().vanguard().gridView().grid();
597 using GridType = std::remove_cv_t<std::remove_reference_t<
decltype(grid)>>;
598 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
600 if (!isCpGrid || (grid.maxLevel() == 0)) {
601 if (this->simulator().episodeIndex() == 0) {
602 eclWriter_->writeInitialFIPReport();
608 unsigned globalDofIdx,
609 unsigned timeIdx)
const override
611 this->aquiferModel_.addToSource(rate, globalDofIdx, timeIdx);
614 const auto& source = this->simulator().vanguard().schedule()[this->episodeIndex()].source();
615 std::array<int,3> ijk;
616 this->simulator().vanguard().cartesianCoordinate(globalDofIdx, ijk);
618 if (source.hasSource(ijk)) {
619 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
620 static std::array<SourceComponent, 3> sc_map = {SourceComponent::WATER, SourceComponent::OIL, SourceComponent::GAS};
621 static std::array<int, 3> phidx_map = {FluidSystem::waterPhaseIdx, FluidSystem::oilPhaseIdx, FluidSystem::gasPhaseIdx};
622 static std::array<int, 3> cidx_map = {waterCompIdx, oilCompIdx, gasCompIdx};
624 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
625 const auto phaseIdx = phidx_map[i];
626 const auto sourceComp = sc_map[i];
627 const auto compIdx = cidx_map[i];
628 if (!FluidSystem::phaseIsActive(phaseIdx)) {
631 Scalar mass_rate = source.rate(ijk, sourceComp) / this->model().dofTotalVolume(globalDofIdx);
632 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
633 mass_rate /= FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
635 rate[FluidSystem::canonicalToActiveCompIdx(compIdx)] += mass_rate;
638 if constexpr (enableSolvent) {
639 Scalar mass_rate = source.rate(ijk, SourceComponent::SOLVENT) / this->model().dofTotalVolume(globalDofIdx);
640 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
641 const auto& solventPvt = SolventModule::solventPvt();
642 mass_rate /= solventPvt.referenceDensity(pvtRegionIdx);
644 rate[Indices::contiSolventEqIdx] += mass_rate;
646 if constexpr (enablePolymer) {
647 rate[Indices::polymerConcentrationIdx] += source.rate(ijk, SourceComponent::POLYMER) / this->model().dofTotalVolume(globalDofIdx);
649 if constexpr (enableMICP) {
650 rate[Indices::microbialConcentrationIdx] += source.rate(ijk, SourceComponent::MICR) / this->model().dofTotalVolume(globalDofIdx);
651 rate[Indices::oxygenConcentrationIdx] += source.rate(ijk, SourceComponent::OXYG) / this->model().dofTotalVolume(globalDofIdx);
652 rate[Indices::ureaConcentrationIdx] += source.rate(ijk, SourceComponent::UREA) / (this->model().dofTotalVolume(globalDofIdx));
654 if constexpr (enableEnergy) {
655 for (
unsigned i = 0; i < phidx_map.size(); ++i) {
656 const auto phaseIdx = phidx_map[i];
657 if (!FluidSystem::phaseIsActive(phaseIdx)) {
660 const auto sourceComp = sc_map[i];
661 const auto source_hrate = source.hrate(ijk, sourceComp);
663 rate[Indices::contiEnergyEqIdx] += source_hrate.value() / this->model().dofTotalVolume(globalDofIdx);
665 const auto& intQuants = this->simulator().model().intensiveQuantities(globalDofIdx, 0);
666 auto fs = intQuants.fluidState();
668 const auto source_temp = source.temperature(ijk, sourceComp);
670 Scalar temperature = source_temp.value();
671 fs.setTemperature(temperature);
673 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, pvtRegionIdx);
674 Scalar mass_rate = source.rate(ijk, sourceComp)/ this->model().dofTotalVolume(globalDofIdx);
675 Scalar energy_rate = getValue(h)*mass_rate;
676 rate[Indices::contiEnergyEqIdx] += energy_rate;
684 if (this->enableDriftCompensation_) {
685 const auto& simulator = this->simulator();
686 const auto& model = this->model();
691 Scalar maxCompensation = model.newtonMethod().tolerance()/10;
692 Scalar poro = this->porosity(globalDofIdx, timeIdx);
693 Scalar dt = simulator.timeStepSize();
694 EqVector dofDriftRate = this->drift_[globalDofIdx];
695 dofDriftRate /= dt*model.dofTotalVolume(globalDofIdx);
698 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
699 Scalar cnv = std::abs(dofDriftRate[eqIdx])*dt*model.eqWeight(globalDofIdx, eqIdx)/poro;
700 if (cnv > maxCompensation) {
701 dofDriftRate[eqIdx] *= maxCompensation/cnv;
705 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx)
706 rate[eqIdx] -= dofDriftRate[eqIdx];
715 template <
class LhsEval>
718 OPM_TIMEBLOCK_LOCAL(permFactTransMultiplier, Subsystem::PvtProps);
719 if constexpr (enableSaltPrecipitation) {
720 const auto& fs = intQuants.fluidState();
721 unsigned tableIdx = this->simulator().problem().satnumRegionIndex(elementIdx);
722 LhsEval porosityFactor = decay<LhsEval>(1. - fs.saltSaturation());
723 porosityFactor = min(porosityFactor, 1.0);
724 const auto& permfactTable = BrineModule::permfactTable(tableIdx);
725 return permfactTable.eval(porosityFactor,
true);
727 else if constexpr (enableBioeffects) {
728 return intQuants.permFactor().value();
737 {
return initialFluidStates_[globalDofIdx]; }
740 {
return initialFluidStates_; }
743 {
return initialFluidStates_; }
746 {
return eclWriter_->eclIO(); }
749 {
return eclWriter_->setSubStepReport(report); }
752 {
return eclWriter_->setSimulationReport(report); }
756 OPM_TIMEBLOCK_LOCAL(boundaryFluidState, Subsystem::Assembly);
757 const auto& bcprop = this->simulator().vanguard().schedule()[this->episodeIndex()].bcprop;
758 if (bcprop.size() > 0) {
759 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
763 if (this->bcindex_(dir)[globalDofIdx] == 0)
764 return initialFluidStates_[globalDofIdx];
766 const auto& bc = bcprop[this->bcindex_(dir)[globalDofIdx]];
767 if (bc.bctype == BCType::DIRICHLET )
769 InitialFluidState fluidState;
770 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
771 fluidState.setPvtRegionIndex(pvtRegionIdx);
773 switch (bc.component) {
774 case BCComponent::OIL:
775 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx))
776 throw std::logic_error(
"oil is not active and you're trying to add oil BC");
778 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0);
780 case BCComponent::GAS:
781 if (!FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
782 throw std::logic_error(
"gas is not active and you're trying to add gas BC");
784 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 1.0);
786 case BCComponent::WATER:
787 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
788 throw std::logic_error(
"water is not active and you're trying to add water BC");
790 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0);
792 case BCComponent::SOLVENT:
793 case BCComponent::POLYMER:
794 case BCComponent::MICR:
795 case BCComponent::OXYG:
796 case BCComponent::UREA:
798 throw std::logic_error(
"you need to specify a valid component (OIL, WATER or GAS) when DIRICHLET type is set in BC");
800 fluidState.setTotalSaturation(1.0);
801 double pressure = initialFluidStates_[globalDofIdx].pressure(this->refPressurePhaseIdx_());
802 const auto pressure_input = bc.pressure;
803 if (pressure_input) {
804 pressure = *pressure_input;
807 std::array<Scalar, numPhases> pc = {0};
808 const auto& matParams = this->materialLawParams(globalDofIdx);
809 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
810 Valgrind::CheckDefined(pressure);
811 Valgrind::CheckDefined(pc);
812 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
813 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
814 if (Indices::oilEnabled)
815 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
816 else if (Indices::gasEnabled)
817 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
818 else if (Indices::waterEnabled)
820 fluidState.setPressure(phaseIdx, pressure);
822 if constexpr (enableEnergy || enableTemperature) {
823 double temperature = initialFluidStates_[globalDofIdx].temperature(0);
824 const auto temperature_input = bc.temperature;
825 if(temperature_input)
826 temperature = *temperature_input;
827 fluidState.setTemperature(temperature);
830 if constexpr (enableDissolvedGas) {
831 if (FluidSystem::enableDissolvedGas()) {
832 fluidState.setRs(0.0);
833 fluidState.setRv(0.0);
836 if constexpr (enableDisgasInWater) {
837 if (FluidSystem::enableDissolvedGasInWater()) {
838 fluidState.setRsw(0.0);
841 if constexpr (enableVapwat) {
842 if (FluidSystem::enableVaporizedWater()) {
843 fluidState.setRvw(0.0);
847 for (
unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
848 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
850 const auto& b = FluidSystem::inverseFormationVolumeFactor(fluidState, phaseIdx, pvtRegionIdx);
851 fluidState.setInvB(phaseIdx, b);
853 const auto& rho = FluidSystem::density(fluidState, phaseIdx, pvtRegionIdx);
854 fluidState.setDensity(phaseIdx, rho);
855 if constexpr (enableEnergy) {
856 const auto& h = FluidSystem::enthalpy(fluidState, phaseIdx, pvtRegionIdx);
857 fluidState.setEnthalpy(phaseIdx, h);
860 fluidState.checkDefined();
864 return initialFluidStates_[globalDofIdx];
869 {
return *eclWriter_; }
872 {
return *eclWriter_; }
880 return this->mixControls_.maxGasDissolutionFactor(timeIdx, globalDofIdx,
881 this->episodeIndex(),
882 this->pvtRegionIndex(globalDofIdx));
891 return this->mixControls_.maxOilVaporizationFactor(timeIdx, globalDofIdx,
892 this->episodeIndex(),
893 this->pvtRegionIndex(globalDofIdx));
906 int episodeIdx = this->episodeIndex();
907 return !this->mixControls_.drsdtActive(episodeIdx) &&
908 !this->mixControls_.drvdtActive(episodeIdx) &&
909 this->rockCompPoroMultWc_.empty() &&
910 this->rockCompPoroMult_.empty();
919 template <
class Context>
922 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
924 values.setPvtRegionIndex(pvtRegionIndex(context, spaceIdx, timeIdx));
925 values.assignNaive(initialFluidStates_[globalDofIdx]);
927 SolventModule::assignPrimaryVars(values,
928 enableSolvent ? this->solventSaturation_[globalDofIdx] : 0.0,
929 enableSolvent ? this->solventRsw_[globalDofIdx] : 0.0);
931 if constexpr (enablePolymer)
932 values[Indices::polymerConcentrationIdx] = this->polymer_.concentration[globalDofIdx];
934 if constexpr (enablePolymerMolarWeight)
935 values[Indices::polymerMoleWeightIdx]= this->polymer_.moleWeight[globalDofIdx];
937 if constexpr (enableBrine) {
938 if (enableSaltPrecipitation && values.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
939 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltSaturation();
942 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltConcentration();
946 if constexpr (enableBioeffects) {
947 values[Indices::microbialConcentrationIdx] = this->bioeffects_.microbialConcentration[globalDofIdx];
948 values[Indices::biofilmVolumeFractionIdx]= this->bioeffects_.biofilmVolumeFraction[globalDofIdx];
949 if constexpr (enableMICP) {
950 values[Indices::oxygenConcentrationIdx]= this->bioeffects_.oxygenConcentration[globalDofIdx];
951 values[Indices::ureaConcentrationIdx]= this->bioeffects_.ureaConcentration[globalDofIdx];
952 values[Indices::calciteVolumeFractionIdx]= this->bioeffects_.calciteVolumeFraction[globalDofIdx];
956 values.checkDefined();
962 return this->mixControls_.drsdtcon(elemIdx, episodeIdx,
963 this->pvtRegionIndex(elemIdx));
968 return this->mixControls_.drsdtConvective(episodeIdx, this->pvtRegionIndex(elemIdx));
976 template <
class Context>
978 const Context& context,
980 unsigned timeIdx)
const
982 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
983 if (!context.intersection(spaceIdx).boundary())
986 if constexpr (!enableEnergy || !enableThermalFluxBoundaries)
994 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
995 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
996 values.setThermalFlow(context, spaceIdx, timeIdx, this->initialFluidStates_[globalDofIdx] );
999 if (this->nonTrivialBoundaryConditions()) {
1000 unsigned indexInInside = context.intersection(spaceIdx).indexInInside();
1001 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1002 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1003 unsigned pvtRegionIdx = pvtRegionIndex(context, spaceIdx, timeIdx);
1004 const auto [type, massrate] = this->boundaryCondition(globalDofIdx, indexInInside);
1005 if (type == BCType::THERMAL)
1006 values.setThermalFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1007 else if (type == BCType::FREE || type == BCType::DIRICHLET)
1008 values.setFreeFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1009 else if (type == BCType::RATE)
1010 values.setMassRate(massrate, pvtRegionIdx);
1020 auto& simulator = this->simulator();
1021 const auto& eclState = simulator.vanguard().eclState();
1023 std::size_t numElems = this->model().numGridDof();
1024 this->initialFluidStates_.resize(numElems);
1025 if constexpr (enableSolvent) {
1026 this->solventSaturation_.resize(numElems, 0.0);
1027 this->solventRsw_.resize(numElems, 0.0);
1030 if constexpr (enablePolymer)
1031 this->polymer_.concentration.resize(numElems, 0.0);
1033 if constexpr (enablePolymerMolarWeight) {
1034 const std::string msg {
"Support of the RESTART for polymer molecular weight "
1035 "is not implemented yet. The polymer weight value will be "
1036 "zero when RESTART begins"};
1037 OpmLog::warning(
"NO_POLYMW_RESTART", msg);
1038 this->polymer_.moleWeight.resize(numElems, 0.0);
1041 if constexpr (enableBioeffects) {
1042 this->bioeffects_.resize(numElems);
1046 this->mixControls_.init(numElems, restart_step, eclState.runspec().tabdims().getNumPVTTables());
1048 if constexpr (enableBioeffects) {
1049 this->bioeffects_ = this->eclWriter_->outputModule().getBioeffects().getSolution();
1052 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1053 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1054 elemFluidState.setPvtRegionIndex(pvtRegionIndex(elemIdx));
1055 this->eclWriter_->outputModule().initHysteresisParams(simulator, elemIdx);
1056 this->eclWriter_->outputModule().assignToFluidState(elemFluidState, elemIdx);
1065 auto ssol = enableSolvent
1066 ? this->eclWriter_->outputModule().getSolventSaturation(elemIdx)
1069 this->processRestartSaturations_(elemFluidState, ssol);
1071 if constexpr (enableSolvent) {
1072 this->solventSaturation_[elemIdx] = ssol;
1073 this->solventRsw_[elemIdx] = this->eclWriter_->outputModule().getSolventRsw(elemIdx);
1078 if constexpr (enableTemperature) {
1079 bool needTemperature = (eclState.runspec().co2Storage() || eclState.runspec().h2Storage());
1080 if (needTemperature) {
1081 const auto& fp = simulator.vanguard().eclState().fieldProps();
1082 elemFluidState.setTemperature(fp.get_double(
"TEMPI")[elemIdx]);
1086 this->mixControls_.updateLastValues(elemIdx, elemFluidState.Rs(), elemFluidState.Rv());
1088 if constexpr (enablePolymer)
1089 this->polymer_.concentration[elemIdx] = this->eclWriter_->outputModule().getPolymerConcentration(elemIdx);
1093 const int episodeIdx = this->episodeIndex();
1094 this->mixControls_.updateMaxValues(episodeIdx, simulator.timeStepSize());
1099 auto& sol = this->model().solution(0);
1100 const auto& gridView = this->gridView();
1102 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1103 elemCtx.updatePrimaryStencil(elem);
1104 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1105 this->initial(sol[elemIdx], elemCtx, 0, 0);
1113 this->model().syncOverlap();
1116 this->updateReferencePorosity_();
1117 this->mixControls_.init(this->model().numGridDof(),
1118 this->episodeIndex(),
1119 eclState.runspec().tabdims().getNumPVTTables());
1127 {
return thresholdPressures_.thresholdPressure(elem1Idx, elem2Idx); }
1130 {
return thresholdPressures_; }
1133 {
return thresholdPressures_; }
1137 return moduleParams_;
1140 template<
class Serializer>
1144 serializer(mixControls_);
1145 serializer(*eclWriter_);
1151 this->updateExplicitQuantities_(first_step_after_restart);
1153 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>())
1154 updateMaxPolymerAdsorption_();
1156 mixControls_.updateExplicitQuantities(episodeIdx, timeStepSize);
1162 this->updateProperty_(
"FlowProblemBlackoil::updateMaxPolymerAdsorption_() failed:",
1165 this->updateMaxPolymerAdsorption_(compressedDofIdx,iq);
1171 const Scalar pa = scalarValue(iq.polymerAdsorption());
1172 auto& mpa = this->polymer_.maxAdsorption;
1173 if (mpa[compressedDofIdx] < pa) {
1174 mpa[compressedDofIdx] = pa;
1183 std::vector<Scalar> sumInvB(numPhases, 0.0);
1184 const auto& gridView = this->gridView();
1186 for(
const auto& elem: elements(gridView, Dune::Partitions::interior)) {
1187 elemCtx.updatePrimaryStencil(elem);
1188 int elemIdx = elemCtx.globalSpaceIndex(0, 0);
1189 const auto& dofFluidState = this->initialFluidStates_[elemIdx];
1190 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1191 if (!FluidSystem::phaseIsActive(phaseIdx))
1194 sumInvB[phaseIdx] += dofFluidState.invB(phaseIdx);
1198 std::size_t numDof = this->model().numGridDof();
1199 const auto& comm = this->simulator().vanguard().grid().comm();
1200 comm.sum(sumInvB.data(),sumInvB.size());
1201 Scalar numTotalDof = comm.sum(numDof);
1203 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1204 if (!FluidSystem::phaseIsActive(phaseIdx))
1207 Scalar avgB = numTotalDof / sumInvB[phaseIdx];
1208 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
1209 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
1210 this->model().setEqWeight(activeSolventCompIdx, avgB);
1217 OPM_TIMEBLOCK(updateCompositionChangeLimits);
1220 int episodeIdx = this->episodeIndex();
1221 std::array<bool,3> active{this->mixControls_.drsdtConvective(episodeIdx),
1222 this->mixControls_.drsdtActive(episodeIdx),
1223 this->mixControls_.drvdtActive(episodeIdx)};
1224 if (!active[0] && !active[1] && !active[2]) {
1228 this->updateProperty_(
"FlowProblemBlackoil::updateCompositionChangeLimits_()) failed:",
1229 [
this,episodeIdx,active](
unsigned compressedDofIdx,
1232 const DimMatrix& perm = this->intrinsicPermeability(compressedDofIdx);
1233 const Scalar distZ = active[0] ? this->simulator().vanguard().cellThickness(compressedDofIdx) : 0.0;
1234 const int pvtRegionIdx = this->pvtRegionIndex(compressedDofIdx);
1235 this->mixControls_.update(compressedDofIdx,
1238 this->gravity_[
dim - 1],
1251 if(this->simulator().vanguard().grid().maxLevel() > 0) {
1252 throw std::invalid_argument(
"Refined grids are not yet supported for restart ");
1256 auto& simulator = this->simulator();
1257 const auto& schedule = simulator.vanguard().schedule();
1258 const auto& eclState = simulator.vanguard().eclState();
1259 const auto& initconfig = eclState.getInitConfig();
1260 const int restart_step = initconfig.getRestartStep();
1262 simulator.setTime(schedule.seconds(restart_step));
1264 simulator.startNextEpisode(simulator.startTime() + simulator.time(),
1265 schedule.stepLength(restart_step));
1266 simulator.setEpisodeIndex(restart_step);
1268 this->eclWriter_->beginRestart();
1270 Scalar dt = std::min(this->eclWriter_->restartTimeStepSize(), simulator.episodeLength());
1271 simulator.setTimeStepSize(dt);
1273 this->readSolutionFromOutputModule(restart_step,
false);
1275 this->eclWriter_->endRestart();
1280 const auto& simulator = this->simulator();
1285 std::size_t numElems = this->model().numGridDof();
1286 this->initialFluidStates_.resize(numElems);
1287 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1288 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1295 const auto& simulator = this->simulator();
1296 const auto& vanguard = simulator.vanguard();
1297 const auto& eclState = vanguard.eclState();
1298 const auto& fp = eclState.fieldProps();
1299 bool has_swat = fp.has_double(
"SWAT");
1300 bool has_sgas = fp.has_double(
"SGAS");
1301 bool has_rs = fp.has_double(
"RS");
1302 bool has_rsw = fp.has_double(
"RSW");
1303 bool has_rv = fp.has_double(
"RV");
1304 bool has_rvw = fp.has_double(
"RVW");
1305 bool has_pressure = fp.has_double(
"PRESSURE");
1306 bool has_salt = fp.has_double(
"SALT");
1307 bool has_saltp = fp.has_double(
"SALTP");
1310 if (Indices::numPhases > 1) {
1311 if (FluidSystem::phaseIsActive(waterPhaseIdx) && !has_swat)
1312 throw std::runtime_error(
"The ECL input file requires the presence of the SWAT keyword if "
1313 "the water phase is active");
1314 if (FluidSystem::phaseIsActive(gasPhaseIdx) && !has_sgas && FluidSystem::phaseIsActive(oilPhaseIdx))
1315 throw std::runtime_error(
"The ECL input file requires the presence of the SGAS keyword if "
1316 "the gas phase is active");
1319 throw std::runtime_error(
"The ECL input file requires the presence of the PRESSURE "
1320 "keyword if the model is initialized explicitly");
1321 if (FluidSystem::enableDissolvedGas() && !has_rs)
1322 throw std::runtime_error(
"The ECL input file requires the RS keyword to be present if"
1323 " dissolved gas is enabled and the model is initialized explicitly");
1324 if (FluidSystem::enableDissolvedGasInWater() && !has_rsw)
1325 OpmLog::warning(
"The model is initialized explicitly and the RSW keyword is not present in the"
1326 " ECL input file. The RSW values are set equal to 0");
1327 if (FluidSystem::enableVaporizedOil() && !has_rv)
1328 throw std::runtime_error(
"The ECL input file requires the RV keyword to be present if"
1329 " vaporized oil is enabled and the model is initialized explicitly");
1330 if (FluidSystem::enableVaporizedWater() && !has_rvw)
1331 throw std::runtime_error(
"The ECL input file requires the RVW keyword to be present if"
1332 " vaporized water is enabled and the model is initialized explicitly");
1333 if (enableBrine && !has_salt)
1334 throw std::runtime_error(
"The ECL input file requires the SALT keyword to be present if"
1335 " brine is enabled and the model is initialized explicitly");
1336 if (enableSaltPrecipitation && !has_saltp)
1337 throw std::runtime_error(
"The ECL input file requires the SALTP keyword to be present if"
1338 " salt precipitation is enabled and the model is initialized explicitly");
1340 std::size_t numDof = this->model().numGridDof();
1342 initialFluidStates_.resize(numDof);
1344 std::vector<double> waterSaturationData;
1345 std::vector<double> gasSaturationData;
1346 std::vector<double> pressureData;
1347 std::vector<double> rsData;
1348 std::vector<double> rswData;
1349 std::vector<double> rvData;
1350 std::vector<double> rvwData;
1351 std::vector<double> tempiData;
1352 std::vector<double> saltData;
1353 std::vector<double> saltpData;
1355 if (FluidSystem::phaseIsActive(waterPhaseIdx) && Indices::numPhases > 1)
1356 waterSaturationData = fp.get_double(
"SWAT");
1358 waterSaturationData.resize(numDof);
1360 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx))
1361 gasSaturationData = fp.get_double(
"SGAS");
1363 gasSaturationData.resize(numDof);
1365 pressureData = fp.get_double(
"PRESSURE");
1366 if (FluidSystem::enableDissolvedGas())
1367 rsData = fp.get_double(
"RS");
1369 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1370 rswData = fp.get_double(
"RSW");
1372 if (FluidSystem::enableVaporizedOil())
1373 rvData = fp.get_double(
"RV");
1375 if (FluidSystem::enableVaporizedWater())
1376 rvwData = fp.get_double(
"RVW");
1379 tempiData = fp.get_double(
"TEMPI");
1382 if constexpr (enableBrine)
1383 saltData = fp.get_double(
"SALT");
1386 if constexpr (enableSaltPrecipitation)
1387 saltpData = fp.get_double(
"SALTP");
1390 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1391 auto& dofFluidState = initialFluidStates_[dofIdx];
1393 dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
1398 if constexpr (enableEnergy || enableTemperature) {
1399 Scalar temperatureLoc = tempiData[dofIdx];
1400 if (!std::isfinite(temperatureLoc) || temperatureLoc <= 0)
1401 temperatureLoc = FluidSystem::surfaceTemperature;
1402 dofFluidState.setTemperature(temperatureLoc);
1408 if constexpr (enableBrine)
1409 dofFluidState.setSaltConcentration(saltData[dofIdx]);
1414 if constexpr (enableSaltPrecipitation)
1415 dofFluidState.setSaltSaturation(saltpData[dofIdx]);
1420 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
1421 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
1422 waterSaturationData[dofIdx]);
1424 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
1425 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
1426 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1428 - waterSaturationData[dofIdx]);
1431 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1432 gasSaturationData[dofIdx]);
1434 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1435 const Scalar soil = 1.0 - waterSaturationData[dofIdx] - gasSaturationData[dofIdx];
1436 if (soil < smallSaturationTolerance_) {
1437 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
1440 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, soil);
1447 Scalar pressure = pressureData[dofIdx];
1451 std::array<Scalar, numPhases> pc = {0};
1452 const auto& matParams = this->materialLawParams(dofIdx);
1453 MaterialLaw::capillaryPressures(pc, matParams, dofFluidState);
1454 Valgrind::CheckDefined(pressure);
1455 Valgrind::CheckDefined(pc);
1456 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1457 if (!FluidSystem::phaseIsActive(phaseIdx))
1460 if (Indices::oilEnabled)
1461 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
1462 else if (Indices::gasEnabled)
1463 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
1464 else if (Indices::waterEnabled)
1466 dofFluidState.setPressure(phaseIdx, pressure);
1469 if constexpr (enableDissolvedGas) {
1470 if (FluidSystem::enableDissolvedGas())
1471 dofFluidState.setRs(rsData[dofIdx]);
1472 else if (Indices::gasEnabled && Indices::oilEnabled)
1473 dofFluidState.setRs(0.0);
1474 if (FluidSystem::enableVaporizedOil())
1475 dofFluidState.setRv(rvData[dofIdx]);
1476 else if (Indices::gasEnabled && Indices::oilEnabled)
1477 dofFluidState.setRv(0.0);
1480 if constexpr (enableDisgasInWater) {
1481 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1482 dofFluidState.setRsw(rswData[dofIdx]);
1485 if constexpr (enableVapwat) {
1486 if (FluidSystem::enableVaporizedWater())
1487 dofFluidState.setRvw(rvwData[dofIdx]);
1493 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1494 if (!FluidSystem::phaseIsActive(phaseIdx))
1497 const auto& b = FluidSystem::inverseFormationVolumeFactor(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1498 dofFluidState.setInvB(phaseIdx, b);
1500 const auto& rho = FluidSystem::density(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1501 dofFluidState.setDensity(phaseIdx, rho);
1512 Scalar sumSaturation = 0.0;
1513 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1514 if (FluidSystem::phaseIsActive(phaseIdx)) {
1515 if (elemFluidState.saturation(phaseIdx) < smallSaturationTolerance_)
1516 elemFluidState.setSaturation(phaseIdx, 0.0);
1518 sumSaturation += elemFluidState.saturation(phaseIdx);
1522 if constexpr (enableSolvent) {
1523 if (solventSaturation < smallSaturationTolerance_)
1524 solventSaturation = 0.0;
1526 sumSaturation += solventSaturation;
1529 assert(sumSaturation > 0.0);
1531 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1532 if (FluidSystem::phaseIsActive(phaseIdx)) {
1533 const Scalar saturation = elemFluidState.saturation(phaseIdx) / sumSaturation;
1534 elemFluidState.setSaturation(phaseIdx, saturation);
1537 if constexpr (enableSolvent) {
1538 solventSaturation = solventSaturation / sumSaturation;
1544 FlowProblemType::readInitialCondition_();
1546 if constexpr (enableSolvent || enablePolymer || enablePolymerMolarWeight || enableBioeffects)
1547 this->readBlackoilExtentionsInitialConditions_(this->model().numGridDof(),
1550 enablePolymerMolarWeight,
1558 if constexpr (!enableSolvent)
1559 throw std::logic_error(
"solvent is disabled and you're trying to add solvent to BC");
1561 rate[Indices::solventSaturationIdx] = bc.rate;
1566 if constexpr (!enablePolymer)
1567 throw std::logic_error(
"polymer is disabled and you're trying to add polymer to BC");
1569 rate[Indices::polymerConcentrationIdx] = bc.rate;
1574 if constexpr (!enableMICP)
1575 throw std::logic_error(
"MICP is disabled and you're trying to add microbes to BC");
1577 rate[Indices::microbialConcentrationIdx] = bc.rate;
1582 if constexpr (!enableMICP)
1583 throw std::logic_error(
"MICP is disabled and you're trying to add oxygen to BC");
1585 rate[Indices::oxygenConcentrationIdx] = bc.rate;
1590 if constexpr (!enableMICP)
1591 throw std::logic_error(
"MICP is disabled and you're trying to add urea to BC");
1593 rate[Indices::ureaConcentrationIdx] = bc.rate;
1595 rate[Indices::ureaConcentrationIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
1600 OPM_TIMEBLOCK(updateExplicitQuantities);
1601 const bool invalidateFromMaxWaterSat = this->updateMaxWaterSaturation_();
1602 const bool invalidateFromMinPressure = this->updateMinPressure_();
1605 const bool invalidateFromHyst = this->updateHysteresis_();
1606 const bool invalidateFromMaxOilSat = this->updateMaxOilSaturation_();
1609 const bool invalidateDRDT = !first_step_after_restart && this->updateCompositionChangeLimits_();
1612 const bool invalidateIntensiveQuantities
1613 = invalidateFromMaxWaterSat || invalidateFromMinPressure || invalidateFromHyst || invalidateFromMaxOilSat || invalidateDRDT;
1614 if (invalidateIntensiveQuantities) {
1615 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
1616 this->model().invalidateAndUpdateIntensiveQuantities(0);
1619 this->updateRockCompTransMultVal_();
1624 if (
const auto nph = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
1625 + FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
1626 + FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1635 const auto numSamplePoints =
static_cast<std::size_t
>
1636 (Parameters::Get<Parameters::NumSatfuncConsistencySamplePoints>());
1638 auto sfuncConsistencyChecks =
1640 numSamplePoints, this->simulator().vanguard().eclState(),
1641 [&cmap = this->simulator().vanguard().cartesianIndexMapper()](
const int elemIdx)
1642 {
return cmap.cartesianIndex(elemIdx); }
1645 const auto ioRank = 0;
1646 const auto isIoRank = this->simulator().vanguard()
1647 .grid().comm().rank() == ioRank;
1653 .
run(this->simulator().vanguard().grid().levelGridView(0),
1654 [&vg = this->simulator().vanguard(),
1655 &emap = this->simulator().model().elementMapper()]
1657 {
return vg.gridIdxToEquilGridIdx(emap.index(elem)); });
1659 using ViolationLevel =
typename Satfunc::PhaseChecks::
1660 SatfuncConsistencyCheckManager<Scalar>::ViolationLevel;
1662 auto reportFailures = [&sfuncConsistencyChecks]
1663 (
const ViolationLevel level)
1665 sfuncConsistencyChecks.reportFailures
1666 (level, [](std::string_view record)
1667 { OpmLog::info(std::string { record }); });
1670 if (sfuncConsistencyChecks.anyFailedStandardChecks()) {
1672 OpmLog::warning(
"Saturation Function "
1673 "End-point Consistency Problems");
1675 reportFailures(ViolationLevel::Standard);
1679 if (sfuncConsistencyChecks.anyFailedCriticalChecks()) {
1681 OpmLog::error(
"Saturation Function "
1682 "End-point Consistency Failures");
1684 reportFailures(ViolationLevel::Critical);
1704 const Scalar smallSaturationTolerance_ = 1.e-6;
1706 bool enableDamarisOutput_ = false ;
1707 std::unique_ptr<DamarisWriterType> damarisWriter_;
1728 bool episodeWillBeOver()
const override
1730 const auto currTime = this->simulator().time()
1731 + this->simulator().timeStepSize();
1733 const auto nextReportStep =
1734 this->simulator().vanguard().schedule()
1735 .seconds(this->simulator().episodeIndex() + 1);
1737 const auto isSubStep = (nextReportStep - currTime)
1738 > (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, enableTemperature, enableEnergy, enableDissolution, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, Indices::numPhases > ScalarFluidState
Definition: EquilInitializer.hpp:99
const ScalarFluidState & initialFluidState(unsigned elemIdx) const
Return the initial thermodynamic state which should be used as the initial condition.
Definition: EquilInitializer.hpp:199
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:1715
void updateExplicitQuantities_(int episodeIdx, int timeStepSize, const bool first_step_after_restart) override
Definition: FlowProblemBlackoil.hpp:1149
bool updateMaxPolymerAdsorption_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblemBlackoil.hpp:1169
void handlePolymerBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1564
void writeOutput(const bool verbose) override
Write the requested quantities of the current solution into the output files.
Definition: FlowProblemBlackoil.hpp:552
void readInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1542
void handleOxygBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1580
void readEquilInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1278
void handleSolventBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1556
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:878
const std::vector< InitialFluidState > & initialFluidStates() const
Definition: FlowProblemBlackoil.hpp:742
void processRestartSaturations_(InitialFluidState &elemFluidState, Scalar &solventSaturation)
Definition: FlowProblemBlackoil.hpp:1508
std::vector< InitialFluidState > & initialFluidStates()
Definition: FlowProblemBlackoil.hpp:739
FlowProblemBlackoil(Simulator &simulator)
Definition: FlowProblemBlackoil.hpp:184
bool enableEclOutput_
Definition: FlowProblemBlackoil.hpp:1701
Scalar drsdtcon(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:960
void endStepApplyAction()
Definition: FlowProblemBlackoil.hpp:473
bool drsdtconIsActive(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:966
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:889
std::vector< InitialFluidState > initialFluidStates_
Definition: FlowProblemBlackoil.hpp:1699
void endTimeStep() override
Called by the simulator after each time integration.
Definition: FlowProblemBlackoil.hpp:467
void updateMaxPolymerAdsorption_()
Definition: FlowProblemBlackoil.hpp:1159
const InitialFluidState & initialFluidState(unsigned globalDofIdx) const
Definition: FlowProblemBlackoil.hpp:736
void endEpisode() override
Called by the simulator after the end of an episode.
Definition: FlowProblemBlackoil.hpp:520
void setSubStepReport(const SimulatorReportSingle &report)
Definition: FlowProblemBlackoil.hpp:748
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: FlowProblemBlackoil.hpp:920
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: FlowProblemBlackoil.hpp:283
void handleMicrBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1572
const FlowThresholdPressure< TypeTag > & thresholdPressure() const
Definition: FlowProblemBlackoil.hpp:1129
void finalizeOutput()
Definition: FlowProblemBlackoil.hpp:573
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: FlowProblemBlackoil.hpp:977
InitialFluidState boundaryFluidState(unsigned globalDofIdx, const int directionId) const
Definition: FlowProblemBlackoil.hpp:754
std::unique_ptr< EclWriterType > eclWriter_
Definition: FlowProblemBlackoil.hpp:1702
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblemBlackoil.hpp:585
const ModuleParams & moduleParams() const
Definition: FlowProblemBlackoil.hpp:1135
void readEclRestartSolution_()
Definition: FlowProblemBlackoil.hpp:1248
FlowThresholdPressure< TypeTag > & thresholdPressure()
Definition: FlowProblemBlackoil.hpp:1132
FlowThresholdPressure< TypeTag > thresholdPressures_
Definition: FlowProblemBlackoil.hpp:1697
void readExplicitInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1293
void beginEpisode() override
Called by the simulator before an episode begins.
Definition: FlowProblemBlackoil.hpp:243
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:904
LhsEval permFactTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to porosity reduction.
Definition: FlowProblemBlackoil.hpp:716
void serializeOp(Serializer &serializer)
Definition: FlowProblemBlackoil.hpp:1141
MixingRateControls< FluidSystem > mixControls_
Definition: FlowProblemBlackoil.hpp:1709
void writeReports(const SimulatorTimer &timer)
Definition: FlowProblemBlackoil.hpp:540
ModuleParams moduleParams_
Definition: FlowProblemBlackoil.hpp:1713
const EclWriterType & eclWriter() const
Definition: FlowProblemBlackoil.hpp:868
void setSimulationReport(const SimulatorReport &report)
Definition: FlowProblemBlackoil.hpp:751
void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const override
Definition: FlowProblemBlackoil.hpp:607
void computeAndSetEqWeights_()
Definition: FlowProblemBlackoil.hpp:1181
void beginTimeStep() override
Called by the simulator before each time integration.
Definition: FlowProblemBlackoil.hpp:274
void updateExplicitQuantities_(const bool first_step_after_restart)
Definition: FlowProblemBlackoil.hpp:1598
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblemBlackoil.hpp:170
ActionHandler< Scalar, IndexTraits > actionHandler_
Definition: FlowProblemBlackoil.hpp:1711
void readSolutionFromOutputModule(const int restart_step, bool fip_init)
Read simulator solution state from the outputmodule (used with restart)
Definition: FlowProblemBlackoil.hpp:1018
const EclipseIO & eclIO() const
Definition: FlowProblemBlackoil.hpp:745
Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
Definition: FlowProblemBlackoil.hpp:1126
void handleUreaBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1588
EclWriterType & eclWriter()
Definition: FlowProblemBlackoil.hpp:871
bool satfuncConsistencyRequirementsMet() const
Definition: FlowProblemBlackoil.hpp:1622
bool updateCompositionChangeLimits_()
Definition: FlowProblemBlackoil.hpp:1215
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:94
@ enableDiffusion
Definition: FlowProblem.hpp:122
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:479
@ gasPhaseIdx
Definition: FlowProblem.hpp:136
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:132
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:826
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:659
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:107
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:131
@ enableTemperature
Definition: FlowProblem.hpp:125
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:101
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:106
@ enableBrine
Definition: FlowProblem.hpp:120
@ enableExtbo
Definition: FlowProblem.hpp:127
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:150
GlobalEqVector drift_
Definition: FlowProblem.hpp:1680
@ numEq
Definition: FlowProblem.hpp:115
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:147
@ waterPhaseIdx
Definition: FlowProblem.hpp:138
@ numComponents
Definition: FlowProblem.hpp:117
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:164
int episodeIndex() const
Definition: FlowProblem.hpp:287
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:108
@ oilPhaseIdx
Definition: FlowProblem.hpp:137
@ gasCompIdx
Definition: FlowProblem.hpp:142
@ enableFoam
Definition: FlowProblem.hpp:128
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:105
@ enableBioeffects
Definition: FlowProblem.hpp:119
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:148
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:134
@ enableSolvent
Definition: FlowProblem.hpp:133
@ enableDispersion
Definition: FlowProblem.hpp:123
TracerModel tracerModel_
Definition: FlowProblem.hpp:1686
WellModel wellModel_
Definition: FlowProblem.hpp:1682
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:295
@ enableConvectiveMixing
Definition: FlowProblem.hpp:121
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:352
@ dimWorld
Definition: FlowProblem.hpp:112
@ numPhases
Definition: FlowProblem.hpp:116
void readThermalParameters_()
Definition: FlowProblem.hpp:1375
@ enablePolymer
Definition: FlowProblem.hpp:130
@ enableEnergy
Definition: FlowProblem.hpp:124
@ waterCompIdx
Definition: FlowProblem.hpp:144
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:159
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:102
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:181
void updatePffDofData_()
Definition: FlowProblem.hpp:1488
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:146
@ enableMICP
Definition: FlowProblem.hpp:129
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1519
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1675
@ oilCompIdx
Definition: FlowProblem.hpp:143
@ dim
Definition: FlowProblem.hpp:111
@ enableExperiments
Definition: FlowProblem.hpp:126
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:104
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:156
void readMaterialParameters_()
Definition: FlowProblem.hpp:1341
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