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/ConditionalStorage.hpp>
35#include <opm/material/common/MathToolbox.hpp>
36#include <opm/material/fluidstates/BlackOilFluidState.hpp>
51#include <opm/common/ErrorMacros.hpp>
52#include <opm/common/utility/gpuDecorators.hpp>
59#include <opm/common/utility/gpuistl_if_available.hpp>
69template <
class TypeTag>
82 using FluidState =
typename IntensiveQuantities::FluidState;
84 enum { conti0EqIdx = Indices::conti0EqIdx };
85 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
86 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
87 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
89 enum { dimWorld = GridView::dimensionworld };
90 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
91 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
92 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
94 enum { gasCompIdx = FluidSystem::gasCompIdx };
95 enum { oilCompIdx = FluidSystem::oilCompIdx };
96 enum { waterCompIdx = FluidSystem::waterCompIdx };
97 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
99 static constexpr bool waterEnabled = Indices::waterEnabled;
100 static constexpr bool gasEnabled = Indices::gasEnabled;
101 static constexpr bool oilEnabled = Indices::oilEnabled;
102 static constexpr bool compositionSwitchEnabled = compositionSwitchIdx >= 0;
104 static constexpr bool blackoilConserveSurfaceVolume
105 = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
107 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
108 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
109 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
110 static constexpr bool enableFullyImplicitThermal
111 = (getPropValue<TypeTag, Properties::EnergyModuleType>()
112 == EnergyModules::FullyImplicitThermal);
113 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
114 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
115 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
116 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
117 static constexpr bool enableConvectiveMixing
118 = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
119 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
120 static constexpr bool enableSaltPrecipitation
121 = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
122 static constexpr bool enableMICP = Indices::enableMICP;
123 static constexpr bool runAssemblyOnGpu = getPropValue<TypeTag, Properties::RunAssemblyOnGpu>();
138 using Toolbox = MathToolbox<Evaluation>;
149 ConditionalStorage<enableFullyImplicitThermal, double>
inAlpha;
150 ConditionalStorage<enableFullyImplicitThermal, double>
outAlpha;
158 template <
class LhsEval>
160 const ElementContext& elemCtx,
162 unsigned timeIdx)
const
164 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
165 computeStorage<LhsEval>(storage, intQuants);
168 template <
class LhsEval,
class StorageType,
class IntensiveQuantitiesType = IntensiveQuantities>
170 const IntensiveQuantitiesType& intQuants)
174 const auto& fs = intQuants.fluidState();
177 const FluidSystem& fsys = intQuants.getFluidSystem();
179 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
180 if (!fsys.phaseIsActive(phaseIdx)) {
183 unsigned activeCompIdx
184 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
185 LhsEval surfaceVolume = Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
186 * Toolbox::template decay<LhsEval>(fs.invB(phaseIdx))
187 * Toolbox::template decay<LhsEval>(intQuants.porosity());
189 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
192 if (phaseIdx == oilPhaseIdx && fsys.enableDissolvedGas()) {
193 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
194 storage[conti0EqIdx + activeGasCompIdx]
195 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs())
200 if (phaseIdx == waterPhaseIdx && fsys.enableDissolvedGasInWater()) {
201 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
202 storage[conti0EqIdx + activeGasCompIdx]
203 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw())
208 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedOil()) {
209 unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
210 storage[conti0EqIdx + activeOilCompIdx]
211 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv())
216 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedWater()) {
217 unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
218 storage[conti0EqIdx + activeWaterCompIdx]
219 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw())
253 template <
class ModuleParamsT,
255 class IntensiveQuantitiesT,
256 class ResidualNBInfoT>
259 const unsigned globalIndexIn,
260 const unsigned globalIndexEx,
261 const IntensiveQuantitiesT& intQuantsIn,
262 const IntensiveQuantitiesT& intQuantsEx,
263 const ResidualNBInfoT& nbInfo,
264 const ModuleParamsT& moduleParams)
266 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
284 computeFlux(RateVector& flux,
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
286 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
287 assert(timeIdx == 0);
290 RateVector darcy = 0.0;
292 const auto& problem = elemCtx.problem();
293 const auto& stencil = elemCtx.stencil(timeIdx);
294 const auto& scvf = stencil.interiorFace(scvfIdx);
296 unsigned interiorDofIdx = scvf.interiorIndex();
297 unsigned exteriorDofIdx = scvf.exteriorIndex();
298 assert(interiorDofIdx != exteriorDofIdx);
302 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, 0);
303 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, 0);
304 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
305 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
306 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
307 Scalar faceArea = scvf.area();
309 Scalar thpres = problem.thresholdPressure(globalIndexIn, globalIndexEx);
314 const Scalar g = problem.gravity()[dimWorld - 1];
315 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
316 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
323 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
324 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
327 const Scalar distZ = zIn - zEx;
329 const Scalar inAlpha = problem.thermalHalfTransmissibility(globalIndexIn, globalIndexEx);
330 const Scalar outAlpha = problem.thermalHalfTransmissibility(globalIndexEx, globalIndexIn);
331 const Scalar diffusivity = problem.diffusivity(globalIndexEx, globalIndexIn);
332 const Scalar dispersivity = problem.dispersivity(globalIndexEx, globalIndexIn);
353 problem.moduleParams());
356 template <
class RateVectorT,
357 class IntensiveQuantitiesT,
358 class ResidualNBInfoT,
362 const IntensiveQuantitiesT& intQuantsIn,
363 const IntensiveQuantitiesT& intQuantsEx,
364 const unsigned& globalIndexIn,
365 const unsigned& globalIndexEx,
366 const ResidualNBInfoT& nbInfo,
367 const ModuleParamsT& moduleParams)
369 OPM_TIMEBLOCK_LOCAL(calculateFluxes, Subsystem::Assembly);
370 const Scalar Vin = nbInfo.Vin;
371 const Scalar Vex = nbInfo.Vex;
372 const Scalar distZg = nbInfo.dZg;
373 const Scalar thpres = nbInfo.thpres;
374 const Scalar trans = nbInfo.trans;
375 const Scalar faceArea = nbInfo.faceArea;
376 FaceDir::DirEnum facedir = nbInfo.faceDir;
378 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
380 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
381 if (!fsys.phaseIsActive(phaseIdx)) {
389 short interiorDofIdx = 0;
390 short exteriorDofIdx = 1;
391 Evaluation pressureDifference;
392 ExtensiveQuantities::calculatePhasePressureDiff_(upIdx,
408 const IntensiveQuantities& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
409 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
412 Evaluation transMult = (intQuantsIn.rockCompTransMultiplier()
413 + Toolbox::value(intQuantsEx.rockCompTransMultiplier()))
415 if constexpr (enableBioeffects || enableSaltPrecipitation) {
417 *= (intQuantsIn.permFactor() + Toolbox::value(intQuantsEx.permFactor())) / 2;
420 Evaluation darcyFlux;
421 if (globalUpIndex == globalIndexIn) {
422 darcyFlux = pressureDifference * up.mobility(phaseIdx, facedir) * transMult
423 * (-trans / faceArea);
425 darcyFlux = pressureDifference
426 * (Toolbox::value(up.mobility(phaseIdx, facedir)) * transMult
427 * (-trans / faceArea));
430 unsigned activeCompIdx
431 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
433 darcy[conti0EqIdx + activeCompIdx] = darcyFlux.value() * faceArea;
435 unsigned pvtRegionIdx = up.pvtRegionIndex();
437 if (globalUpIndex == globalIndexIn) {
438 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
439 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
440 const auto& surfaceVolumeFlux = invB * darcyFlux;
442 evalPhaseFluxes_<Evaluation>(
443 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
444 if constexpr (enableFullyImplicitThermal) {
445 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
446 flux, phaseIdx, darcyFlux, up.fluidState());
448 if constexpr (enableBioeffects) {
449 BioeffectsModule::template addBioeffectsFluxes_<Evaluation>(
450 flux, phaseIdx, darcyFlux, up);
452 if constexpr (enableBrine) {
453 BrineModule::template addBrineFluxes_<Evaluation, FluidState>(
454 flux, phaseIdx, darcyFlux, up.fluidState());
457 const auto& invB = getInvB_<FluidSystem, FluidState, Scalar>(
458 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
459 const auto& surfaceVolumeFlux = invB * darcyFlux;
460 evalPhaseFluxes_<Scalar>(
461 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
462 if constexpr (enableFullyImplicitThermal) {
463 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
464 flux, phaseIdx, darcyFlux, up.fluidState());
466 if constexpr (enableBioeffects) {
467 BioeffectsModule::template addBioeffectsFluxes_<Scalar>(
468 flux, phaseIdx, darcyFlux, up);
470 if constexpr (enableBrine) {
471 BrineModule::template addBrineFluxes_<Scalar, FluidState>(
472 flux, phaseIdx, darcyFlux, up.fluidState());
480 "Relevant computeFlux() method must be implemented for this module before enabling.");
486 "Relevant computeFlux() method must be implemented for this module before enabling.");
492 "Relevant computeFlux() method must be implemented for this module before enabling.");
496 if constexpr (enableConvectiveMixing) {
497 ConvectiveMixingModule::addConvectiveMixingFlux(
506 moduleParams.convectiveMixingModuleParam);
510 if constexpr (enableFullyImplicitThermal) {
511 const Scalar inAlpha = nbInfo.inAlpha;
512 const Scalar outAlpha = nbInfo.outAlpha;
515 short interiorDofIdx = 0;
516 short exteriorDofIdx = 1;
518 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
524 intQuantsIn.fluidState(),
525 intQuantsEx.fluidState(),
537 "Relevant computeFlux() method must be implemented for this module before enabling.");
542 if constexpr (enableDiffusion) {
543 typename DiffusionModule::ExtensiveQuantities::EvaluationArray
544 effectiveDiffusionCoefficient;
545 DiffusionModule::ExtensiveQuantities::update(
546 effectiveDiffusionCoefficient, intQuantsIn, intQuantsEx);
547 const Scalar diffusivity = nbInfo.diffusivity;
548 const Scalar tmpdiffusivity = diffusivity / faceArea;
549 DiffusionModule::addDiffusiveFlux(
550 flux, intQuantsIn, intQuantsEx, tmpdiffusivity, effectiveDiffusionCoefficient);
555 if constexpr (enableDispersion) {
556 typename DispersionModule::ExtensiveQuantities::ScalarArray normVelocityAvg;
557 DispersionModule::ExtensiveQuantities::update(
558 normVelocityAvg, intQuantsIn, intQuantsEx);
559 const Scalar dispersivity = nbInfo.dispersivity;
560 const Scalar tmpdispersivity = dispersivity / faceArea;
561 DispersionModule::addDispersiveFlux(
562 flux, intQuantsIn, intQuantsEx, tmpdispersivity, normVelocityAvg);
566 if constexpr (enableMICP) {
571 template <
class BoundaryConditionData,
class RateVectorLocal,
class LocalProblem>
573 const LocalProblem& problem,
574 const BoundaryConditionData& bdyInfo,
575 const IntensiveQuantities& insideIntQuants,
576 unsigned globalSpaceIdx)
578#if OPM_IS_INSIDE_HOST_FUNCTION
579 switch (bdyInfo.type) {
587 case BCType::DIRICHLET:
590 case BCType::THERMAL:
594 throw std::logic_error(
"Unknown boundary condition type "
596 +
" in computeBoundaryFlux().");
599 switch (bdyInfo.type) {
603 case BCType::THERMAL:
607 OPM_THROW(std::logic_error,
608 "Boundary condition type " +
std::to_string(
static_cast<int>(bdyInfo.type))
609 +
" is not supported for GPU fluid systems in computeBoundaryFlux().");
614 template <
class BoundaryConditionData>
617 bdyFlux.setMassRate(bdyInfo.massRate, bdyInfo.pvtRegionIdx);
620 template <
class BoundaryConditionData>
623 const BoundaryConditionData& bdyInfo,
624 const IntensiveQuantities& insideIntQuants,
625 unsigned globalSpaceIdx)
628 std::array<short, numPhases> upIdx;
629 std::array<short, numPhases> dnIdx;
630 std::array<Evaluation, numPhases> volumeFlux;
631 std::array<Evaluation, numPhases> pressureDifference;
633 ExtensiveQuantities::calculateBoundaryGradients_(problem,
636 bdyInfo.boundaryFaceIndex,
639 bdyInfo.exFluidState,
649 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
650 if (!FluidSystem::phaseIsActive(phaseIdx)) {
653 const auto& pBoundary = bdyInfo.exFluidState.pressure(phaseIdx);
654 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
655 const unsigned pvtRegionIdx = insideIntQuants.pvtRegionIndex();
658 const auto& darcyFlux = volumeFlux[phaseIdx];
660 if (pBoundary < pInside) {
662 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
663 insideIntQuants.fluidState(), phaseIdx, pvtRegionIdx);
664 Evaluation surfaceVolumeFlux = invB * darcyFlux;
665 evalPhaseFluxes_<Evaluation>(tmp,
667 insideIntQuants.pvtRegionIndex(),
669 insideIntQuants.fluidState());
670 if constexpr (enableFullyImplicitThermal) {
671 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
672 tmp, phaseIdx, darcyFlux, insideIntQuants.fluidState());
674 }
else if (pBoundary > pInside) {
676 using ScalarFluidState =
decltype(bdyInfo.exFluidState);
677 const auto& invB = getInvB_<FluidSystem, ScalarFluidState, Scalar>(
678 bdyInfo.exFluidState, phaseIdx, pvtRegionIdx);
679 Evaluation surfaceVolumeFlux = invB * darcyFlux;
680 evalPhaseFluxes_<Scalar>(tmp,
682 insideIntQuants.pvtRegionIndex(),
684 bdyInfo.exFluidState);
685 if constexpr (enableFullyImplicitThermal) {
686 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
687 tmp, phaseIdx, darcyFlux, bdyInfo.exFluidState);
691 for (
unsigned i = 0; i < tmp.size(); ++i) {
692 bdyFlux[i] += tmp[i];
697 if constexpr (enableFullyImplicitThermal) {
700 Scalar alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
701 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
704 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
709 bdyInfo.exFluidState);
715 "Relevant treatment of boundary conditions must be implemented before enabling.");
718 "Relevant treatment of boundary conditions must be implemented before enabling.");
720 const FluidSystem& fsys = insideIntQuants.getFluidSystem();
726 for (
unsigned i = 0; i < numEq; ++i) {
727 Valgrind::CheckDefined(bdyFlux[i]);
729 Valgrind::CheckDefined(bdyFlux);
733 template <
class ProblemLocal,
class BoundaryConditionData,
class RateVectorLocal>
735 RateVectorLocal& bdyFlux,
736 const BoundaryConditionData& bdyInfo,
737 const IntensiveQuantities& insideIntQuants,
738 [[maybe_unused]]
unsigned globalSpaceIdx)
745 if constexpr (enableFullyImplicitThermal) {
750 if constexpr (runAssemblyOnGpu) {
753 alpha = problem.getAlpha(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
755 alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
756 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
761 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
766 bdyInfo.exFluidState);
771 for (
unsigned i = 0; i < numEq; ++i) {
772 Valgrind::CheckDefined(bdyFlux[i]);
774 Valgrind::CheckDefined(bdyFlux);
779 const Problem& problem,
780 const IntensiveQuantities& insideIntQuants,
781 unsigned globalSpaceIdex,
786 problem.source(source, globalSpaceIdex, timeIdx);
792 if constexpr (enableFullyImplicitThermal) {
793 source[Indices::contiEnergyEqIdx]
794 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
799 const Problem& problem,
800 const IntensiveQuantities& insideIntQuants,
801 unsigned globalSpaceIdex,
805 problem.addToSourceDense(source, globalSpaceIdex, timeIdx);
811 if constexpr (enableFullyImplicitThermal) {
812 source[Indices::contiEnergyEqIdx]
813 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
821 const ElementContext& elemCtx,
823 unsigned timeIdx)
const
827 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
833 if constexpr (enableFullyImplicitThermal) {
834 source[Indices::contiEnergyEqIdx]
835 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
839 template <
class UpEval,
class Flu
idState>
842 unsigned pvtRegionIdx,
843 const ExtensiveQuantities& extQuants,
844 const FluidState& upFs)
846 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
847 const auto& surfaceVolumeFlux = invB * extQuants.volumeFlux(phaseIdx);
848 evalPhaseFluxes_<UpEval>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, upFs);
855 template <
class UpEval,
class Eval,
class Flu
idState,
class RateVectorT = RateVector>
858 unsigned pvtRegionIdx,
859 const Eval& surfaceVolumeFlux,
860 const FluidState& upFs)
862 const FluidSystem& fsys = upFs.fluidSystem();
864 unsigned activeCompIdx
865 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
867 if constexpr (blackoilConserveSurfaceVolume) {
868 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
870 flux[conti0EqIdx + activeCompIdx]
871 += surfaceVolumeFlux * fsys.referenceDensity(phaseIdx, pvtRegionIdx);
874 if (phaseIdx == oilPhaseIdx) {
876 if (fsys.enableDissolvedGas()) {
878 = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
880 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
881 if constexpr (blackoilConserveSurfaceVolume) {
882 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
884 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux
885 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
888 }
else if (phaseIdx == waterPhaseIdx) {
890 if (fsys.enableDissolvedGasInWater()) {
892 = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
894 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
895 if constexpr (blackoilConserveSurfaceVolume) {
896 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
898 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux
899 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
902 }
else if (phaseIdx == gasPhaseIdx) {
904 if (fsys.enableVaporizedOil()) {
906 = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
908 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
909 if constexpr (blackoilConserveSurfaceVolume) {
910 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
912 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux
913 * fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
917 if (fsys.enableVaporizedWater()) {
919 = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
921 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
922 if constexpr (blackoilConserveSurfaceVolume) {
923 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
925 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux
926 * fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
940 template <
class Scalar>
942 unsigned pvtRegionIdx)
961 template <
class ScalarVector,
class FsysType>
963 unsigned pvtRegionIdx,
964 const FsysType& fsys)
966 if constexpr (!blackoilConserveSurfaceVolume) {
971 if constexpr (waterEnabled) {
972 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
973 container[conti0EqIdx + activeWaterCompIdx]
974 *= fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
977 if constexpr (gasEnabled) {
978 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
979 container[conti0EqIdx + activeGasCompIdx]
980 *= fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
983 if constexpr (oilEnabled) {
984 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
985 container[conti0EqIdx + activeOilCompIdx]
986 *= fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
994 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:95
static void applyScaling(RateVector &flux)
Definition: blackoilbioeffectsmodules.hh:238
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbioeffectsmodules.hh:179
static void addSource(RateVector &source, const Problem &problem, const IntensiveQuantities &intQuants, unsigned globalSpaceIdex)
Definition: blackoilbioeffectsmodules.hh:279
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:58
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:165
Definition: blackoilconvectivemixingmodule.hh:99
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: blackoildiffusionmodule.hh:54
Provides the auxiliary methods required for consideration of the dispersion equation.
Definition: blackoildispersionmodule.hh:56
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:63
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:158
static OPM_HOST_DEVICE void addHeatFlux(RateVectorT &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:221
Contains the high level supplements required to extend the black oil model.
Definition: blackoilextbomodules.hh:64
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilextbomodules.hh:158
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition: blackoilfoammodules.hh:59
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilfoammodules.hh:165
Calculates the local residual of the black oil model.
Definition: blackoillocalresidualtpfa.hh:71
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidualtpfa.hh:820
static void computeBoundaryFluxFree(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:621
static OPM_HOST_DEVICE void computeStorage(StorageType &storage, const IntensiveQuantitiesType &intQuants)
Definition: blackoillocalresidualtpfa.hh:169
static void computeBoundaryFluxRate(RateVector &bdyFlux, const BoundaryConditionData &bdyInfo)
Definition: blackoillocalresidualtpfa.hh:615
static FaceDir::DirEnum faceDirFromDirId(const int dirId)
Definition: blackoillocalresidualtpfa.hh:992
static OPM_HOST_DEVICE void computeBoundaryThermal(const ProblemLocal &problem, RateVectorLocal &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:734
static OPM_HOST_DEVICE void calculateFluxes_(RateVectorT &flux, RateVectorT &darcy, const IntensiveQuantitiesT &intQuantsIn, const IntensiveQuantitiesT &intQuantsEx, const unsigned &globalIndexIn, const unsigned &globalIndexEx, const ResidualNBInfoT &nbInfo, const ModuleParamsT &moduleParams)
Definition: blackoillocalresidualtpfa.hh:360
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Definition: blackoillocalresidualtpfa.hh:840
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:284
static OPM_HOST_DEVICE void evalPhaseFluxes_(RateVectorT &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:856
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:941
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:159
static OPM_HOST_DEVICE void computeFlux(RateVectorT &flux, RateVectorT &darcy, const unsigned globalIndexIn, const unsigned globalIndexEx, const IntensiveQuantitiesT &intQuantsIn, const IntensiveQuantitiesT &intQuantsEx, const ResidualNBInfoT &nbInfo, const ModuleParamsT &moduleParams)
Definition: blackoillocalresidualtpfa.hh:257
static OPM_HOST_DEVICE void computeBoundaryFlux(RateVectorLocal &bdyFlux, const LocalProblem &problem, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:572
static void computeSource(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:778
static void computeSourceDense(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:798
static OPM_HOST_DEVICE void adaptMassConservationQuantities_(ScalarVector &container, unsigned pvtRegionIdx, const FsysType &fsys)
Helper function to convert the mass-related parts of a vector that stores conservation quantities in ...
Definition: blackoillocalresidualtpfa.hh:962
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:65
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilpolymermodules.hh:237
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:185
@ NONE
Definition: DeferredLogger.hpp:46
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: blackoillocalresidualtpfa.hh:141
FaceDir::DirEnum faceDir
Definition: blackoillocalresidualtpfa.hh:146
double faceArea
Definition: blackoillocalresidualtpfa.hh:143
ConditionalStorage< enableFullyImplicitThermal, double > outAlpha
Definition: blackoillocalresidualtpfa.hh:150
ConditionalStorage< enableFullyImplicitThermal, double > inAlpha
Definition: blackoillocalresidualtpfa.hh:149
ConditionalStorage< enableDiffusion, double > diffusivity
Definition: blackoillocalresidualtpfa.hh:151
double dZg
Definition: blackoillocalresidualtpfa.hh:145
double Vin
Definition: blackoillocalresidualtpfa.hh:147
ConditionalStorage< enableDispersion, double > dispersivity
Definition: blackoillocalresidualtpfa.hh:152
double thpres
Definition: blackoillocalresidualtpfa.hh:144
double trans
Definition: blackoillocalresidualtpfa.hh:142
double Vex
Definition: blackoillocalresidualtpfa.hh:148