28#ifndef EWOMS_BLACK_OIL_SOLVENT_MODULE_HH
29#define EWOMS_BLACK_OIL_SOLVENT_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/common/Exceptions.hpp>
34#include <opm/common/utility/gpuDecorators.hpp>
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/common/Valgrind.hpp>
38#include <opm/material/fluidsystems/blackoilpvt/SolventPvt.hpp>
67template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
81 using Toolbox = MathToolbox<Evaluation>;
90 static constexpr unsigned solventSaturationIdx = Indices::solventSaturationIdx;
91 static constexpr unsigned contiSolventEqIdx = Indices::contiSolventEqIdx;
92 static constexpr unsigned enableSolvent = enableSolventV;
93 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
94 static constexpr bool blackoilConserveSurfaceVolume =
95 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
96 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
99 static constexpr double cutOff = 1e-12;
102 { params_ = params; }
108 { params_.solventPvt_ = value; }
118 if constexpr (enableSolvent) {
127 Simulator& simulator)
129 if constexpr (enableSolvent) {
136 if constexpr (enableSolvent) {
137 return pvIdx == solventSaturationIdx;
148 return "saturation_solvent";
156 return static_cast<Scalar
>(1.0);
161 if constexpr (enableSolvent) {
162 return eqIdx == contiSolventEqIdx;
169 static std::string
eqName([[maybe_unused]]
unsigned eqIdx)
173 return "conti^solvent";
176 static Scalar
eqWeight([[maybe_unused]]
unsigned eqIdx)
181 return static_cast<Scalar
>(1.0);
184 template <
class StorageType>
186 const IntensiveQuantities& intQuants)
188 using LhsEval =
typename StorageType::value_type;
190 if constexpr (enableSolvent) {
191 if constexpr (blackoilConserveSurfaceVolume) {
192 storage[contiSolventEqIdx] +=
193 Toolbox::template decay<LhsEval>(intQuants.porosity()) *
194 Toolbox::template decay<LhsEval>(intQuants.solventSaturation()) *
195 Toolbox::template decay<LhsEval>(intQuants.solventInverseFormationVolumeFactor());
197 storage[contiSolventEqIdx] +=
198 Toolbox::template decay<LhsEval>(intQuants.porosity()) *
199 Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx)) *
200 Toolbox::template decay<LhsEval>(intQuants.fluidState().invB(waterPhaseIdx)) *
201 Toolbox::template decay<LhsEval>(intQuants.rsSolw());
205 storage[contiSolventEqIdx] +=
206 Toolbox::template decay<LhsEval>(intQuants.porosity()) *
207 Toolbox::template decay<LhsEval>(intQuants.solventSaturation()) *
208 Toolbox::template decay<LhsEval>(intQuants.solventDensity());
210 storage[contiSolventEqIdx] +=
211 Toolbox::template decay<LhsEval>(intQuants.porosity()) *
212 Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx)) *
213 Toolbox::template decay<LhsEval>(intQuants.fluidState().density(waterPhaseIdx)) *
214 Toolbox::template decay<LhsEval>(intQuants.rsSolw());
221 [[maybe_unused]]
const ElementContext& elemCtx,
222 [[maybe_unused]]
unsigned scvfIdx,
223 [[maybe_unused]]
unsigned timeIdx)
226 if constexpr (enableSolvent) {
227 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
229 const unsigned upIdx = extQuants.solventUpstreamIndex();
230 const unsigned inIdx = extQuants.interiorIndex();
231 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
233 if constexpr (blackoilConserveSurfaceVolume) {
234 if (upIdx == inIdx) {
235 flux[contiSolventEqIdx] = extQuants.solventVolumeFlux() *
236 up.solventInverseFormationVolumeFactor();
239 flux[contiSolventEqIdx] = extQuants.solventVolumeFlux() *
240 decay<Scalar>(up.solventInverseFormationVolumeFactor());
245 if (upIdx == inIdx) {
246 flux[contiSolventEqIdx] +=
247 extQuants.volumeFlux(waterPhaseIdx) *
248 up.fluidState().invB(waterPhaseIdx) *
252 flux[contiSolventEqIdx] +=
253 extQuants.volumeFlux(waterPhaseIdx) *
254 decay<Scalar>(up.fluidState().invB(waterPhaseIdx)) *
255 decay<Scalar>(up.rsSolw());
260 if (upIdx == inIdx) {
261 flux[contiSolventEqIdx] = extQuants.solventVolumeFlux() * up.solventDensity();
264 flux[contiSolventEqIdx] = extQuants.solventVolumeFlux() * decay<Scalar>(up.solventDensity());
268 if (upIdx == inIdx) {
269 flux[contiSolventEqIdx] +=
270 extQuants.volumeFlux(waterPhaseIdx) *
271 up.fluidState().density(waterPhaseIdx) *
275 flux[contiSolventEqIdx] +=
276 extQuants.volumeFlux(waterPhaseIdx) *
277 decay<Scalar>(up.fluidState().density(waterPhaseIdx)) *
278 decay<Scalar>(up.rsSolw());
289 Scalar solventSaturation,
292 if constexpr (!enableSolvent) {
293 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Disabled);
298 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Ss);
299 priVars[solventSaturationIdx] = solventSaturation;
301 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Rsolw);
302 priVars[solventSaturationIdx] = solventRsw;
310 const PrimaryVariables& oldPv,
311 const EqVector& delta)
313 if constexpr (enableSolvent) {
315 newPv[solventSaturationIdx] = oldPv[solventSaturationIdx] - delta[solventSaturationIdx];
328 return static_cast<Scalar
>(0.0);
337 return std::abs(Toolbox::scalarValue(resid[contiSolventEqIdx]));
340 template <
class DofEntity>
341 static void serializeEntity(
const Model& model, std::ostream& outstream,
const DofEntity& dof)
343 if constexpr (enableSolvent) {
344 const unsigned dofIdx = model.dofMapper().index(dof);
346 const PrimaryVariables& priVars = model.solution(0)[dofIdx];
347 outstream << priVars[solventSaturationIdx];
351 template <
class DofEntity>
354 if constexpr (enableSolvent) {
355 const unsigned dofIdx = model.dofMapper().index(dof);
357 PrimaryVariables& priVars0 = model.solution(0)[dofIdx];
358 PrimaryVariables& priVars1 = model.solution(1)[dofIdx];
360 instream >> priVars0[solventSaturationIdx];
363 priVars1 = priVars0[solventSaturationIdx];
368 {
return params_.solventPvt_; }
371 {
return params_.co2GasPvt_; }
374 {
return params_.h2GasPvt_; }
377 {
return params_.brineCo2Pvt_; }
380 {
return params_.brineH2Pvt_; }
382 template <
class ElemContext>
383 static const TabulatedFunction&
ssfnKrg(
const ElemContext& elemCtx,
387 const unsigned satnumRegionIdx =
388 elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
389 return params_.ssfnKrg_[satnumRegionIdx];
392 template <
class ElemContext>
393 static const TabulatedFunction&
ssfnKrs(
const ElemContext& elemCtx,
397 const unsigned satnumRegionIdx =
398 elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
399 return params_.ssfnKrs_[satnumRegionIdx];
402 static const TabulatedFunction&
ssfnKrg(
const unsigned satnumRegionIdx)
404 return params_.ssfnKrg_[satnumRegionIdx];
407 static const TabulatedFunction&
ssfnKrs(
const unsigned satnumRegionIdx)
409 return params_.ssfnKrs_[satnumRegionIdx];
412 template <
class ElemContext>
413 static const TabulatedFunction&
sof2Krn(
const ElemContext& elemCtx,
417 const unsigned satnumRegionIdx =
418 elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
419 return params_.sof2Krn_[satnumRegionIdx];
422 template <
class ElemContext>
423 static const TabulatedFunction&
misc(
const ElemContext& elemCtx,
427 const unsigned miscnumRegionIdx =
428 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
429 return params_.misc_[miscnumRegionIdx];
432 template <
class ElemContext>
433 static const TabulatedFunction&
pmisc(
const ElemContext& elemCtx,
437 const unsigned miscnumRegionIdx =
438 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
439 return params_.pmisc_[miscnumRegionIdx];
442 template <
class ElemContext>
443 static const TabulatedFunction&
msfnKrsg(
const ElemContext& elemCtx,
447 const unsigned satnumRegionIdx =
448 elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
449 return params_.msfnKrsg_[satnumRegionIdx];
452 template <
class ElemContext>
453 static const TabulatedFunction&
msfnKro(
const ElemContext& elemCtx,
457 const unsigned satnumRegionIdx =
458 elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
459 return params_.msfnKro_[satnumRegionIdx];
462 template <
class ElemContext>
463 static const TabulatedFunction&
sorwmis(
const ElemContext& elemCtx,
467 const unsigned miscnumRegionIdx =
468 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
469 return params_.sorwmis_[miscnumRegionIdx];
472 template <
class ElemContext>
473 static const TabulatedFunction&
sgcwmis(
const ElemContext& elemCtx,
477 const unsigned miscnumRegionIdx =
478 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
479 return params_.sgcwmis_[miscnumRegionIdx];
482 template <
class ElemContext>
483 static const TabulatedFunction&
tlPMixTable(
const ElemContext& elemCtx,
487 const unsigned miscnumRegionIdx =
488 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
489 return params_.tlPMixTable_[miscnumRegionIdx];
492 template <
class ElemContext>
497 const unsigned miscnumRegionIdx =
498 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
499 return params_.tlMixParamViscosity_[miscnumRegionIdx];
503 template <
class ElemContext>
508 const unsigned miscnumRegionIdx =
509 elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
510 return params_.tlMixParamDensity_[miscnumRegionIdx];
514 {
return params_.isMiscible_; }
516 template <
class Value>
518 const Value& pressure,
const Value& saltConcentration)
526 return brineCo2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature,
527 pressure, saltConcentration);
530 return brineH2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature,
531 pressure, saltConcentration);
536 {
return params_.rsSolw_active_; }
539 {
return params_.co2sol_; }
542 {
return params_.h2sol_; }
548template <
class TypeTag,
bool enableSolventV>
549BlackOilSolventParams<typename BlackOilSolventModule<TypeTag, enableSolventV>::Scalar>
550BlackOilSolventModule<TypeTag, enableSolventV>::params_;
552template <
class TypeTag,
bool enableSolventV>
562template <
class TypeTag>
578 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
579 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
580 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
581 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
582 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
583 static constexpr double cutOff = SolventModule::cutOff;
597 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
598 this->solventPreSatFuncUpdate_(priVars, timeIdx, elemCtx.linearizationType());
602 const unsigned timeIdx,
605 auto& fs = asImp_().fluidState_;
606 Evaluation solventSat{0.0};
607 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
608 solventSat = priVars.makeEvaluation(solventSaturationIdx, timeIdx, linearizationType);
611 fs.setSolventSaturation(solventSat);
613 hydrocarbonSaturation_ = fs.saturation(gasPhaseIdx);
616 if (solventSaturation().value() < cutOff) {
622 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_ + solventSat);
636 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
637 this->solventPostSatFuncUpdate_(elemCtx.problem(), priVars, elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx, elemCtx.linearizationType());
642 template <
class Problem>
643 struct ProblemAndCellIndexOnlyContext
645 const Problem& problem_;
647 const Problem& problem()
const {
return problem_; }
648 unsigned int globalSpaceIndex([[maybe_unused]]
const unsigned int spaceIdx,
649 [[maybe_unused]]
const unsigned int timeIdx)
const
657 const PrimaryVariables& priVars,
658 const unsigned globalSpaceIdx,
659 const unsigned timeIdx,
664 auto& fs = asImp_().fluidState_;
665 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_);
668 Evaluation rsSolw{0.0};
669 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
670 rsSolw = SolventModule::solubilityLimit(asImp_().pvtRegionIndex(),
671 fs.temperature(waterPhaseIdx),
672 fs.pressure(waterPhaseIdx),
673 fs.saltConcentration());
674 }
else if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Rsolw) {
675 rsSolw = priVars.makeEvaluation(solventSaturationIdx, timeIdx, linearizationType);
678 fs.setRsSolw(rsSolw);
680 solventMobility_ = 0.0;
683 if (solventSaturation().value() < cutOff) {
687 ProblemAndCellIndexOnlyContext<Problem> elemCtx{problem, globalSpaceIdx};
691 if (SolventModule::isMiscible()) {
692 const Evaluation& p =
693 FluidSystem::phaseIsActive(oilPhaseIdx)
694 ? fs.pressure(oilPhaseIdx)
695 : fs.pressure(gasPhaseIdx);
696 const Evaluation pmisc =
697 SolventModule::pmisc(elemCtx, dofIdx, timeIdx).eval(p,
true);
698 const Evaluation& pgImisc = fs.pressure(gasPhaseIdx);
701 Evaluation pgMisc = 0.0;
702 std::array<Evaluation, numPhases> pC;
703 const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
704 MaterialLaw::capillaryPressures(pC, materialParams, fs);
707 if (priVars.primaryVarsMeaningPressure() == PrimaryVariables::PressureMeaning::Pg) {
708 pgMisc = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx,
712 const Evaluation& po =
713 priVars.makeEvaluation(Indices::pressureSwitchIdx,
714 timeIdx, linearizationType);
715 pgMisc = po + (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
718 fs.setPressure(gasPhaseIdx, pmisc * pgMisc + (1.0 - pmisc) * pgImisc);
721 const Evaluation gasSolventSat = hydrocarbonSaturation_ + solventSaturation();
723 if (gasSolventSat.value() < cutOff) {
727 const Evaluation Fhydgas = hydrocarbonSaturation_ / gasSolventSat;
728 const Evaluation Fsolgas = solventSaturation() / gasSolventSat;
731 if (SolventModule::isMiscible() && FluidSystem::phaseIsActive(oilPhaseIdx)) {
732 const auto& misc = SolventModule::misc(elemCtx, dofIdx, timeIdx);
733 const auto& pmisc = SolventModule::pmisc(elemCtx, dofIdx, timeIdx);
734 const Evaluation& p =
735 FluidSystem::phaseIsActive(oilPhaseIdx)
736 ? fs.pressure(oilPhaseIdx)
737 : fs.pressure(gasPhaseIdx);
738 const Evaluation miscibility =
739 misc.eval(Fsolgas,
true) *
743 const unsigned cellIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
744 const auto& materialLawManager = elemCtx.problem().materialLawManager();
745 const auto& scaledDrainageInfo =
746 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
748 const Scalar sogcr = scaledDrainageInfo.Sogcr;
749 Evaluation sor = sogcr;
750 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
751 const Evaluation& sw = fs.saturation(waterPhaseIdx);
752 const auto& sorwmis = SolventModule::sorwmis(elemCtx, dofIdx, timeIdx);
754 sorwmis.eval(sw,
true) + (1.0 - miscibility) * sogcr;
756 const Scalar sgcr = scaledDrainageInfo.Sgcr;
757 Evaluation sgc = sgcr;
758 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
759 const Evaluation& sw = fs.saturation(waterPhaseIdx);
760 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, dofIdx, timeIdx);
761 sgc = miscibility * sgcwmis.eval(sw,
true) + (1.0 - miscibility) * sgcr;
764 Evaluation oilGasSolventSat = gasSolventSat;
765 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
766 oilGasSolventSat += fs.saturation(oilPhaseIdx);
768 const Evaluation zero = 0.0;
769 const Evaluation oilGasSolventEffSat = std::max(oilGasSolventSat - sor - sgc, zero);
771 Evaluation F_totalGas = 0.0;
772 if (oilGasSolventEffSat.value() > cutOff) {
773 const Evaluation gasSolventEffSat = std::max(gasSolventSat - sgc, zero);
774 F_totalGas = gasSolventEffSat / oilGasSolventEffSat;
776 const auto& msfnKro = SolventModule::msfnKro(elemCtx, dofIdx, timeIdx);
777 const auto& msfnKrsg = SolventModule::msfnKrsg(elemCtx, dofIdx, timeIdx);
778 const auto& sof2Krn = SolventModule::sof2Krn(elemCtx, dofIdx, timeIdx);
780 const Evaluation mkrgt = msfnKrsg.eval(F_totalGas,
true) *
781 sof2Krn.eval(oilGasSolventSat,
true);
782 const Evaluation mkro = msfnKro.eval(F_totalGas,
true) *
783 sof2Krn.eval(oilGasSolventSat,
true);
785 Evaluation& kro = asImp_().mobility_[oilPhaseIdx];
786 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
789 krg *= 1.0 - miscibility;
790 krg += miscibility * mkrgt;
791 kro *= 1.0 - miscibility;
792 kro += miscibility * mkro;
796 const auto& ssfnKrg = SolventModule::ssfnKrg(elemCtx, dofIdx, timeIdx);
797 const auto& ssfnKrs = SolventModule::ssfnKrs(elemCtx, dofIdx, timeIdx);
799 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
800 solventMobility_ = krg * ssfnKrs.eval(Fsolgas,
true);
801 krg *= ssfnKrg.eval(Fhydgas,
true);
814 const auto& iq = asImp_();
815 auto& fs = asImp_().fluidState_;
816 const unsigned pvtRegionIdx = iq.pvtRegionIndex();
817 const Evaluation& T = iq.fluidState().temperature(gasPhaseIdx);
818 const Evaluation& p = iq.fluidState().pressure(gasPhaseIdx);
820 Evaluation solventInvB;
821 const Evaluation rv = 0.0;
822 const Evaluation rvw = 0.0;
823 if (SolventModule::isCO2Sol() || SolventModule::isH2Sol() ){
824 if (SolventModule::isCO2Sol()) {
825 const auto& co2gasPvt = SolventModule::co2GasPvt();
827 co2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
828 solventRefDensity_ = co2gasPvt.gasReferenceDensity(pvtRegionIdx);
829 solventViscosity_ = co2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
831 const auto& brineCo2Pvt = SolventModule::brineCo2Pvt();
832 const auto& bw = brineCo2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
834 const auto denw = bw * brineCo2Pvt.waterReferenceDensity(pvtRegionIdx) +
835 rsSolw() * bw * brineCo2Pvt.gasReferenceDensity(pvtRegionIdx);
836 fs.setDensity(waterPhaseIdx, denw);
837 fs.setInvB(waterPhaseIdx, bw);
838 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
839 const auto& muWat = fs.viscosity(waterPhaseIdx);
840 const auto& muWatEff = brineCo2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
841 mobw *= muWat / muWatEff;
843 const auto& h2gasPvt = SolventModule::h2GasPvt();
845 h2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
846 solventRefDensity_ = h2gasPvt.gasReferenceDensity(pvtRegionIdx);
847 solventViscosity_ = h2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
849 const auto& brineH2Pvt = SolventModule::brineH2Pvt();
850 const auto& bw = brineH2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
852 const auto denw = bw * brineH2Pvt.waterReferenceDensity(pvtRegionIdx) +
853 rsSolw() * bw * brineH2Pvt.gasReferenceDensity(pvtRegionIdx);
854 fs.setDensity(waterPhaseIdx, denw);
855 fs.setInvB(waterPhaseIdx, bw);
856 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
857 const auto& muWat = fs.viscosity(waterPhaseIdx);
858 const auto& muWatEff = brineH2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
859 mobw *= muWat / muWatEff;
862 const auto& solventPvt = SolventModule::solventPvt();
863 solventInvB = solventPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p);
864 solventRefDensity_ = solventPvt.referenceDensity(pvtRegionIdx);
865 solventViscosity_ = solventPvt.viscosity(pvtRegionIdx, T, p);
869 fs.setSolventInvB(solventInvB);
870 fs.setSolventDensity(solventInvB * solventRefDensity_);
872 effectiveProperties(elemCtx, scvIdx, timeIdx);
874 solventMobility_ /= solventViscosity_;
878 {
return asImp_().fluidState_.solventSaturation(); }
881 {
return asImp_().fluidState_.rsSolw(); }
884 {
return asImp_().fluidState_.solventDensity(); }
887 {
return solventViscosity_; }
890 {
return solventMobility_; }
893 {
return asImp_().fluidState_.solventInvB(); }
897 {
return solventRefDensity_; }
902 void effectiveProperties(
const ElementContext& elemCtx,
906 if (!SolventModule::isMiscible()) {
912 if (solventSaturation() < cutOff) {
918 auto& fs = asImp_().fluidState_;
921 Evaluation oilEffSat = 0.0;
922 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
923 oilEffSat = fs.saturation(oilPhaseIdx);
925 Evaluation gasEffSat = fs.saturation(gasPhaseIdx);
926 Evaluation solventEffSat = solventSaturation();
927 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
928 const auto& sorwmis = SolventModule::sorwmis(elemCtx, scvIdx, timeIdx);
929 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, scvIdx, timeIdx);
930 const Evaluation zero = 0.0;
931 const Evaluation& sw = fs.saturation(waterPhaseIdx);
932 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
933 oilEffSat = std::max(oilEffSat - sorwmis.eval(sw,
true), zero);
935 gasEffSat = std::max(gasEffSat - sgcwmis.eval(sw,
true), zero);
936 solventEffSat = std::max(solventEffSat - sgcwmis.eval(sw,
true), zero);
938 const Evaluation oilGasSolventEffSat = oilEffSat + gasEffSat + solventEffSat;
939 const Evaluation oilSolventEffSat = oilEffSat + solventEffSat;
940 const Evaluation solventGasEffSat = solventEffSat + gasEffSat;
947 const Evaluation& p =
948 FluidSystem::phaseIsActive(oilPhaseIdx)
949 ? fs.pressure(oilPhaseIdx)
950 : fs.pressure(gasPhaseIdx);
952 const auto& pmiscTable = SolventModule::pmisc(elemCtx, scvIdx, timeIdx);
953 const Evaluation pmisc = pmiscTable.eval(p,
true);
954 const auto& tlPMixTable = SolventModule::tlPMixTable(elemCtx, scvIdx, timeIdx);
955 const Evaluation tlMixParamMu = SolventModule::tlMixParamViscosity(elemCtx, scvIdx, timeIdx) *
956 tlPMixTable.eval(p,
true);
958 const Evaluation& muGas = fs.viscosity(gasPhaseIdx);
959 const Evaluation& muSolvent = solventViscosity_;
961 assert(muGas.value() > 0);
962 assert(muSolvent.value() > 0);
963 const Evaluation muGasPow = pow(muGas, 0.25);
964 const Evaluation muSolventPow = pow(muSolvent, 0.25);
966 Evaluation muMixSolventGas = muGas;
967 if (solventGasEffSat > cutOff) {
968 muMixSolventGas *= muSolvent / pow(((gasEffSat / solventGasEffSat) * muSolventPow) +
969 ((solventEffSat / solventGasEffSat) * muGasPow), 4.0);
972 Evaluation muOil = 1.0;
973 Evaluation muOilPow = 1.0;
974 Evaluation muMixOilSolvent = 1.0;
975 Evaluation muOilEff = 1.0;
976 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
977 muOil = fs.viscosity(oilPhaseIdx);
978 assert(muOil.value() > 0);
979 muOilPow = pow(muOil, 0.25);
980 muMixOilSolvent = muOil;
981 if (oilSolventEffSat > cutOff) {
982 muMixOilSolvent *= muSolvent / pow(((oilEffSat / oilSolventEffSat) * muSolventPow) +
983 ((solventEffSat / oilSolventEffSat) * muOilPow), 4.0);
986 muOilEff = pow(muOil,1.0 - tlMixParamMu) * pow(muMixOilSolvent, tlMixParamMu);
988 Evaluation muMixSolventGasOil = muOil;
989 if (oilGasSolventEffSat > cutOff) {
990 muMixSolventGasOil *= muSolvent * muGas / pow(((oilEffSat / oilGasSolventEffSat) *
991 muSolventPow * muGasPow) +
992 ((solventEffSat / oilGasSolventEffSat) *
993 muOilPow * muGasPow) +
994 ((gasEffSat / oilGasSolventEffSat) *
995 muSolventPow * muOilPow), 4.0);
998 Evaluation muGasEff = pow(muGas,1.0 - tlMixParamMu) * pow(muMixSolventGas, tlMixParamMu);
999 const Evaluation muSolventEff = pow(muSolvent,1.0 - tlMixParamMu) * pow(muMixSolventGasOil, tlMixParamMu);
1002 const Evaluation& rhoGas = fs.density(gasPhaseIdx);
1003 Evaluation rhoOil = 0.0;
1004 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1005 rhoOil = fs.density(oilPhaseIdx);
1008 Evaluation rhoSolvent = solventDensity();
1014 const Evaluation tlMixParamRho = SolventModule::tlMixParamDensity(elemCtx, scvIdx, timeIdx) *
1015 tlPMixTable.eval(p,
true);
1019 const Evaluation muOilEffPow = pow(pow(muOil, 1.0 - tlMixParamRho) *
1020 pow(muMixOilSolvent, tlMixParamRho), 0.25);
1021 const Evaluation muGasEffPow = pow(pow(muGas, 1.0 - tlMixParamRho) *
1022 pow(muMixSolventGas, tlMixParamRho), 0.25);
1023 const Evaluation muSolventEffPow = pow(pow(muSolvent, 1.0 - tlMixParamRho) *
1024 pow(muMixSolventGasOil, tlMixParamRho), 0.25);
1026 const Evaluation oilGasEffSaturation = oilEffSat + gasEffSat;
1027 Evaluation sof = 0.0;
1028 Evaluation sgf = 0.0;
1029 if (oilGasEffSaturation.value() > cutOff) {
1030 sof = oilEffSat / oilGasEffSaturation;
1031 sgf = gasEffSat / oilGasEffSaturation;
1034 const Evaluation muSolventOilGasPow = muSolventPow * ((sgf * muOilPow) + (sof * muGasPow));
1036 Evaluation rhoMixSolventGasOil = 0.0;
1037 if (oilGasSolventEffSat.value() > cutOff) {
1038 rhoMixSolventGasOil = (rhoOil * oilEffSat / oilGasSolventEffSat) +
1039 (rhoGas * gasEffSat / oilGasSolventEffSat) +
1040 (rhoSolvent * solventEffSat / oilGasSolventEffSat);
1043 Evaluation rhoGasEff = 0.0;
1044 if (std::abs(muSolventPow.value() - muGasPow.value()) < cutOff) {
1045 rhoGasEff = ((1.0 - tlMixParamRho) * rhoGas) +
1046 (tlMixParamRho * rhoMixSolventGasOil);
1049 const Evaluation solventGasEffFraction = (muGasPow * (muSolventPow - muGasEffPow)) /
1050 (muGasEffPow * (muSolventPow - muGasPow));
1051 rhoGasEff = (rhoGas * solventGasEffFraction) + (rhoSolvent * (1.0 - solventGasEffFraction));
1054 Evaluation rhoSolventEff = 0.0;
1055 if (std::abs((muSolventOilGasPow.value() - (muOilPow.value() * muGasPow.value()))) < cutOff) {
1056 rhoSolventEff = ((1.0 - tlMixParamRho) * rhoSolvent) + (tlMixParamRho * rhoMixSolventGasOil);
1059 const Evaluation sfraction_se = (muSolventOilGasPow -
1060 muOilPow * muGasPow * muSolventPow / muSolventEffPow) /
1061 (muSolventOilGasPow - (muOilPow * muGasPow));
1062 rhoSolventEff = (rhoSolvent * sfraction_se) +
1063 (rhoGas * sgf * (1.0 - sfraction_se)) +
1064 (rhoOil * sof * (1.0 - sfraction_se));
1068 const unsigned pvtRegionIdx = asImp_().pvtRegionIndex();
1069 Evaluation bGasEff = rhoGasEff;
1070 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1071 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) +
1072 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) * fs.Rv());
1074 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1076 const Evaluation bSolventEff = rhoSolventEff / solventRefDensity();
1079 const Evaluation bg = fs.invB(gasPhaseIdx);
1080 const Evaluation bs = solventInverseFormationVolumeFactor();
1081 const Evaluation bg_eff = pmisc * bGasEff + (1.0 - pmisc) * bg;
1082 const Evaluation bs_eff = pmisc * bSolventEff + (1.0 - pmisc) * bs;
1085 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1086 fs.setDensity(gasPhaseIdx,
1088 (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) +
1089 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)*fs.Rv()));
1091 fs.setDensity(gasPhaseIdx,
1092 bg_eff * FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1094 rhoSolvent = bs_eff * solventRefDensity();
1095 fs.setSolventDensity(rhoSolvent);
1098 Evaluation& mobg = asImp_().mobility_[gasPhaseIdx];
1099 muGasEff = bg_eff / (pmisc * bGasEff / muGasEff + (1.0 - pmisc) * bg / muGas);
1100 mobg *= muGas / muGasEff;
1103 solventViscosity_ = bs_eff / (pmisc * bSolventEff / muSolventEff +
1104 (1.0 - pmisc) * bs / muSolvent);
1106 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1107 Evaluation rhoOilEff = 0.0;
1108 if (std::abs(muOilPow.value() - muSolventPow.value()) < cutOff) {
1109 rhoOilEff = ((1.0 - tlMixParamRho) * rhoOil) + (tlMixParamRho * rhoMixSolventGasOil);
1112 const Evaluation solventOilEffFraction = (muOilPow * (muOilEffPow - muSolventPow)) /
1113 (muOilEffPow * (muOilPow - muSolventPow));
1114 rhoOilEff = (rhoOil * solventOilEffFraction) + (rhoSolvent * (1.0 - solventOilEffFraction));
1116 const Evaluation bOilEff =
1117 rhoOilEff / (FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) +
1118 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) * fs.Rs());
1119 const Evaluation bo = fs.invB(oilPhaseIdx);
1120 const Evaluation bo_eff = pmisc * bOilEff + (1.0 - pmisc) * bo;
1121 fs.setDensity(oilPhaseIdx,
1122 bo_eff * (FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) +
1123 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) * fs.Rs()));
1126 Evaluation& mobo = asImp_().mobility_[oilPhaseIdx];
1127 muOilEff = bo_eff / (pmisc * bOilEff / muOilEff + (1.0 - pmisc) * bo / muOil);
1128 mobo *= muOil / muOilEff;
1134 {
return *
static_cast<Implementation*
>(
this); }
1137 {
return *
static_cast<const Implementation*
>(
this); }
1146template <
class TypeTag>
1170 {
throw std::runtime_error(
"solventSaturation() called but solvents are disabled"); }
1173 {
throw std::runtime_error(
"rsSolw() called but solvents are disabled"); }
1176 {
throw std::runtime_error(
"solventDensity() called but solvents are disabled"); }
1179 {
throw std::runtime_error(
"solventViscosity() called but solvents are disabled"); }
1182 {
throw std::runtime_error(
"solventMobility() called but solvents are disabled"); }
1185 {
throw std::runtime_error(
"solventInverseFormationVolumeFactor() called but solvents are disabled"); }
1188 {
throw std::runtime_error(
"solventRefDensity() called but solvents are disabled"); }
1198template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
1211 using Toolbox = MathToolbox<Evaluation>;
1213 static constexpr unsigned gasPhaseIdx = FluidSystem::gasPhaseIdx;
1214 static constexpr int dimWorld = GridView::dimensionworld;
1216 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
1217 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
1224 template <
class Dummy =
bool>
1230 const auto& gradCalc = elemCtx.gradientCalculator();
1233 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
1234 const auto& faceNormal = scvf.normal();
1236 const unsigned i = scvf.interiorIndex();
1237 const unsigned j = scvf.exteriorIndex();
1240 DimEvalVector solventPGrad;
1242 gradCalc.calculateGradient(solventPGrad,
1246 Valgrind::CheckDefined(solventPGrad);
1249 if (Parameters::Get<Parameters::EnableGravity>()) {
1252 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
1253 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
1255 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
1256 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
1258 const auto& posIn = elemCtx.pos(i, timeIdx);
1259 const auto& posEx = elemCtx.pos(j, timeIdx);
1260 const auto& posFace = scvf.integrationPos();
1263 DimVector distVecIn(posIn);
1264 DimVector distVecEx(posEx);
1265 DimVector distVecTotal(posEx);
1267 distVecIn -= posFace;
1268 distVecEx -= posFace;
1269 distVecTotal -= posIn;
1270 const Scalar absDistTotalSquared = distVecTotal.two_norm2();
1273 auto rhoIn = intQuantsIn.solventDensity();
1274 auto pStatIn = - rhoIn * (gIn * distVecIn);
1278 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1279 Scalar pStatEx = - rhoEx * (gEx * distVecEx);
1285 DimEvalVector f(distVecTotal);
1286 f *= (pStatEx - pStatIn) / absDistTotalSquared;
1289 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1290 solventPGrad[dimIdx] += f[dimIdx];
1292 if (!isfinite(solventPGrad[dimIdx])) {
1293 throw NumericalProblem(
"Non-finite potential gradient for solvent 'phase'");
1299 Evaluation solventPGradNormal = 0.0;
1300 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
1301 solventPGradNormal += solventPGrad[dimIdx]*faceNormal[dimIdx];
1304 if (solventPGradNormal > 0) {
1305 solventUpstreamDofIdx_ = j;
1306 solventDownstreamDofIdx_ = i;
1309 solventUpstreamDofIdx_ = i;
1310 solventDownstreamDofIdx_ = j;
1313 const auto& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1319 if (solventUpstreamDofIdx_ == i) {
1320 solventVolumeFlux_ = solventPGradNormal * up.solventMobility();
1323 solventVolumeFlux_ = solventPGradNormal * scalarValue(up.solventMobility());
1331 template <
class Dummy =
bool>
1337 const ExtensiveQuantities& extQuants = asImp_();
1339 const unsigned interiorDofIdx = extQuants.interiorIndex();
1340 const unsigned exteriorDofIdx = extQuants.exteriorIndex();
1341 assert(interiorDofIdx != exteriorDofIdx);
1343 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
1344 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
1346 const unsigned I = elemCtx.globalSpaceIndex(interiorDofIdx, timeIdx);
1347 const unsigned J = elemCtx.globalSpaceIndex(exteriorDofIdx, timeIdx);
1349 const Scalar thpres = elemCtx.problem().thresholdPressure(I, J);
1350 const Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
1351 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
1353 const Scalar zIn = elemCtx.problem().dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
1354 const Scalar zEx = elemCtx.problem().dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
1355 const Scalar distZ = zIn - zEx;
1357 const Evaluation rhoIn = intQuantsIn.solventDensity();
1358 const Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1359 const Evaluation rhoAvg = rhoIn * 0.5 + rhoEx * 0.5;
1361 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(gasPhaseIdx);
1362 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(gasPhaseIdx));
1363 pressureExterior += distZ * g * rhoAvg;
1365 Evaluation pressureDiffSolvent = pressureExterior - pressureInterior;
1366 if (std::abs(scalarValue(pressureDiffSolvent)) > thpres) {
1367 if (pressureDiffSolvent < 0.0) {
1368 pressureDiffSolvent += thpres;
1371 pressureDiffSolvent -= thpres;
1375 pressureDiffSolvent = 0.0;
1378 if (pressureDiffSolvent > 0.0) {
1379 solventUpstreamDofIdx_ = exteriorDofIdx;
1380 solventDownstreamDofIdx_ = interiorDofIdx;
1382 else if (pressureDiffSolvent < 0.0) {
1383 solventUpstreamDofIdx_ = interiorDofIdx;
1384 solventDownstreamDofIdx_ = exteriorDofIdx;
1390 solventUpstreamDofIdx_ = std::min(interiorDofIdx, exteriorDofIdx);
1391 solventDownstreamDofIdx_ = std::max(interiorDofIdx, exteriorDofIdx);
1392 solventVolumeFlux_ = 0.0;
1396 const Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
1397 const IntensiveQuantities& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1398 if (solventUpstreamDofIdx_ == interiorDofIdx) {
1399 solventVolumeFlux_ = up.solventMobility() * (-trans / faceArea) * pressureDiffSolvent;
1402 solventVolumeFlux_ =
1403 scalarValue(up.solventMobility()) * (-trans / faceArea) * pressureDiffSolvent;
1408 {
return solventUpstreamDofIdx_; }
1411 {
return solventDownstreamDofIdx_; }
1414 {
return solventVolumeFlux_; }
1420 Implementation& asImp_()
1421 {
return *
static_cast<Implementation*
>(
this); }
1423 Evaluation solventVolumeFlux_;
1424 unsigned solventUpstreamDofIdx_;
1425 unsigned solventDownstreamDofIdx_;
1428template <
class TypeTag>
1446 {
throw std::runtime_error(
"solventUpstreamIndex() called but solvents are disabled"); }
1449 {
throw std::runtime_error(
"solventDownstreamIndex() called but solvents are disabled"); }
1452 {
throw std::runtime_error(
"solventVolumeFlux() called but solvents are disabled"); }
1455 {
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:1445
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1448
void setSolventVolumeFlux(const Evaluation &)
Definition: blackoilsolventmodules.hh:1454
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1451
void updateVolumeFluxPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1435
void updateVolumeFluxTrans(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1440
Provides the solvent specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilsolventmodules.hh:1200
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1407
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:1226
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:1333
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1410
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1413
void setSolventVolumeFlux(const Evaluation &solventVolumeFlux)
Definition: blackoilsolventmodules.hh:1416
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1175
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1184
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1172
void solventPreSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1154
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1178
Scalar solventRefDensity() const
Definition: blackoilsolventmodules.hh:1187
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1181
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1169
void solventPvtUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1164
void solventPostSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1159
void solventPreSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called before the saturation functions are doing their magic.
Definition: blackoilsolventmodules.hh:593
Implementation & asImp_()
Definition: blackoilsolventmodules.hh:1133
Evaluation solventViscosity_
Definition: blackoilsolventmodules.hh:1140
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:886
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:810
const Implementation & asImp_() const
Definition: blackoilsolventmodules.hh:1136
void solventPostSatFuncUpdate_(const Problem &problem, const PrimaryVariables &priVars, const unsigned globalSpaceIdx, const unsigned timeIdx, const LinearizationType linearizationType)
Definition: blackoilsolventmodules.hh:656
Evaluation solventDensity() const
Definition: blackoilsolventmodules.hh:883
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:889
Evaluation hydrocarbonSaturation_
Definition: blackoilsolventmodules.hh:1139
void solventPreSatFuncUpdate_(const PrimaryVariables &priVars, const unsigned timeIdx, const LinearizationType linearizationType)
Definition: blackoilsolventmodules.hh:601
void solventPostSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called after the saturation functions have been doing their magic.
Definition: blackoilsolventmodules.hh:632
Evaluation solventSaturation() const
Definition: blackoilsolventmodules.hh:877
Scalar solventRefDensity_
Definition: blackoilsolventmodules.hh:1143
Scalar solventRefDensity() const
Definition: blackoilsolventmodules.hh:896
Evaluation rsSolw() const
Definition: blackoilsolventmodules.hh:880
Evaluation solventMobility_
Definition: blackoilsolventmodules.hh:1141
Evaluation solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:892
Provides the volumetric quantities required for the equations needed by the solvents extension of the...
Definition: blackoilsolventmodules.hh:553
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:107
static const TabulatedFunction & misc(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:423
static const TabulatedFunction & sof2Krn(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:413
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition: blackoilsolventmodules.hh:334
static Scalar tlMixParamDensity(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:504
static const TabulatedFunction & ssfnKrg(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:383
static const BrineH2Pvt & brineH2Pvt()
Definition: blackoilsolventmodules.hh:379
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:151
static void setParams(BlackOilSolventParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilsolventmodules.hh:101
static std::string eqName(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:169
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:176
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:341
static bool isSolubleInWater()
Definition: blackoilsolventmodules.hh:535
static bool isCO2Sol()
Definition: blackoilsolventmodules.hh:538
static const TabulatedFunction & sgcwmis(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:473
static Scalar tlMixParamViscosity(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:493
static void registerOutputModules(Model &model, Simulator &simulator)
Register all solvent specific VTK and ECL output modules.
Definition: blackoilsolventmodules.hh:126
static const TabulatedFunction & tlPMixTable(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:483
static const TabulatedFunction & sorwmis(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:463
static Value solubilityLimit(unsigned pvtIdx, const Value &temperature, const Value &pressure, const Value &saltConcentration)
Definition: blackoilsolventmodules.hh:517
static const TabulatedFunction & pmisc(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:433
static void setIsMiscible(const bool isMiscible)
Definition: blackoilsolventmodules.hh:110
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilsolventmodules.hh:116
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilsolventmodules.hh:322
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:367
static bool isH2Sol()
Definition: blackoilsolventmodules.hh:541
static constexpr double cutOff
Definition: blackoilsolventmodules.hh:99
static bool isMiscible()
Definition: blackoilsolventmodules.hh:513
static const H2GasPvt & h2GasPvt()
Definition: blackoilsolventmodules.hh:373
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:309
static void assignPrimaryVars(PrimaryVariables &priVars, Scalar solventSaturation, Scalar solventRsw)
Assign the solvent specific primary variables to a PrimaryVariables object.
Definition: blackoilsolventmodules.hh:288
static const TabulatedFunction & msfnKrsg(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:443
static const TabulatedFunction & ssfnKrs(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:393
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:144
static bool eqApplies(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:159
static const TabulatedFunction & msfnKro(const ElemContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:453
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:352
static const TabulatedFunction & ssfnKrg(const unsigned satnumRegionIdx)
Definition: blackoilsolventmodules.hh:402
static const Co2GasPvt & co2GasPvt()
Definition: blackoilsolventmodules.hh:370
static const BrineCo2Pvt & brineCo2Pvt()
Definition: blackoilsolventmodules.hh:376
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:185
static const TabulatedFunction & ssfnKrs(const unsigned satnumRegionIdx)
Definition: blackoilsolventmodules.hh:407
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:220
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:134
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:85
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:109
VTK output module for the black oil model's solvent related quantities.
Definition: vtkblackoilsolventmodule.hpp:54
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilsolventmodule.hpp:84
Defines the common parameters for the porous medium multi-phase models.
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
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.hpp:47
::Opm::SolventPvt< Scalar > SolventPvt
Definition: blackoilsolventparams.hpp:52
::Opm::H2GasPvt< Scalar > H2GasPvt
Definition: blackoilsolventparams.hpp:51
::Opm::BrineCo2Pvt< Scalar > BrineCo2Pvt
Definition: blackoilsolventparams.hpp:48
::Opm::Co2GasPvt< Scalar > Co2GasPvt
Definition: blackoilsolventparams.hpp:50
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilsolventparams.hpp:53
::Opm::BrineH2Pvt< Scalar > BrineH2Pvt
Definition: blackoilsolventparams.hpp:49
Definition: linearizationtype.hh:34