28#ifndef EWOMS_BLACK_OIL_LOCAL_TPFA_RESIDUAL_HH
29#define EWOMS_BLACK_OIL_LOCAL_TPFA_RESIDUAL_HH
31#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
32#include <opm/input/eclipse/Schedule/BCProp.hpp>
34#include <opm/material/common/MathToolbox.hpp>
35#include <opm/material/fluidstates/BlackOilFluidState.hpp>
36#include <opm/material/common/ConditionalStorage.hpp>
62template <
class TypeTag>
75 using FluidState =
typename IntensiveQuantities::FluidState;
77 enum { conti0EqIdx = Indices::conti0EqIdx };
78 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
79 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
80 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
82 enum { dimWorld = GridView::dimensionworld };
83 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
84 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
85 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
87 enum { gasCompIdx = FluidSystem::gasCompIdx };
88 enum { oilCompIdx = FluidSystem::oilCompIdx };
89 enum { waterCompIdx = FluidSystem::waterCompIdx };
90 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
92 static constexpr bool waterEnabled = Indices::waterEnabled;
93 static constexpr bool gasEnabled = Indices::gasEnabled;
94 static constexpr bool oilEnabled = Indices::oilEnabled;
95 static constexpr bool compositionSwitchEnabled = compositionSwitchIdx >= 0;
97 static constexpr bool blackoilConserveSurfaceVolume =
98 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
100 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
101 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
102 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
103 static constexpr bool enableFullyImplicitThermal = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
104 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
105 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
106 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
107 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
108 static constexpr bool enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
109 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
110 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
111 static constexpr bool enableMICP = Indices::enableMICP;
126 using Toolbox = MathToolbox<Evaluation>;
138 ConditionalStorage<enableFullyImplicitThermal, double>
inAlpha;
139 ConditionalStorage<enableFullyImplicitThermal, double>
outAlpha;
147 template <
class LhsEval>
149 const ElementContext& elemCtx,
151 unsigned timeIdx)
const
153 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
157 template <
class LhsEval>
159 const IntensiveQuantities& intQuants)
163 const auto& fs = intQuants.fluidState();
166 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
167 if (!FluidSystem::phaseIsActive(phaseIdx)) {
170 unsigned activeCompIdx =
171 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
172 LhsEval surfaceVolume =
173 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
174 Toolbox::template decay<LhsEval>(fs.invB(phaseIdx)) *
175 Toolbox::template decay<LhsEval>(intQuants.porosity());
177 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
180 if (phaseIdx == oilPhaseIdx && FluidSystem::enableDissolvedGas()) {
181 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
182 storage[conti0EqIdx + activeGasCompIdx] +=
183 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs()) *
188 if (phaseIdx == waterPhaseIdx && FluidSystem::enableDissolvedGasInWater()) {
189 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
190 storage[conti0EqIdx + activeGasCompIdx] +=
191 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw()) *
196 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedOil()) {
197 unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
198 storage[conti0EqIdx + activeOilCompIdx] +=
199 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv()) *
204 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedWater()) {
205 unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
206 storage[conti0EqIdx + activeWaterCompIdx] +=
207 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw()) *
243 const unsigned globalIndexIn,
244 const unsigned globalIndexEx,
245 const IntensiveQuantities& intQuantsIn,
246 const IntensiveQuantities& intQuantsEx,
250 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
268 const ElementContext& elemCtx,
272 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
273 assert(timeIdx == 0);
276 RateVector darcy = 0.0;
278 const auto& problem = elemCtx.problem();
279 const auto& stencil = elemCtx.stencil(timeIdx);
280 const auto& scvf = stencil.interiorFace(scvfIdx);
282 unsigned interiorDofIdx = scvf.interiorIndex();
283 unsigned exteriorDofIdx = scvf.exteriorIndex();
284 assert(interiorDofIdx != exteriorDofIdx);
288 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, 0);
289 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, 0);
290 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
291 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
292 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
293 Scalar faceArea = scvf.area();
295 Scalar thpres = problem.thresholdPressure(globalIndexIn, globalIndexEx);
300 const Scalar g = problem.gravity()[dimWorld - 1];
301 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
302 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
309 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
310 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
313 const Scalar distZ = zIn - zEx;
315 const Scalar inAlpha = problem.thermalHalfTransmissibility(globalIndexIn, globalIndexEx);
316 const Scalar outAlpha = problem.thermalHalfTransmissibility(globalIndexEx, globalIndexIn);
317 const Scalar diffusivity = problem.diffusivity(globalIndexEx, globalIndexIn);
318 const Scalar dispersivity = problem.dispersivity(globalIndexEx, globalIndexIn);
321 trans, faceArea, thpres, distZ * g, faceDir, Vin, Vex,
322 inAlpha, outAlpha, diffusivity, dispersivity
332 problem.moduleParams());
337 const IntensiveQuantities& intQuantsIn,
338 const IntensiveQuantities& intQuantsEx,
339 const unsigned& globalIndexIn,
340 const unsigned& globalIndexEx,
344 OPM_TIMEBLOCK_LOCAL(calculateFluxes, Subsystem::Assembly);
345 const Scalar Vin = nbInfo.
Vin;
346 const Scalar Vex = nbInfo.
Vex;
347 const Scalar distZg = nbInfo.
dZg;
348 const Scalar thpres = nbInfo.
thpres;
349 const Scalar trans = nbInfo.
trans;
350 const Scalar faceArea = nbInfo.
faceArea;
351 FaceDir::DirEnum facedir = nbInfo.
faceDir;
353 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
354 if (!FluidSystem::phaseIsActive(phaseIdx)) {
362 short interiorDofIdx = 0;
363 short exteriorDofIdx = 1;
364 Evaluation pressureDifference;
365 ExtensiveQuantities::calculatePhasePressureDiff_(upIdx,
381 const IntensiveQuantities& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
382 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
384 Evaluation transMult = (intQuantsIn.rockCompTransMultiplier() +
385 Toolbox::value(intQuantsEx.rockCompTransMultiplier())) / 2;
386 if constexpr (enableBioeffects || enableSaltPrecipitation) {
387 transMult *= (intQuantsIn.permFactor() + Toolbox::value(intQuantsEx.permFactor())) / 2;
389 Evaluation darcyFlux;
390 if (globalUpIndex == globalIndexIn) {
391 darcyFlux = pressureDifference * up.mobility(phaseIdx, facedir) * transMult * (-trans / faceArea);
393 darcyFlux = pressureDifference *
394 (Toolbox::value(up.mobility(phaseIdx, facedir)) * transMult * (-trans / faceArea));
397 unsigned activeCompIdx =
398 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
400 darcy[conti0EqIdx + activeCompIdx] = darcyFlux.value() * faceArea;
402 unsigned pvtRegionIdx = up.pvtRegionIndex();
404 if (globalUpIndex == globalIndexIn) {
406 = getInvB_<FluidSystem, FluidState, Evaluation>(up.fluidState(), phaseIdx, pvtRegionIdx);
407 const auto& surfaceVolumeFlux = invB * darcyFlux;
408 evalPhaseFluxes_<Evaluation>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
409 if constexpr (enableFullyImplicitThermal) {
410 EnergyModule::template
411 addPhaseEnthalpyFluxes_<Evaluation>(flux, phaseIdx, darcyFlux, up.fluidState());
413 if constexpr (enableBioeffects) {
414 BioeffectsModule::template
415 addBioeffectsFluxes_<Evaluation>(flux, phaseIdx, darcyFlux, up);
417 if constexpr (enableBrine) {
418 BrineModule::template
419 addBrineFluxes_<Evaluation, FluidState>(flux, phaseIdx, darcyFlux, up.fluidState());
422 const auto& invB = getInvB_<FluidSystem, FluidState, Scalar>(up.fluidState(), phaseIdx, pvtRegionIdx);
423 const auto& surfaceVolumeFlux = invB * darcyFlux;
424 evalPhaseFluxes_<Scalar>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
425 if constexpr (enableFullyImplicitThermal) {
426 EnergyModule::template
427 addPhaseEnthalpyFluxes_<Scalar>(flux, phaseIdx, darcyFlux, up.fluidState());
429 if constexpr (enableBioeffects) {
430 BioeffectsModule::template
431 addBioeffectsFluxes_<Scalar>(flux, phaseIdx, darcyFlux, up);
433 if constexpr (enableBrine) {
434 BrineModule::template
435 addBrineFluxes_<Scalar, FluidState>(flux, phaseIdx, darcyFlux, up.fluidState());
441 static_assert(!enableSolvent,
442 "Relevant computeFlux() method must be implemented for this module before enabling.");
446 static_assert(!enableExtbo,
447 "Relevant computeFlux() method must be implemented for this module before enabling.");
451 static_assert(!enablePolymer,
452 "Relevant computeFlux() method must be implemented for this module before enabling.");
456 if constexpr (enableConvectiveMixing) {
457 ConvectiveMixingModule::addConvectiveMixingFlux(flux,
469 if constexpr (enableFullyImplicitThermal) {
470 const Scalar inAlpha = nbInfo.
inAlpha;
471 const Scalar outAlpha = nbInfo.
outAlpha;
474 short interiorDofIdx = 0;
475 short exteriorDofIdx = 1;
477 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
483 intQuantsIn.fluidState(),
484 intQuantsEx.fluidState(),
494 static_assert(!enableFoam,
495 "Relevant computeFlux() method must be implemented for this module before enabling.");
499 if constexpr (enableDiffusion) {
500 typename DiffusionModule::ExtensiveQuantities::EvaluationArray effectiveDiffusionCoefficient;
501 DiffusionModule::ExtensiveQuantities::update(effectiveDiffusionCoefficient, intQuantsIn, intQuantsEx);
503 const Scalar tmpdiffusivity = diffusivity / faceArea;
504 DiffusionModule::addDiffusiveFlux(flux,
508 effectiveDiffusionCoefficient);
512 if constexpr (enableDispersion) {
513 typename DispersionModule::ExtensiveQuantities::ScalarArray normVelocityAvg;
514 DispersionModule::ExtensiveQuantities::update(normVelocityAvg, intQuantsIn, intQuantsEx);
516 const Scalar tmpdispersivity = dispersivity / faceArea;
517 DispersionModule::addDispersiveFlux(flux,
525 if constexpr (enableMICP) {
530 template <
class BoundaryConditionData>
532 const Problem& problem,
533 const BoundaryConditionData& bdyInfo,
534 const IntensiveQuantities& insideIntQuants,
535 unsigned globalSpaceIdx)
537 switch (bdyInfo.type) {
545 case BCType::DIRICHLET:
548 case BCType::THERMAL:
552 throw std::logic_error(
"Unknown boundary condition type " +
554 " in computeBoundaryFlux()." );
558 template <
class BoundaryConditionData>
560 const BoundaryConditionData& bdyInfo)
562 bdyFlux.setMassRate(bdyInfo.massRate, bdyInfo.pvtRegionIdx);
565 template <
class BoundaryConditionData>
568 const BoundaryConditionData& bdyInfo,
569 const IntensiveQuantities& insideIntQuants,
570 unsigned globalSpaceIdx)
573 std::array<short, numPhases> upIdx;
574 std::array<short, numPhases> dnIdx;
575 std::array<Evaluation, numPhases> volumeFlux;
576 std::array<Evaluation, numPhases> pressureDifference;
578 ExtensiveQuantities::calculateBoundaryGradients_(problem,
581 bdyInfo.boundaryFaceIndex,
584 bdyInfo.exFluidState,
594 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
595 if (!FluidSystem::phaseIsActive(phaseIdx)) {
598 const auto& pBoundary = bdyInfo.exFluidState.pressure(phaseIdx);
599 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
600 const unsigned pvtRegionIdx = insideIntQuants.pvtRegionIndex();
603 const auto& darcyFlux = volumeFlux[phaseIdx];
605 if (pBoundary < pInside) {
608 getInvB_<FluidSystem, FluidState, Evaluation>(insideIntQuants.fluidState(), phaseIdx, pvtRegionIdx);
609 Evaluation surfaceVolumeFlux = invB * darcyFlux;
610 evalPhaseFluxes_<Evaluation>(tmp,
612 insideIntQuants.pvtRegionIndex(),
614 insideIntQuants.fluidState());
615 if constexpr (enableFullyImplicitThermal) {
616 EnergyModule::template
617 addPhaseEnthalpyFluxes_<Evaluation>(tmp, phaseIdx, darcyFlux, insideIntQuants.fluidState());
619 }
else if (pBoundary > pInside) {
621 using ScalarFluidState =
decltype(bdyInfo.exFluidState);
623 getInvB_<FluidSystem, ScalarFluidState, Scalar>(bdyInfo.exFluidState, phaseIdx, pvtRegionIdx);
624 Evaluation surfaceVolumeFlux = invB * darcyFlux;
625 evalPhaseFluxes_<Scalar>(tmp,
627 insideIntQuants.pvtRegionIndex(),
629 bdyInfo.exFluidState);
630 if constexpr (enableFullyImplicitThermal) {
631 EnergyModule::template
632 addPhaseEnthalpyFluxes_<Scalar>(tmp, phaseIdx, darcyFlux, bdyInfo.exFluidState);
636 for (
unsigned i = 0; i < tmp.size(); ++i) {
637 bdyFlux[i] += tmp[i];
642 if constexpr (enableFullyImplicitThermal) {
646 problem.eclTransmissibilities().thermalHalfTransBoundary(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
649 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
654 bdyInfo.exFluidState);
658 static_assert(!enableSolvent,
659 "Relevant treatment of boundary conditions must be implemented before enabling.");
660 static_assert(!enablePolymer,
661 "Relevant treatment of boundary conditions must be implemented before enabling.");
667 for (
unsigned i = 0; i < numEq; ++i) {
668 Valgrind::CheckDefined(bdyFlux[i]);
670 Valgrind::CheckDefined(bdyFlux);
674 template <
class BoundaryConditionData>
677 const BoundaryConditionData& bdyInfo,
678 const IntensiveQuantities& insideIntQuants,
679 [[maybe_unused]]
unsigned globalSpaceIdx)
686 if constexpr (enableFullyImplicitThermal) {
690 problem.eclTransmissibilities().thermalHalfTransBoundary(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
693 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
698 bdyInfo.exFluidState);
703 for (
unsigned i = 0; i < numEq; ++i) {
704 Valgrind::CheckDefined(bdyFlux[i]);
706 Valgrind::CheckDefined(bdyFlux);
711 const Problem& problem,
712 const IntensiveQuantities& insideIntQuants,
713 unsigned globalSpaceIdex,
718 problem.source(source, globalSpaceIdex, timeIdx);
724 if constexpr (enableFullyImplicitThermal) {
725 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
730 const Problem& problem,
731 const IntensiveQuantities& insideIntQuants,
732 unsigned globalSpaceIdex,
736 problem.addToSourceDense(source, globalSpaceIdex, timeIdx);
742 if constexpr (enableFullyImplicitThermal) {
743 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
751 const ElementContext& elemCtx,
753 unsigned timeIdx)
const
757 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
763 if constexpr (enableFullyImplicitThermal) {
764 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
768 template <
class UpEval,
class Flu
idState>
771 unsigned pvtRegionIdx,
772 const ExtensiveQuantities& extQuants,
773 const FluidState& upFs)
775 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
776 const auto& surfaceVolumeFlux = invB * extQuants.volumeFlux(phaseIdx);
777 evalPhaseFluxes_<UpEval>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, upFs);
784 template <
class UpEval,
class Eval,
class Flu
idState>
787 unsigned pvtRegionIdx,
788 const Eval& surfaceVolumeFlux,
789 const FluidState& upFs)
791 unsigned activeCompIdx =
792 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
794 if constexpr (blackoilConserveSurfaceVolume) {
795 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
798 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux *
799 FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
802 if (phaseIdx == oilPhaseIdx) {
804 if (FluidSystem::enableDissolvedGas()) {
805 const auto& Rs = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
807 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
808 if constexpr (blackoilConserveSurfaceVolume) {
809 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
812 flux[conti0EqIdx + activeGasCompIdx] +=
813 Rs * surfaceVolumeFlux *
814 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
818 else if (phaseIdx == waterPhaseIdx) {
820 if (FluidSystem::enableDissolvedGasInWater()) {
821 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
823 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
824 if constexpr (blackoilConserveSurfaceVolume) {
825 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
828 flux[conti0EqIdx + activeGasCompIdx] +=
829 Rsw * surfaceVolumeFlux *
830 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
834 else if (phaseIdx == gasPhaseIdx) {
836 if (FluidSystem::enableVaporizedOil()) {
837 const auto& Rv = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
839 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
840 if constexpr (blackoilConserveSurfaceVolume) {
841 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
844 flux[conti0EqIdx + activeOilCompIdx] +=
845 Rv * surfaceVolumeFlux *
846 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
850 if (FluidSystem::enableVaporizedWater()) {
851 const auto& Rvw = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
853 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
854 if constexpr (blackoilConserveSurfaceVolume) {
855 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
858 flux[conti0EqIdx + activeWaterCompIdx] +=
859 Rvw * surfaceVolumeFlux *
860 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
877 template <
class Scalar>
879 unsigned pvtRegionIdx)
881 if constexpr (!blackoilConserveSurfaceVolume) {
886 if constexpr (waterEnabled) {
887 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
888 container[conti0EqIdx + activeWaterCompIdx] *=
889 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
892 if constexpr (gasEnabled) {
893 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
894 container[conti0EqIdx + activeGasCompIdx] *=
895 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
898 if constexpr (oilEnabled) {
899 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
900 container[conti0EqIdx + activeOilCompIdx] *=
901 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
908 {
return dirId < 0 ? FaceDir::DirEnum::Unknown : FaceDir::FromIntersectionIndex(dirId); }
Contains the classes required to extend the black-oil model by bioeffects.
Contains the classes required to extend the black-oil model by brine.
Classes required for dynamic convective mixing.
Classes required for molecular diffusion.
Classes required for mechanical dispersion.
Contains the classes required to extend the black-oil model by energy.
Contains the classes required to extend the black-oil model by solvent component. For details,...
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by polymer.
Declares the properties required by the black oil model.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:93
static void applyScaling(RateVector &flux)
Definition: blackoilbioeffectsmodules.hh:235
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbioeffectsmodules.hh:177
static void addSource(RateVector &source, const Problem &problem, const IntensiveQuantities &intQuants, unsigned globalSpaceIdex)
Definition: blackoilbioeffectsmodules.hh:276
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:56
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:163
Definition: blackoilconvectivemixingmodule.hh:103
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: blackoildiffusionmodule.hh:50
Provides the auxiliary methods required for consideration of the dispersion equation.
Definition: blackoildispersionmodule.hh:58
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:61
static void addHeatFlux(RateVector &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:214
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:156
Contains the high level supplements required to extend the black oil model.
Definition: blackoilextbomodules.hh:62
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilextbomodules.hh:156
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition: blackoilfoammodules.hh:58
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilfoammodules.hh:164
Calculates the local residual of the black oil model.
Definition: blackoillocalresidualtpfa.hh:64
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidualtpfa.hh:750
static void computeBoundaryFlux(RateVector &bdyFlux, const Problem &problem, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:531
static void computeBoundaryFluxFree(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:566
static void computeBoundaryFluxRate(RateVector &bdyFlux, const BoundaryConditionData &bdyInfo)
Definition: blackoillocalresidualtpfa.hh:559
static FaceDir::DirEnum faceDirFromDirId(const int dirId)
Definition: blackoillocalresidualtpfa.hh:907
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Definition: blackoillocalresidualtpfa.hh:769
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:267
static void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoillocalresidualtpfa.hh:158
static void computeFlux(RateVector &flux, RateVector &darcy, const unsigned globalIndexIn, const unsigned globalIndexEx, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const ResidualNBInfo &nbInfo, const ModuleParams &moduleParams)
Definition: blackoillocalresidualtpfa.hh:241
static void computeBoundaryThermal(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:675
static void adaptMassConservationQuantities_(Dune::FieldVector< Scalar, numEq > &container, unsigned pvtRegionIdx)
Helper function to convert the mass-related parts of a Dune::FieldVector that stores conservation qua...
Definition: blackoillocalresidualtpfa.hh:878
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const Eval &surfaceVolumeFlux, const FluidState &upFs)
Helper function to calculate the flux of mass in terms of conservation quantities via specific fluid ...
Definition: blackoillocalresidualtpfa.hh:785
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: blackoillocalresidualtpfa.hh:148
static void calculateFluxes_(RateVector &flux, RateVector &darcy, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned &globalIndexIn, const unsigned &globalIndexEx, const ResidualNBInfo &nbInfo, const ModuleParams &moduleParams)
Definition: blackoillocalresidualtpfa.hh:335
static void computeSource(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:710
static void computeSourceDense(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:729
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilpolymermodules.hh:236
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:68
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:183
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:43
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: blackoillocalresidualtpfa.hh:130
FaceDir::DirEnum faceDir
Definition: blackoillocalresidualtpfa.hh:135
double faceArea
Definition: blackoillocalresidualtpfa.hh:132
ConditionalStorage< enableFullyImplicitThermal, double > outAlpha
Definition: blackoillocalresidualtpfa.hh:139
ConditionalStorage< enableFullyImplicitThermal, double > inAlpha
Definition: blackoillocalresidualtpfa.hh:138
ConditionalStorage< enableDiffusion, double > diffusivity
Definition: blackoillocalresidualtpfa.hh:140
double dZg
Definition: blackoillocalresidualtpfa.hh:134
double Vin
Definition: blackoillocalresidualtpfa.hh:136
ConditionalStorage< enableDispersion, double > dispersivity
Definition: blackoillocalresidualtpfa.hh:141
double thpres
Definition: blackoillocalresidualtpfa.hh:133
double trans
Definition: blackoillocalresidualtpfa.hh:131
double Vex
Definition: blackoillocalresidualtpfa.hh:137
ConvectiveMixingModuleParamT convectiveMixingModuleParam
Definition: blackoilmoduleparams.hh:23