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>
61template <
class TypeTag>
74 using FluidState =
typename IntensiveQuantities::FluidState;
76 enum { conti0EqIdx = Indices::conti0EqIdx };
77 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
78 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
79 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
81 enum { dimWorld = GridView::dimensionworld };
82 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
83 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
84 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
86 enum { gasCompIdx = FluidSystem::gasCompIdx };
87 enum { oilCompIdx = FluidSystem::oilCompIdx };
88 enum { waterCompIdx = FluidSystem::waterCompIdx };
89 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
91 static constexpr bool waterEnabled = Indices::waterEnabled;
92 static constexpr bool gasEnabled = Indices::gasEnabled;
93 static constexpr bool oilEnabled = Indices::oilEnabled;
94 static constexpr bool compositionSwitchEnabled = compositionSwitchIdx >= 0;
96 static constexpr bool blackoilConserveSurfaceVolume =
97 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
99 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
100 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
101 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
102 static constexpr bool enableEnergy = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
103 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
104 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
105 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
106 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
107 static constexpr bool enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
108 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
109 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
110 static constexpr bool enableMICP = Indices::enableMICP;
120 using ConvectiveMixingModuleParam =
typename ConvectiveMixingModule::ConvectiveMixingModuleParam;
125 using Toolbox = MathToolbox<Evaluation>;
137 ConditionalStorage<enableEnergy, double>
inAlpha;
151 template <
class LhsEval>
153 const ElementContext& elemCtx,
155 unsigned timeIdx)
const
157 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
161 template <
class LhsEval>
163 const IntensiveQuantities& intQuants)
167 const auto& fs = intQuants.fluidState();
170 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
171 if (!FluidSystem::phaseIsActive(phaseIdx)) {
174 unsigned activeCompIdx =
175 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
176 LhsEval surfaceVolume =
177 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
178 Toolbox::template decay<LhsEval>(fs.invB(phaseIdx)) *
179 Toolbox::template decay<LhsEval>(intQuants.porosity());
181 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
184 if (phaseIdx == oilPhaseIdx && FluidSystem::enableDissolvedGas()) {
185 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
186 storage[conti0EqIdx + activeGasCompIdx] +=
187 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs()) *
192 if (phaseIdx == waterPhaseIdx && FluidSystem::enableDissolvedGasInWater()) {
193 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
194 storage[conti0EqIdx + activeGasCompIdx] +=
195 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw()) *
200 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedOil()) {
201 unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
202 storage[conti0EqIdx + activeOilCompIdx] +=
203 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv()) *
208 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedWater()) {
209 unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
210 storage[conti0EqIdx + activeWaterCompIdx] +=
211 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw()) *
247 const unsigned globalIndexIn,
248 const unsigned globalIndexEx,
249 const IntensiveQuantities& intQuantsIn,
250 const IntensiveQuantities& intQuantsEx,
254 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
272 const ElementContext& elemCtx,
276 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
277 assert(timeIdx == 0);
280 RateVector darcy = 0.0;
282 const auto& problem = elemCtx.problem();
283 const auto& stencil = elemCtx.stencil(timeIdx);
284 const auto& scvf = stencil.interiorFace(scvfIdx);
286 unsigned interiorDofIdx = scvf.interiorIndex();
287 unsigned exteriorDofIdx = scvf.exteriorIndex();
288 assert(interiorDofIdx != exteriorDofIdx);
292 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, 0);
293 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, 0);
294 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
295 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
296 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
297 Scalar faceArea = scvf.area();
299 Scalar thpres = problem.thresholdPressure(globalIndexIn, globalIndexEx);
304 const Scalar g = problem.gravity()[dimWorld - 1];
305 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
306 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
313 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
314 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
317 const Scalar distZ = zIn - zEx;
319 const Scalar inAlpha = problem.thermalHalfTransmissibility(globalIndexIn, globalIndexEx);
320 const Scalar outAlpha = problem.thermalHalfTransmissibility(globalIndexEx, globalIndexIn);
321 const Scalar diffusivity = problem.diffusivity(globalIndexEx, globalIndexIn);
322 const Scalar dispersivity = problem.dispersivity(globalIndexEx, globalIndexIn);
325 trans, faceArea, thpres, distZ * g, faceDir, Vin, Vex,
326 inAlpha, outAlpha, diffusivity, dispersivity
336 problem.moduleParams());
341 const IntensiveQuantities& intQuantsIn,
342 const IntensiveQuantities& intQuantsEx,
343 const unsigned& globalIndexIn,
344 const unsigned& globalIndexEx,
348 OPM_TIMEBLOCK_LOCAL(calculateFluxes, Subsystem::Assembly);
349 const Scalar Vin = nbInfo.
Vin;
350 const Scalar Vex = nbInfo.
Vex;
351 const Scalar distZg = nbInfo.
dZg;
352 const Scalar thpres = nbInfo.
thpres;
353 const Scalar trans = nbInfo.
trans;
354 const Scalar faceArea = nbInfo.
faceArea;
355 FaceDir::DirEnum facedir = nbInfo.
faceDir;
357 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
358 if (!FluidSystem::phaseIsActive(phaseIdx)) {
366 short interiorDofIdx = 0;
367 short exteriorDofIdx = 1;
368 Evaluation pressureDifference;
369 ExtensiveQuantities::calculatePhasePressureDiff_(upIdx,
385 const IntensiveQuantities& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
386 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
388 Evaluation transMult = (intQuantsIn.rockCompTransMultiplier() +
389 Toolbox::value(intQuantsEx.rockCompTransMultiplier())) / 2;
390 if constexpr (enableBioeffects || enableSaltPrecipitation) {
391 transMult *= (intQuantsIn.permFactor() + Toolbox::value(intQuantsEx.permFactor())) / 2;
393 Evaluation darcyFlux;
394 if (globalUpIndex == globalIndexIn) {
395 darcyFlux = pressureDifference * up.mobility(phaseIdx, facedir) * transMult * (-trans / faceArea);
397 darcyFlux = pressureDifference *
398 (Toolbox::value(up.mobility(phaseIdx, facedir)) * transMult * (-trans / faceArea));
401 unsigned activeCompIdx =
402 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
404 darcy[conti0EqIdx + activeCompIdx] = darcyFlux.value() * faceArea;
406 unsigned pvtRegionIdx = up.pvtRegionIndex();
408 if (globalUpIndex == globalIndexIn) {
410 = getInvB_<FluidSystem, FluidState, Evaluation>(up.fluidState(), phaseIdx, pvtRegionIdx);
411 const auto& surfaceVolumeFlux = invB * darcyFlux;
412 evalPhaseFluxes_<Evaluation>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
413 if constexpr (enableEnergy) {
414 EnergyModule::template
415 addPhaseEnthalpyFluxes_<Evaluation>(flux, phaseIdx, darcyFlux, up.fluidState());
417 if constexpr (enableBioeffects) {
418 BioeffectsModule::template
419 addBioeffectsFluxes_<Evaluation>(flux, phaseIdx, darcyFlux, up);
421 if constexpr (enableBrine) {
422 BrineModule::template
423 addBrineFluxes_<Evaluation, FluidState>(flux, phaseIdx, darcyFlux, up.fluidState());
426 const auto& invB = getInvB_<FluidSystem, FluidState, Scalar>(up.fluidState(), phaseIdx, pvtRegionIdx);
427 const auto& surfaceVolumeFlux = invB * darcyFlux;
428 evalPhaseFluxes_<Scalar>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
429 if constexpr (enableEnergy) {
430 EnergyModule::template
431 addPhaseEnthalpyFluxes_<Scalar>(flux, phaseIdx, darcyFlux, up.fluidState());
433 if constexpr (enableBioeffects) {
434 BioeffectsModule::template
435 addBioeffectsFluxes_<Scalar>(flux, phaseIdx, darcyFlux, up);
437 if constexpr (enableBrine) {
438 BrineModule::template
439 addBrineFluxes_<Scalar, FluidState>(flux, phaseIdx, darcyFlux, up.fluidState());
445 static_assert(!enableSolvent,
446 "Relevant computeFlux() method must be implemented for this module before enabling.");
450 static_assert(!enableExtbo,
451 "Relevant computeFlux() method must be implemented for this module before enabling.");
455 static_assert(!enablePolymer,
456 "Relevant computeFlux() method must be implemented for this module before enabling.");
460 if constexpr (enableConvectiveMixing) {
461 ConvectiveMixingModule::addConvectiveMixingFlux(flux,
473 if constexpr (enableEnergy) {
474 const Scalar inAlpha = nbInfo.
inAlpha;
475 const Scalar outAlpha = nbInfo.
outAlpha;
478 short interiorDofIdx = 0;
479 short exteriorDofIdx = 1;
481 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
487 intQuantsIn.fluidState(),
488 intQuantsEx.fluidState(),
498 static_assert(!enableFoam,
499 "Relevant computeFlux() method must be implemented for this module before enabling.");
503 if constexpr (enableDiffusion) {
504 typename DiffusionModule::ExtensiveQuantities::EvaluationArray effectiveDiffusionCoefficient;
505 DiffusionModule::ExtensiveQuantities::update(effectiveDiffusionCoefficient, intQuantsIn, intQuantsEx);
507 const Scalar tmpdiffusivity = diffusivity / faceArea;
508 DiffusionModule::addDiffusiveFlux(flux,
512 effectiveDiffusionCoefficient);
516 if constexpr (enableDispersion) {
517 typename DispersionModule::ExtensiveQuantities::ScalarArray normVelocityAvg;
518 DispersionModule::ExtensiveQuantities::update(normVelocityAvg, intQuantsIn, intQuantsEx);
520 const Scalar tmpdispersivity = dispersivity / faceArea;
521 DispersionModule::addDispersiveFlux(flux,
529 if constexpr (enableMICP) {
534 template <
class BoundaryConditionData>
536 const Problem& problem,
537 const BoundaryConditionData& bdyInfo,
538 const IntensiveQuantities& insideIntQuants,
539 unsigned globalSpaceIdx)
541 switch (bdyInfo.type) {
549 case BCType::DIRICHLET:
552 case BCType::THERMAL:
556 throw std::logic_error(
"Unknown boundary condition type " +
558 " in computeBoundaryFlux()." );
562 template <
class BoundaryConditionData>
564 const BoundaryConditionData& bdyInfo)
566 bdyFlux.setMassRate(bdyInfo.massRate, bdyInfo.pvtRegionIdx);
569 template <
class BoundaryConditionData>
572 const BoundaryConditionData& bdyInfo,
573 const IntensiveQuantities& insideIntQuants,
574 unsigned globalSpaceIdx)
577 std::array<short, numPhases> upIdx;
578 std::array<short, numPhases> dnIdx;
579 std::array<Evaluation, numPhases> volumeFlux;
580 std::array<Evaluation, numPhases> pressureDifference;
582 ExtensiveQuantities::calculateBoundaryGradients_(problem,
585 bdyInfo.boundaryFaceIndex,
588 bdyInfo.exFluidState,
598 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
599 if (!FluidSystem::phaseIsActive(phaseIdx)) {
602 const auto& pBoundary = bdyInfo.exFluidState.pressure(phaseIdx);
603 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
604 const unsigned pvtRegionIdx = insideIntQuants.pvtRegionIndex();
607 const auto& darcyFlux = volumeFlux[phaseIdx];
609 if (pBoundary < pInside) {
612 getInvB_<FluidSystem, FluidState, Evaluation>(insideIntQuants.fluidState(), phaseIdx, pvtRegionIdx);
613 Evaluation surfaceVolumeFlux = invB * darcyFlux;
614 evalPhaseFluxes_<Evaluation>(tmp,
616 insideIntQuants.pvtRegionIndex(),
618 insideIntQuants.fluidState());
619 if constexpr (enableEnergy) {
620 EnergyModule::template
621 addPhaseEnthalpyFluxes_<Evaluation>(tmp, phaseIdx, darcyFlux, insideIntQuants.fluidState());
623 }
else if (pBoundary > pInside) {
625 using ScalarFluidState =
decltype(bdyInfo.exFluidState);
627 getInvB_<FluidSystem, ScalarFluidState, Scalar>(bdyInfo.exFluidState, phaseIdx, pvtRegionIdx);
628 Evaluation surfaceVolumeFlux = invB * darcyFlux;
629 evalPhaseFluxes_<Scalar>(tmp,
631 insideIntQuants.pvtRegionIndex(),
633 bdyInfo.exFluidState);
634 if constexpr (enableEnergy) {
635 EnergyModule::template
636 addPhaseEnthalpyFluxes_<Scalar>(tmp, phaseIdx, darcyFlux, bdyInfo.exFluidState);
640 for (
unsigned i = 0; i < tmp.size(); ++i) {
641 bdyFlux[i] += tmp[i];
646 if constexpr (enableEnergy) {
650 problem.eclTransmissibilities().thermalHalfTransBoundary(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
653 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
658 bdyInfo.exFluidState);
662 static_assert(!enableSolvent,
663 "Relevant treatment of boundary conditions must be implemented before enabling.");
664 static_assert(!enablePolymer,
665 "Relevant treatment of boundary conditions must be implemented before enabling.");
671 for (
unsigned i = 0; i < numEq; ++i) {
672 Valgrind::CheckDefined(bdyFlux[i]);
674 Valgrind::CheckDefined(bdyFlux);
678 template <
class BoundaryConditionData>
681 const BoundaryConditionData& bdyInfo,
682 const IntensiveQuantities& insideIntQuants,
683 [[maybe_unused]]
unsigned globalSpaceIdx)
690 if constexpr (enableEnergy) {
694 problem.eclTransmissibilities().thermalHalfTransBoundary(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
697 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
702 bdyInfo.exFluidState);
707 for (
unsigned i = 0; i < numEq; ++i) {
708 Valgrind::CheckDefined(bdyFlux[i]);
710 Valgrind::CheckDefined(bdyFlux);
715 const Problem& problem,
716 const IntensiveQuantities& insideIntQuants,
717 unsigned globalSpaceIdex,
722 problem.source(source, globalSpaceIdex, timeIdx);
728 if constexpr (enableEnergy) {
729 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
734 const Problem& problem,
735 const IntensiveQuantities& insideIntQuants,
736 unsigned globalSpaceIdex,
740 problem.addToSourceDense(source, globalSpaceIdex, timeIdx);
746 if constexpr (enableEnergy) {
747 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
755 const ElementContext& elemCtx,
757 unsigned timeIdx)
const
761 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
767 if constexpr (enableEnergy) {
768 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
772 template <
class UpEval,
class Flu
idState>
775 unsigned pvtRegionIdx,
776 const ExtensiveQuantities& extQuants,
777 const FluidState& upFs)
779 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
780 const auto& surfaceVolumeFlux = invB * extQuants.volumeFlux(phaseIdx);
781 evalPhaseFluxes_<UpEval>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, upFs);
788 template <
class UpEval,
class Eval,
class Flu
idState>
791 unsigned pvtRegionIdx,
792 const Eval& surfaceVolumeFlux,
793 const FluidState& upFs)
795 unsigned activeCompIdx =
796 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
798 if constexpr (blackoilConserveSurfaceVolume) {
799 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
802 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux *
803 FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
806 if (phaseIdx == oilPhaseIdx) {
808 if (FluidSystem::enableDissolvedGas()) {
809 const auto& Rs = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
811 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
812 if constexpr (blackoilConserveSurfaceVolume) {
813 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
816 flux[conti0EqIdx + activeGasCompIdx] +=
817 Rs * surfaceVolumeFlux *
818 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
822 else if (phaseIdx == waterPhaseIdx) {
824 if (FluidSystem::enableDissolvedGasInWater()) {
825 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
827 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
828 if constexpr (blackoilConserveSurfaceVolume) {
829 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
832 flux[conti0EqIdx + activeGasCompIdx] +=
833 Rsw * surfaceVolumeFlux *
834 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
838 else if (phaseIdx == gasPhaseIdx) {
840 if (FluidSystem::enableVaporizedOil()) {
841 const auto& Rv = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
843 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
844 if constexpr (blackoilConserveSurfaceVolume) {
845 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
848 flux[conti0EqIdx + activeOilCompIdx] +=
849 Rv * surfaceVolumeFlux *
850 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
854 if (FluidSystem::enableVaporizedWater()) {
855 const auto& Rvw = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
857 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
858 if constexpr (blackoilConserveSurfaceVolume) {
859 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
862 flux[conti0EqIdx + activeWaterCompIdx] +=
863 Rvw * surfaceVolumeFlux *
864 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
881 template <
class Scalar>
883 unsigned pvtRegionIdx)
885 if constexpr (!blackoilConserveSurfaceVolume) {
890 if constexpr (waterEnabled) {
891 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
892 container[conti0EqIdx + activeWaterCompIdx] *=
893 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
896 if constexpr (gasEnabled) {
897 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
898 container[conti0EqIdx + activeGasCompIdx] *=
899 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
902 if constexpr (oilEnabled) {
903 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
904 container[conti0EqIdx + activeOilCompIdx] *=
905 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
912 {
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:65
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:60
static void addHeatFlux(RateVector &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:213
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:155
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:63
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidualtpfa.hh:754
static void computeBoundaryFlux(RateVector &bdyFlux, const Problem &problem, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:535
static void computeBoundaryFluxFree(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:570
static void computeBoundaryFluxRate(RateVector &bdyFlux, const BoundaryConditionData &bdyInfo)
Definition: blackoillocalresidualtpfa.hh:563
static FaceDir::DirEnum faceDirFromDirId(const int dirId)
Definition: blackoillocalresidualtpfa.hh:911
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Definition: blackoillocalresidualtpfa.hh:773
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:271
static void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoillocalresidualtpfa.hh:162
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:245
static void computeBoundaryThermal(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:679
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:882
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:789
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:152
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:339
static void computeSource(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:714
static void computeSourceDense(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:733
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:144
ConvectiveMixingModuleParam convectiveMixingModuleParam
Definition: blackoillocalresidualtpfa.hh:145
Definition: blackoillocalresidualtpfa.hh:129
FaceDir::DirEnum faceDir
Definition: blackoillocalresidualtpfa.hh:134
double faceArea
Definition: blackoillocalresidualtpfa.hh:131
ConditionalStorage< enableEnergy, double > inAlpha
Definition: blackoillocalresidualtpfa.hh:137
ConditionalStorage< enableDiffusion, double > diffusivity
Definition: blackoillocalresidualtpfa.hh:139
double dZg
Definition: blackoillocalresidualtpfa.hh:133
double Vin
Definition: blackoillocalresidualtpfa.hh:135
ConditionalStorage< enableDispersion, double > dispersivity
Definition: blackoillocalresidualtpfa.hh:140
double thpres
Definition: blackoillocalresidualtpfa.hh:132
ConditionalStorage< enableEnergy, double > outAlpha
Definition: blackoillocalresidualtpfa.hh:138
double trans
Definition: blackoillocalresidualtpfa.hh:130
double Vex
Definition: blackoillocalresidualtpfa.hh:136