28#ifndef EWOMS_BLACK_OIL_SOLVENT_MODULE_HH
29#define EWOMS_BLACK_OIL_SOLVENT_MODULE_HH
33#include <opm/common/Exceptions.hpp>
40#include <opm/material/fluidsystems/blackoilpvt/SolventPvt.hpp>
43#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
44#include <opm/input/eclipse/EclipseState/Tables/SsfnTable.hpp>
45#include <opm/input/eclipse/EclipseState/Tables/Sof2Table.hpp>
46#include <opm/input/eclipse/EclipseState/Tables/MsfnTable.hpp>
47#include <opm/input/eclipse/EclipseState/Tables/PmiscTable.hpp>
48#include <opm/input/eclipse/EclipseState/Tables/MiscTable.hpp>
49#include <opm/input/eclipse/EclipseState/Tables/SorwmisTable.hpp>
50#include <opm/input/eclipse/EclipseState/Tables/SgcwmisTable.hpp>
51#include <opm/input/eclipse/EclipseState/Tables/TlpmixpaTable.hpp>
54#include <opm/material/common/Valgrind.hpp>
56#include <dune/common/fvector.hh>
67template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
82 using Toolbox = MathToolbox<Evaluation>;
91 static constexpr unsigned solventSaturationIdx = Indices::solventSaturationIdx;
92 static constexpr unsigned contiSolventEqIdx = Indices::contiSolventEqIdx;
93 static constexpr unsigned enableSolvent = enableSolventV;
94 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
95 static constexpr unsigned numPhases = FluidSystem::numPhases;
96 static constexpr bool blackoilConserveSurfaceVolume = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
97 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
105 static void initFromState(
const EclipseState& eclState,
const Schedule& schedule)
109 if (enableSolvent && !eclState.runspec().phases().active(Phase::SOLVENT))
110 throw std::runtime_error(
"Non-trivial solvent treatment requested at compile "
111 "time, but the deck does not contain the SOLVENT keyword");
112 else if (!enableSolvent && eclState.runspec().phases().active(Phase::SOLVENT))
113 throw std::runtime_error(
"Solvent treatment disabled at compile time, but the deck "
114 "contains the SOLVENT keyword");
116 if (!eclState.runspec().phases().active(Phase::SOLVENT))
119 params_.co2sol_ = eclState.runspec().co2Sol();
120 params_.h2sol_ = eclState.runspec().h2Sol();
123 throw std::runtime_error(
"CO2SOL and H2SOL can not be used together");
128 params_.co2GasPvt_.initFromState(eclState, schedule);
129 params_.brineCo2Pvt_.initFromState(eclState, schedule);
131 params_.h2GasPvt_.initFromState(eclState, schedule);
132 params_.brineH2Pvt_.initFromState(eclState, schedule);
134 if (eclState.getSimulationConfig().hasDISGASW()) {
135 params_.rsSolw_active_ =
true;
138 params_.solventPvt_.initFromState(eclState, schedule);
140 const auto& tableManager = eclState.getTableManager();
142 const auto& ssfnTables = tableManager.getSsfnTables();
143 unsigned numSatRegions = tableManager.getTabdims().getNumSatTables();
144 params_.setNumSatRegions(numSatRegions);
145 for (
unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
146 const auto& ssfnTable = ssfnTables.template getTable<SsfnTable>(satRegionIdx);
147 params_.ssfnKrg_[satRegionIdx].setXYContainers(ssfnTable.getSolventFractionColumn(),
148 ssfnTable.getGasRelPermMultiplierColumn(),
150 params_.ssfnKrs_[satRegionIdx].setXYContainers(ssfnTable.getSolventFractionColumn(),
151 ssfnTable.getSolventRelPermMultiplierColumn(),
156 params_.isMiscible_ =
false;
157 if (!eclState.getTableManager().getMiscTables().empty()) {
158 params_.isMiscible_ =
true;
160 unsigned numMiscRegions = 1;
163 const auto& sof2Tables = tableManager.getSof2Tables();
164 if (!sof2Tables.empty()) {
166 params_.sof2Krn_.resize(numSatRegions);
167 for (
unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
168 const auto& sof2Table = sof2Tables.template getTable<Sof2Table>(satRegionIdx);
169 params_.sof2Krn_[satRegionIdx].setXYContainers(sof2Table.getSoColumn(),
170 sof2Table.getKroColumn(),
174 else if(eclState.runspec().phases().active(Phase::OIL))
175 throw std::runtime_error(
"SOF2 must be specified in MISCIBLE (SOLVENT and OIL) runs\n");
177 const auto& miscTables = tableManager.getMiscTables();
178 if (!miscTables.empty()) {
179 assert(numMiscRegions == miscTables.size());
182 params_.misc_.resize(numMiscRegions);
183 for (
unsigned miscRegionIdx = 0; miscRegionIdx < numMiscRegions; ++miscRegionIdx) {
184 const auto& miscTable = miscTables.template getTable<MiscTable>(miscRegionIdx);
187 const auto& solventFraction = miscTable.getSolventFractionColumn();
188 const auto&
misc = miscTable.getMiscibilityColumn();
189 params_.misc_[miscRegionIdx].setXYContainers(solventFraction,
misc);
193 throw std::runtime_error(
"MISC must be specified in MISCIBLE (SOLVENT) runs\n");
196 params_.pmisc_.resize(numMiscRegions);
197 const auto& pmiscTables = tableManager.getPmiscTables();
198 if (!pmiscTables.empty()) {
199 assert(numMiscRegions == pmiscTables.size());
201 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
202 const auto& pmiscTable = pmiscTables.template getTable<PmiscTable>(regionIdx);
205 const auto& po = pmiscTable.getOilPhasePressureColumn();
206 const auto&
pmisc = pmiscTable.getMiscibilityColumn();
208 params_.pmisc_[regionIdx].setXYContainers(po,
pmisc);
212 std::vector<double> x = {0.0,1.0e20};
213 std::vector<double> y = {1.0,1.0};
214 TabulatedFunction constant = TabulatedFunction(2, x, y);
215 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
216 params_.pmisc_[regionIdx] = constant;
221 params_.msfnKrsg_.resize(numSatRegions);
222 params_.msfnKro_.resize(numSatRegions);
223 const auto& msfnTables = tableManager.getMsfnTables();
224 if (!msfnTables.empty()) {
225 assert(numSatRegions == msfnTables.size());
227 for (
unsigned regionIdx = 0; regionIdx < numSatRegions; ++regionIdx) {
228 const MsfnTable& msfnTable = msfnTables.template getTable<MsfnTable>(regionIdx);
232 const auto& Ssg = msfnTable.getGasPhaseFractionColumn();
233 const auto& krsg = msfnTable.getGasSolventRelpermMultiplierColumn();
234 const auto& kro = msfnTable.getOilRelpermMultiplierColumn();
236 params_.msfnKrsg_[regionIdx].setXYContainers(Ssg, krsg);
237 params_.msfnKro_[regionIdx].setXYContainers(Ssg, kro);
241 std::vector<double> x = {0.0,1.0};
242 std::vector<double> y = {1.0,0.0};
243 TabulatedFunction unit = TabulatedFunction(2, x, x);
244 TabulatedFunction invUnit = TabulatedFunction(2, x, y);
246 for (
unsigned regionIdx = 0; regionIdx < numSatRegions; ++regionIdx) {
247 params_.setMsfn(regionIdx, unit, invUnit);
251 params_.sorwmis_.resize(numMiscRegions);
252 const auto& sorwmisTables = tableManager.getSorwmisTables();
253 if (!sorwmisTables.empty()) {
254 assert(numMiscRegions == sorwmisTables.size());
256 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
257 const auto& sorwmisTable = sorwmisTables.template getTable<SorwmisTable>(regionIdx);
260 const auto& sw = sorwmisTable.getWaterSaturationColumn();
261 const auto&
sorwmis = sorwmisTable.getMiscibleResidualOilColumn();
263 params_.sorwmis_[regionIdx].setXYContainers(sw,
sorwmis);
268 std::vector<double> x = {0.0,1.0};
269 std::vector<double> y = {0.0,0.0};
270 TabulatedFunction zero = TabulatedFunction(2, x, y);
271 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
272 params_.sorwmis_[regionIdx] = zero;
277 params_.sgcwmis_.resize(numMiscRegions);
278 const auto& sgcwmisTables = tableManager.getSgcwmisTables();
279 if (!sgcwmisTables.empty()) {
280 assert(numMiscRegions == sgcwmisTables.size());
282 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
283 const auto& sgcwmisTable = sgcwmisTables.template getTable<SgcwmisTable>(regionIdx);
286 const auto& sw = sgcwmisTable.getWaterSaturationColumn();
287 const auto&
sgcwmis = sgcwmisTable.getMiscibleResidualGasColumn();
289 params_.sgcwmis_[regionIdx].setXYContainers(sw,
sgcwmis);
294 std::vector<double> x = {0.0,1.0};
295 std::vector<double> y = {0.0,0.0};
296 TabulatedFunction zero = TabulatedFunction(2, x, y);
297 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx)
298 params_.sgcwmis_[regionIdx] = zero;
301 const auto& tlmixpar = eclState.getTableManager().getTLMixpar();
302 if (!tlmixpar.empty()) {
304 params_.tlMixParamViscosity_.resize(numMiscRegions);
305 params_.tlMixParamDensity_.resize(numMiscRegions);
307 assert(numMiscRegions == tlmixpar.size());
308 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
309 const auto& tlp = tlmixpar[regionIdx];
310 params_.tlMixParamViscosity_[regionIdx] = tlp.viscosity_parameter;
311 params_.tlMixParamDensity_[regionIdx] = tlp.density_parameter;
315 throw std::runtime_error(
"TLMIXPAR must be specified in MISCIBLE (SOLVENT) runs\n");
318 params_.tlPMixTable_.resize(numMiscRegions);
319 if (!eclState.getTableManager().getTlpmixpaTables().empty()) {
320 const auto& tlpmixparTables = tableManager.getTlpmixpaTables();
321 if (!tlpmixparTables.empty()) {
322 assert(numMiscRegions == tlpmixparTables.size());
323 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
324 const auto& tlpmixparTable = tlpmixparTables.template getTable<TlpmixpaTable>(regionIdx);
327 const auto& po = tlpmixparTable.getOilPhasePressureColumn();
328 const auto& tlpmixpa = tlpmixparTable.getMiscibilityColumn();
330 params_.tlPMixTable_[regionIdx].setXYContainers(po, tlpmixpa);
335 if (params_.pmisc_.size() > 0)
336 params_.tlPMixTable_ = params_.pmisc_;
338 throw std::invalid_argument(
"If the pressure dependent TL values in "
339 "TLPMIXPA is defaulted (no entries), then "
340 "the PMISC tables must be specified.");
345 std::vector<double> x = {0.0,1.0e20};
346 std::vector<double> y = {1.0,1.0};
347 TabulatedFunction ones = TabulatedFunction(2, x, y);
348 for (
unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx)
349 params_.tlPMixTable_[regionIdx] = ones;
359 { params_.solventPvt_ = value; }
370 if constexpr (enableSolvent)
378 Simulator& simulator)
380 if constexpr (enableSolvent)
386 if constexpr (enableSolvent)
387 return pvIdx == solventSaturationIdx;
396 return "saturation_solvent";
404 return static_cast<Scalar
>(1.0);
409 if constexpr (enableSolvent)
410 return eqIdx == contiSolventEqIdx;
415 static std::string
eqName([[maybe_unused]]
unsigned eqIdx)
419 return "conti^solvent";
422 static Scalar
eqWeight([[maybe_unused]]
unsigned eqIdx)
427 return static_cast<Scalar
>(1.0);
430 template <
class LhsEval>
431 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
432 const IntensiveQuantities& intQuants)
434 if constexpr (enableSolvent) {
435 if constexpr (blackoilConserveSurfaceVolume) {
436 storage[contiSolventEqIdx] +=
437 Toolbox::template decay<LhsEval>(intQuants.porosity())
438 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
439 * Toolbox::template decay<LhsEval>(intQuants.solventInverseFormationVolumeFactor());
441 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
442 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
443 * Toolbox::template decay<LhsEval>(intQuants.fluidState().invB(waterPhaseIdx))
444 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
448 storage[contiSolventEqIdx] +=
449 Toolbox::template decay<LhsEval>(intQuants.porosity())
450 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
451 * Toolbox::template decay<LhsEval>(intQuants.solventDensity());
453 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
454 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
455 * Toolbox::template decay<LhsEval>(intQuants.fluidState().density(waterPhaseIdx))
456 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
464 [[maybe_unused]]
const ElementContext& elemCtx,
465 [[maybe_unused]]
unsigned scvfIdx,
466 [[maybe_unused]]
unsigned timeIdx)
469 if constexpr (enableSolvent) {
470 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
472 unsigned upIdx = extQuants.solventUpstreamIndex();
473 unsigned inIdx = extQuants.interiorIndex();
474 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
476 if constexpr (blackoilConserveSurfaceVolume) {
478 flux[contiSolventEqIdx] =
479 extQuants.solventVolumeFlux()
480 *up.solventInverseFormationVolumeFactor();
482 flux[contiSolventEqIdx] =
483 extQuants.solventVolumeFlux()
484 *decay<Scalar>(up.solventInverseFormationVolumeFactor());
489 flux[contiSolventEqIdx] +=
490 extQuants.volumeFlux(waterPhaseIdx)
491 * up.fluidState().invB(waterPhaseIdx)
494 flux[contiSolventEqIdx] +=
495 extQuants.volumeFlux(waterPhaseIdx)
496 *decay<Scalar>(up.fluidState().invB(waterPhaseIdx))
497 *decay<Scalar>(up.rsSolw());
502 flux[contiSolventEqIdx] =
503 extQuants.solventVolumeFlux()
504 *up.solventDensity();
506 flux[contiSolventEqIdx] =
507 extQuants.solventVolumeFlux()
508 *decay<Scalar>(up.solventDensity());
513 flux[contiSolventEqIdx] +=
514 extQuants.volumeFlux(waterPhaseIdx)
515 * up.fluidState().density(waterPhaseIdx)
518 flux[contiSolventEqIdx] +=
519 extQuants.volumeFlux(waterPhaseIdx)
520 *decay<Scalar>(up.fluidState().density(waterPhaseIdx))
521 *decay<Scalar>(up.rsSolw());
531 Scalar solventSaturation,
534 if constexpr (!enableSolvent) {
535 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Disabled);
540 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Ss);
541 priVars[solventSaturationIdx] = solventSaturation;
543 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Rsolw);
544 priVars[solventSaturationIdx] = solventRsw;
552 const PrimaryVariables& oldPv,
553 const EqVector& delta)
555 if constexpr (enableSolvent)
557 newPv[solventSaturationIdx] = oldPv[solventSaturationIdx] - delta[solventSaturationIdx];
569 return static_cast<Scalar
>(0.0);
578 return std::abs(Toolbox::scalarValue(resid[contiSolventEqIdx]));
581 template <
class DofEntity>
582 static void serializeEntity(
const Model& model, std::ostream& outstream,
const DofEntity& dof)
584 if constexpr (enableSolvent) {
585 unsigned dofIdx = model.dofMapper().index(dof);
587 const PrimaryVariables& priVars = model.solution(0)[dofIdx];
588 outstream << priVars[solventSaturationIdx];
592 template <
class DofEntity>
595 if constexpr (enableSolvent) {
596 unsigned dofIdx = model.dofMapper().index(dof);
598 PrimaryVariables& priVars0 = model.solution(0)[dofIdx];
599 PrimaryVariables& priVars1 = model.solution(1)[dofIdx];
601 instream >> priVars0[solventSaturationIdx];
604 priVars1 = priVars0[solventSaturationIdx];
610 return params_.solventPvt_;
616 return params_.co2GasPvt_;
621 return params_.h2GasPvt_;
626 return params_.brineCo2Pvt_;
631 return params_.brineH2Pvt_;
634 static const TabulatedFunction&
ssfnKrg(
const ElementContext& elemCtx,
638 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
639 return params_.ssfnKrg_[satnumRegionIdx];
642 static const TabulatedFunction&
ssfnKrs(
const ElementContext& elemCtx,
646 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
647 return params_.ssfnKrs_[satnumRegionIdx];
650 static const TabulatedFunction&
sof2Krn(
const ElementContext& elemCtx,
654 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
655 return params_.sof2Krn_[satnumRegionIdx];
658 static const TabulatedFunction&
misc(
const ElementContext& elemCtx,
662 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
663 return params_.misc_[miscnumRegionIdx];
666 static const TabulatedFunction&
pmisc(
const ElementContext& elemCtx,
670 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
671 return params_.pmisc_[miscnumRegionIdx];
674 static const TabulatedFunction&
msfnKrsg(
const ElementContext& elemCtx,
678 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
679 return params_.msfnKrsg_[satnumRegionIdx];
682 static const TabulatedFunction&
msfnKro(
const ElementContext& elemCtx,
686 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
687 return params_.msfnKro_[satnumRegionIdx];
690 static const TabulatedFunction&
sorwmis(
const ElementContext& elemCtx,
694 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
695 return params_.sorwmis_[miscnumRegionIdx];
698 static const TabulatedFunction&
sgcwmis(
const ElementContext& elemCtx,
702 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
703 return params_.sgcwmis_[miscnumRegionIdx];
706 static const TabulatedFunction&
tlPMixTable(
const ElementContext& elemCtx,
710 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
711 return params_.tlPMixTable_[miscnumRegionIdx];
718 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
719 return params_.tlMixParamViscosity_[miscnumRegionIdx];
726 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
727 return params_.tlMixParamDensity_[miscnumRegionIdx];
732 return params_.isMiscible_;
735 template <
class Value>
736 static const Value
solubilityLimit(
unsigned pvtIdx,
const Value& temperature,
const Value& pressure,
const Value& saltConcentration)
743 return brineCo2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
745 return brineH2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
750 return params_.rsSolw_active_;
755 return params_.co2sol_;
760 return params_.h2sol_;
767template <
class TypeTag,
bool enableSolventV>
768BlackOilSolventParams<typename BlackOilSolventModule<TypeTag, enableSolventV>::Scalar>
769BlackOilSolventModule<TypeTag, enableSolventV>::params_;
778template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
793 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
794 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
795 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
796 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
797 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
798 static constexpr double cutOff = 1e-12;
812 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
814 auto& fs =
asImp_().fluidState_;
816 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
817 solventSaturation_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
844 auto& fs =
asImp_().fluidState_;
850 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
851 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
853 }
else if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Rsolw) {
854 rsSolw_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
865 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
867 const Evaluation& pgImisc = fs.pressure(gasPhaseIdx);
870 const auto& problem = elemCtx.problem();
871 Evaluation pgMisc = 0.0;
872 std::array<Evaluation, numPhases> pC;
873 const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
874 MaterialLaw::capillaryPressures(pC, materialParams, fs);
877 const auto linearizationType = elemCtx.linearizationType();
878 if (priVars.primaryVarsMeaningPressure() == PrimaryVariables::PressureMeaning::Pg)
879 pgMisc = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
881 const Evaluation& po = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
882 pgMisc = po + (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
885 fs.setPressure(gasPhaseIdx, pmisc * pgMisc + (1.0 - pmisc) * pgImisc);
891 if (gasSolventSat.value() < cutOff)
901 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
902 const Evaluation miscibility = misc.eval(Fsolgas,
true) * pmisc.eval(p,
true);
905 unsigned cellIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
906 const auto& materialLawManager = elemCtx.problem().materialLawManager();
907 const auto& scaledDrainageInfo =
908 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
910 const Scalar& sogcr = scaledDrainageInfo.Sogcr;
911 Evaluation sor = sogcr;
912 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
913 const Evaluation& sw = fs.saturation(waterPhaseIdx);
915 sor = miscibility * sorwmis.eval(sw,
true) + (1.0 - miscibility) * sogcr;
917 const Scalar& sgcr = scaledDrainageInfo.Sgcr;
918 Evaluation sgc = sgcr;
919 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
920 const Evaluation& sw = fs.saturation(waterPhaseIdx);
922 sgc = miscibility * sgcwmis.eval(sw,
true) + (1.0 - miscibility) * sgcr;
925 Evaluation oilGasSolventSat = gasSolventSat;
926 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
927 oilGasSolventSat += fs.saturation(oilPhaseIdx);
929 const Evaluation zero = 0.0;
930 const Evaluation oilGasSolventEffSat = std::max(oilGasSolventSat - sor - sgc, zero);
932 Evaluation F_totalGas = 0.0;
933 if (oilGasSolventEffSat.value() > cutOff) {
934 const Evaluation gasSolventEffSat = std::max(gasSolventSat - sgc, zero);
935 F_totalGas = gasSolventEffSat / oilGasSolventEffSat;
941 const Evaluation mkrgt = msfnKrsg.eval(F_totalGas,
true) * sof2Krn.eval(oilGasSolventSat,
true);
942 const Evaluation mkro = msfnKro.eval(F_totalGas,
true) * sof2Krn.eval(oilGasSolventSat,
true);
944 Evaluation& kro =
asImp_().mobility_[oilPhaseIdx];
945 Evaluation& krg =
asImp_().mobility_[gasPhaseIdx];
948 krg *= (1.0 - miscibility);
949 krg += miscibility * mkrgt;
950 kro *= (1.0 - miscibility);
951 kro += miscibility * mkro;
958 Evaluation& krg =
asImp_().mobility_[gasPhaseIdx];
960 krg *= ssfnKrg.eval(Fhydgas,
true);
974 const auto& iq =
asImp_();
975 unsigned pvtRegionIdx = iq.pvtRegionIndex();
976 const Evaluation& T = iq.fluidState().temperature(gasPhaseIdx);
977 const Evaluation& p = iq.fluidState().pressure(gasPhaseIdx);
979 const Evaluation rv = 0.0;
980 const Evaluation rvw = 0.0;
989 auto& fs =
asImp_().fluidState_;
990 const auto& bw = brineCo2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p,
rsSolw());
992 const auto denw = bw*brineCo2Pvt.waterReferenceDensity(pvtRegionIdx)
993 +
rsSolw()*bw*brineCo2Pvt.gasReferenceDensity(pvtRegionIdx);
994 fs.setDensity(waterPhaseIdx, denw);
995 fs.setInvB(waterPhaseIdx, bw);
996 Evaluation& mobw =
asImp_().mobility_[waterPhaseIdx];
997 const auto& muWat = fs.viscosity(waterPhaseIdx);
998 const auto& muWatEff = brineCo2Pvt.viscosity(pvtRegionIdx, T, p,
rsSolw());
999 mobw *= muWat / muWatEff;
1007 auto& fs =
asImp_().fluidState_;
1008 const auto& bw = brineH2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p,
rsSolw());
1010 const auto denw = bw*brineH2Pvt.waterReferenceDensity(pvtRegionIdx)
1011 +
rsSolw()*bw*brineH2Pvt.gasReferenceDensity(pvtRegionIdx);
1012 fs.setDensity(waterPhaseIdx, denw);
1013 fs.setInvB(waterPhaseIdx, bw);
1014 Evaluation& mobw =
asImp_().mobility_[waterPhaseIdx];
1015 const auto& muWat = fs.viscosity(waterPhaseIdx);
1016 const auto& muWatEff = brineH2Pvt.viscosity(pvtRegionIdx, T, p,
rsSolw());
1017 mobw *= muWat / muWatEff;
1027 effectiveProperties(elemCtx, scvIdx, timeIdx);
1059 void effectiveProperties(
const ElementContext& elemCtx,
1073 auto& fs =
asImp_().fluidState_;
1076 Evaluation oilEffSat = 0.0;
1077 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1078 oilEffSat = fs.saturation(oilPhaseIdx);
1080 Evaluation gasEffSat = fs.saturation(gasPhaseIdx);
1082 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
1085 const Evaluation zero = 0.0;
1086 const Evaluation& sw = fs.saturation(waterPhaseIdx);
1087 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1088 oilEffSat = std::max(oilEffSat - sorwmis.eval(sw,
true), zero);
1090 gasEffSat = std::max(gasEffSat - sgcwmis.eval(sw,
true), zero);
1091 solventEffSat = std::max(solventEffSat - sgcwmis.eval(sw,
true), zero);
1093 const Evaluation oilGasSolventEffSat = oilEffSat + gasEffSat + solventEffSat;
1094 const Evaluation oilSolventEffSat = oilEffSat + solventEffSat;
1095 const Evaluation solventGasEffSat = solventEffSat + gasEffSat;
1102 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
1105 const Evaluation pmisc = pmiscTable.eval(p,
true);
1109 const Evaluation& muGas = fs.viscosity(gasPhaseIdx);
1112 assert(muGas.value() > 0);
1113 assert(muSolvent.value() > 0);
1114 const Evaluation muGasPow = pow(muGas, 0.25);
1115 const Evaluation muSolventPow = pow(muSolvent, 0.25);
1117 Evaluation muMixSolventGas = muGas;
1118 if (solventGasEffSat > cutOff)
1119 muMixSolventGas *= muSolvent / pow(((gasEffSat / solventGasEffSat) * muSolventPow) + ((solventEffSat / solventGasEffSat) * muGasPow) , 4.0);
1121 Evaluation muOil = 1.0;
1122 Evaluation muOilPow = 1.0;
1123 Evaluation muMixOilSolvent = 1.0;
1124 Evaluation muOilEff = 1.0;
1125 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1126 muOil = fs.viscosity(oilPhaseIdx);
1127 assert(muOil.value() > 0);
1128 muOilPow = pow(muOil, 0.25);
1129 muMixOilSolvent = muOil;
1130 if (oilSolventEffSat > cutOff)
1131 muMixOilSolvent *= muSolvent / pow(((oilEffSat / oilSolventEffSat) * muSolventPow) + ((solventEffSat / oilSolventEffSat) * muOilPow) , 4.0);
1133 muOilEff = pow(muOil,1.0 - tlMixParamMu) * pow(muMixOilSolvent, tlMixParamMu);
1135 Evaluation muMixSolventGasOil = muOil;
1136 if (oilGasSolventEffSat > cutOff)
1137 muMixSolventGasOil *= muSolvent * muGas / pow(((oilEffSat / oilGasSolventEffSat) * muSolventPow * muGasPow)
1138 + ((solventEffSat / oilGasSolventEffSat) * muOilPow * muGasPow) + ((gasEffSat / oilGasSolventEffSat) * muSolventPow * muOilPow), 4.0);
1140 Evaluation muGasEff = pow(muGas,1.0 - tlMixParamMu) * pow(muMixSolventGas, tlMixParamMu);
1141 Evaluation muSolventEff = pow(muSolvent,1.0 - tlMixParamMu) * pow(muMixSolventGasOil, tlMixParamMu);
1144 const Evaluation& rhoGas = fs.density(gasPhaseIdx);
1145 Evaluation rhoOil = 0.0;
1146 if (FluidSystem::phaseIsActive(oilPhaseIdx))
1147 rhoOil = fs.density(oilPhaseIdx);
1158 const Evaluation muOilEffPow = pow(pow(muOil, 1.0 - tlMixParamRho) * pow(muMixOilSolvent, tlMixParamRho), 0.25);
1159 const Evaluation muGasEffPow = pow(pow(muGas, 1.0 - tlMixParamRho) * pow(muMixSolventGas, tlMixParamRho), 0.25);
1160 const Evaluation muSolventEffPow = pow(pow(muSolvent, 1.0 - tlMixParamRho) * pow(muMixSolventGasOil, tlMixParamRho), 0.25);
1162 const Evaluation oilGasEffSaturation = oilEffSat + gasEffSat;
1163 Evaluation sof = 0.0;
1164 Evaluation sgf = 0.0;
1165 if (oilGasEffSaturation.value() > cutOff) {
1166 sof = oilEffSat / oilGasEffSaturation;
1167 sgf = gasEffSat / oilGasEffSaturation;
1170 const Evaluation muSolventOilGasPow = muSolventPow * ((sgf * muOilPow) + (sof * muGasPow));
1172 Evaluation rhoMixSolventGasOil = 0.0;
1173 if (oilGasSolventEffSat.value() > cutOff)
1174 rhoMixSolventGasOil = (rhoOil * oilEffSat / oilGasSolventEffSat) + (rhoGas * gasEffSat / oilGasSolventEffSat) + (rhoSolvent * solventEffSat / oilGasSolventEffSat);
1176 Evaluation rhoGasEff = 0.0;
1177 if (std::abs(muSolventPow.value() - muGasPow.value()) < cutOff)
1178 rhoGasEff = ((1.0 - tlMixParamRho) * rhoGas) + (tlMixParamRho * rhoMixSolventGasOil);
1180 const Evaluation solventGasEffFraction = (muGasPow * (muSolventPow - muGasEffPow)) / (muGasEffPow * (muSolventPow - muGasPow));
1181 rhoGasEff = (rhoGas * solventGasEffFraction) + (rhoSolvent * (1.0 - solventGasEffFraction));
1184 Evaluation rhoSolventEff = 0.0;
1185 if (std::abs((muSolventOilGasPow.value() - (muOilPow.value() * muGasPow.value()))) < cutOff)
1186 rhoSolventEff = ((1.0 - tlMixParamRho) * rhoSolvent) + (tlMixParamRho * rhoMixSolventGasOil);
1188 const Evaluation sfraction_se = (muSolventOilGasPow - (muOilPow * muGasPow * muSolventPow / muSolventEffPow)) / (muSolventOilGasPow - (muOilPow * muGasPow));
1189 rhoSolventEff = (rhoSolvent * sfraction_se) + (rhoGas * sgf * (1.0 - sfraction_se)) + (rhoOil * sof * (1.0 - sfraction_se));
1193 unsigned pvtRegionIdx =
asImp_().pvtRegionIndex();
1194 Evaluation bGasEff = rhoGasEff;
1195 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1196 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) * fs.Rv());
1198 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1203 const Evaluation bg = fs.invB(gasPhaseIdx);
1205 const Evaluation bg_eff = pmisc * bGasEff + (1.0 - pmisc) * bg;
1206 const Evaluation bs_eff = pmisc * bSolventEff + (1.0 - pmisc) * bs;
1209 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1210 fs.setDensity(gasPhaseIdx,
1212 *(FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)
1213 + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)*fs.Rv()));
1215 fs.setDensity(gasPhaseIdx,
1217 *FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1222 Evaluation& mobg =
asImp_().mobility_[gasPhaseIdx];
1223 muGasEff = bg_eff / (pmisc * bGasEff / muGasEff + (1.0 - pmisc) * bg / muGas);
1224 mobg *= muGas / muGasEff;
1227 solventViscosity_ = bs_eff / (pmisc * bSolventEff / muSolventEff + (1.0 - pmisc) * bs / muSolvent);
1229 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1230 Evaluation rhoOilEff = 0.0;
1231 if (std::abs(muOilPow.value() - muSolventPow.value()) < cutOff) {
1232 rhoOilEff = ((1.0 - tlMixParamRho) * rhoOil) + (tlMixParamRho * rhoMixSolventGasOil);
1235 const Evaluation solventOilEffFraction = (muOilPow * (muOilEffPow - muSolventPow)) / (muOilEffPow * (muOilPow - muSolventPow));
1236 rhoOilEff = (rhoOil * solventOilEffFraction) + (rhoSolvent * (1.0 - solventOilEffFraction));
1238 const Evaluation bOilEff = rhoOilEff / (FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) * fs.Rs());
1239 const Evaluation bo = fs.invB(oilPhaseIdx);
1240 const Evaluation bo_eff = pmisc * bOilEff + (1.0 - pmisc) * bo;
1241 fs.setDensity(oilPhaseIdx,
1243 *(FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)
1244 + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)*fs.Rs()));
1247 Evaluation& mobo =
asImp_().mobility_[oilPhaseIdx];
1248 muOilEff = bo_eff / (pmisc * bOilEff / muOilEff + (1.0 - pmisc) * bo / muOil);
1249 mobo *= muOil / muOilEff;
1255 {
return *
static_cast<Implementation*
>(
this); }
1268template <
class TypeTag>
1293 {
throw std::runtime_error(
"solventSaturation() called but solvents are disabled"); }
1296 {
throw std::runtime_error(
"rsSolw() called but solvents are disabled"); }
1299 {
throw std::runtime_error(
"solventDensity() called but solvents are disabled"); }
1302 {
throw std::runtime_error(
"solventViscosity() called but solvents are disabled"); }
1305 {
throw std::runtime_error(
"solventMobility() called but solvents are disabled"); }
1308 {
throw std::runtime_error(
"solventInverseFormationVolumeFactor() called but solvents are disabled"); }
1311 {
throw std::runtime_error(
"solventRefDensity() called but solvents are disabled"); }
1321template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
1334 using Toolbox = MathToolbox<Evaluation>;
1336 static constexpr unsigned gasPhaseIdx = FluidSystem::gasPhaseIdx;
1337 static constexpr int dimWorld = GridView::dimensionworld;
1339 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
1340 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
1347 template <
class Dummy =
bool>
1353 const auto& gradCalc = elemCtx.gradientCalculator();
1356 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
1357 const auto& faceNormal = scvf.normal();
1359 unsigned i = scvf.interiorIndex();
1360 unsigned j = scvf.exteriorIndex();
1363 DimEvalVector solventPGrad;
1365 gradCalc.calculateGradient(solventPGrad,
1369 Valgrind::CheckDefined(solventPGrad);
1372 if (Parameters::Get<Parameters::EnableGravity>()) {
1375 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
1376 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
1378 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
1379 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
1381 const auto& posIn = elemCtx.pos(i, timeIdx);
1382 const auto& posEx = elemCtx.pos(j, timeIdx);
1383 const auto& posFace = scvf.integrationPos();
1386 DimVector distVecIn(posIn);
1387 DimVector distVecEx(posEx);
1388 DimVector distVecTotal(posEx);
1390 distVecIn -= posFace;
1391 distVecEx -= posFace;
1392 distVecTotal -= posIn;
1393 Scalar absDistTotalSquared = distVecTotal.two_norm2();
1396 auto rhoIn = intQuantsIn.solventDensity();
1397 auto pStatIn = - rhoIn*(gIn*distVecIn);
1401 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1402 Scalar pStatEx = - rhoEx*(gEx*distVecEx);
1408 DimEvalVector f(distVecTotal);
1409 f *= (pStatEx - pStatIn)/absDistTotalSquared;
1412 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1413 solventPGrad[dimIdx] += f[dimIdx];
1415 if (!isfinite(solventPGrad[dimIdx]))
1416 throw NumericalProblem(
"Non-finite potential gradient for solvent 'phase'");
1421 Evaluation solventPGradNormal = 0.0;
1422 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
1423 solventPGradNormal += solventPGrad[dimIdx]*faceNormal[dimIdx];
1425 if (solventPGradNormal > 0) {
1426 solventUpstreamDofIdx_ = j;
1427 solventDownstreamDofIdx_ = i;
1430 solventUpstreamDofIdx_ = i;
1431 solventDownstreamDofIdx_ = j;
1434 const auto& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1440 if (solventUpstreamDofIdx_ == i)
1441 solventVolumeFlux_ = solventPGradNormal*up.solventMobility();
1443 solventVolumeFlux_ = solventPGradNormal*scalarValue(up.solventMobility());
1450 template <
class Dummy =
bool>
1456 const ExtensiveQuantities& extQuants = asImp_();
1458 unsigned interiorDofIdx = extQuants.interiorIndex();
1459 unsigned exteriorDofIdx = extQuants.exteriorIndex();
1460 assert(interiorDofIdx != exteriorDofIdx);
1462 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
1463 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
1465 unsigned I = elemCtx.globalSpaceIndex(interiorDofIdx, timeIdx);
1466 unsigned J = elemCtx.globalSpaceIndex(exteriorDofIdx, timeIdx);
1468 Scalar thpres = elemCtx.problem().thresholdPressure(I, J);
1469 Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
1470 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
1472 Scalar zIn = elemCtx.problem().dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
1473 Scalar zEx = elemCtx.problem().dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
1474 Scalar distZ = zIn - zEx;
1476 const Evaluation& rhoIn = intQuantsIn.solventDensity();
1477 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1478 const Evaluation& rhoAvg = rhoIn*0.5 + rhoEx*0.5;
1480 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(gasPhaseIdx);
1481 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(gasPhaseIdx));
1482 pressureExterior += distZ*g*rhoAvg;
1484 Evaluation pressureDiffSolvent = pressureExterior - pressureInterior;
1485 if (std::abs(scalarValue(pressureDiffSolvent)) > thpres) {
1486 if (pressureDiffSolvent < 0.0)
1487 pressureDiffSolvent += thpres;
1489 pressureDiffSolvent -= thpres;
1492 pressureDiffSolvent = 0.0;
1494 if (pressureDiffSolvent > 0.0) {
1495 solventUpstreamDofIdx_ = exteriorDofIdx;
1496 solventDownstreamDofIdx_ = interiorDofIdx;
1498 else if (pressureDiffSolvent < 0.0) {
1499 solventUpstreamDofIdx_ = interiorDofIdx;
1500 solventDownstreamDofIdx_ = exteriorDofIdx;
1506 solventUpstreamDofIdx_ = std::min(interiorDofIdx, exteriorDofIdx);
1507 solventDownstreamDofIdx_ = std::max(interiorDofIdx, exteriorDofIdx);
1508 solventVolumeFlux_ = 0.0;
1512 Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
1513 const IntensiveQuantities& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1514 if (solventUpstreamDofIdx_ == interiorDofIdx)
1515 solventVolumeFlux_ =
1516 up.solventMobility()
1518 *pressureDiffSolvent;
1520 solventVolumeFlux_ =
1521 scalarValue(up.solventMobility())
1523 *pressureDiffSolvent;
1527 {
return solventUpstreamDofIdx_; }
1530 {
return solventDownstreamDofIdx_; }
1533 {
return solventVolumeFlux_; }
1539 Implementation& asImp_()
1540 {
return *
static_cast<Implementation*
>(
this); }
1542 Evaluation solventVolumeFlux_;
1543 unsigned solventUpstreamDofIdx_;
1544 unsigned solventDownstreamDofIdx_;
1547template <
class TypeTag>
1565 {
throw std::runtime_error(
"solventUpstreamIndex() called but solvents are disabled"); }
1568 {
throw std::runtime_error(
"solventDownstreamIndex() called but solvents are disabled"); }
1571 {
throw std::runtime_error(
"solventVolumeFlux() called but solvents are disabled"); }
1574 {
throw std::runtime_error(
"setSolventVolumeFlux() called but solvents are disabled"); }
Declares the properties required by the black oil model.
Contains the parameters required to extend the black-oil model by solvents.
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1564
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1567
void setSolventVolumeFlux(const Evaluation &)
Definition: blackoilsolventmodules.hh:1573
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1570
void updateVolumeFluxPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1554
void updateVolumeFluxTrans(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1559
Provides the solvent specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilsolventmodules.hh:1323
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1526
void updateVolumeFluxPerm(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the volume flux of the polymer "phase" using the pressure potential gradient ...
Definition: blackoilsolventmodules.hh:1349
void updateVolumeFluxTrans(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the volume flux of the polymer "phase" using the gas pressure potential diffe...
Definition: blackoilsolventmodules.hh:1452
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1529
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1532
void setSolventVolumeFlux(const Evaluation &solventVolumeFlux)
Definition: blackoilsolventmodules.hh:1535
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:1310
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1298
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1307
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1295
void solventPreSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1277
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1301
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1304
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1292
void solventPvtUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1287
void solventPostSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1282
Provides the volumetric quantities required for the equations needed by the solvents extension of the...
Definition: blackoilsolventmodules.hh:780
void solventPreSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called before the saturation functions are doing their magic.
Definition: blackoilsolventmodules.hh:808
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1046
Evaluation solventViscosity_
Definition: blackoilsolventmodules.hh:1261
Evaluation solventMobility_
Definition: blackoilsolventmodules.hh:1262
Scalar solventRefDensity_
Definition: blackoilsolventmodules.hh:1265
Evaluation solventDensity_
Definition: blackoilsolventmodules.hh:1260
Evaluation rsSolw_
Definition: blackoilsolventmodules.hh:1259
Implementation & asImp_()
Definition: blackoilsolventmodules.hh:1254
Evaluation solventSaturation_
Definition: blackoilsolventmodules.hh:1258
void solventPvtUpdate_(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Update the intensive PVT properties needed to handle solvents from the primary variables.
Definition: blackoilsolventmodules.hh:970
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1040
Evaluation solventInvFormationVolumeFactor_
Definition: blackoilsolventmodules.hh:1263
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1037
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1043
void solventPostSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called after the saturation functions have been doing their magic.
Definition: blackoilsolventmodules.hh:838
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1049
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1034
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:1053
Evaluation hydrocarbonSaturation_
Definition: blackoilsolventmodules.hh:1257
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
static void setSolventPvt(const SolventPvt &value)
Specify the solvent PVT of a all PVT regions.
Definition: blackoilsolventmodules.hh:358
static const TabulatedFunction & sof2Krn(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:650
static const TabulatedFunction & ssfnKrs(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:642
static const TabulatedFunction & pmisc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:666
static const TabulatedFunction & msfnKrsg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:674
static const TabulatedFunction & sgcwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:698
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition: blackoilsolventmodules.hh:575
static const BrineH2Pvt & brineH2Pvt()
Definition: blackoilsolventmodules.hh:629
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:399
static std::string eqName(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:415
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:422
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:582
static bool isSolubleInWater()
Definition: blackoilsolventmodules.hh:748
static bool isCO2Sol()
Definition: blackoilsolventmodules.hh:753
static const TabulatedFunction & msfnKro(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:682
static void registerOutputModules(Model &model, Simulator &simulator)
Register all solvent specific VTK and ECL output modules.
Definition: blackoilsolventmodules.hh:377
static const TabulatedFunction & ssfnKrg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:634
static void setIsMiscible(const bool isMiscible)
Definition: blackoilsolventmodules.hh:362
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilsolventmodules.hh:368
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilsolventmodules.hh:563
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:608
static bool isH2Sol()
Definition: blackoilsolventmodules.hh:758
static bool isMiscible()
Definition: blackoilsolventmodules.hh:730
static const H2GasPvt & h2GasPvt()
Definition: blackoilsolventmodules.hh:619
static const Scalar & tlMixParamDensity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:722
static void updatePrimaryVars(PrimaryVariables &newPv, const PrimaryVariables &oldPv, const EqVector &delta)
Do a Newton-Raphson update the primary variables of the solvents.
Definition: blackoilsolventmodules.hh:551
static const TabulatedFunction & tlPMixTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:706
static const Scalar & tlMixParamViscosity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:714
static void assignPrimaryVars(PrimaryVariables &priVars, Scalar solventSaturation, Scalar solventRsw)
Assign the solvent specific primary variables to a PrimaryVariables object.
Definition: blackoilsolventmodules.hh:530
static const TabulatedFunction & sorwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:690
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:392
static bool eqApplies(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:407
static const Value solubilityLimit(unsigned pvtIdx, const Value &temperature, const Value &pressure, const Value &saltConcentration)
Definition: blackoilsolventmodules.hh:736
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:431
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:593
static const Co2GasPvt & co2GasPvt()
Definition: blackoilsolventmodules.hh:614
static const BrineCo2Pvt & brineCo2Pvt()
Definition: blackoilsolventmodules.hh:624
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:463
static const TabulatedFunction & misc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:658
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:384
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:84
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:108
VTK output module for the black oil model's solvent related quantities.
Definition: vtkblackoilsolventmodule.hh:63
Defines the common parameters for the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:37
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:235
This method contains all callback classes for quantities that are required by some extensive quantiti...
Struct holding the parameters for the BlackOilSolventModule class.
Definition: blackoilsolventparams.hh:43
::Opm::SolventPvt< Scalar > SolventPvt
Definition: blackoilsolventparams.hh:46
::Opm::H2GasPvt< Scalar > H2GasPvt
Definition: blackoilsolventparams.hh:52
::Opm::BrineCo2Pvt< Scalar > BrineCo2Pvt
Definition: blackoilsolventparams.hh:55
::Opm::Co2GasPvt< Scalar > Co2GasPvt
Definition: blackoilsolventparams.hh:49
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilsolventparams.hh:44
::Opm::BrineH2Pvt< Scalar > BrineH2Pvt
Definition: blackoilsolventparams.hh:58