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>
46#include <opm/common/ErrorMacros.hpp>
47#include <opm/common/utility/gpuDecorators.hpp>
54#include <opm/common/utility/gpuistl_if_available.hpp>
64template <
class TypeTag>
77 using FluidState =
typename IntensiveQuantities::FluidState;
79 enum { conti0EqIdx = Indices::conti0EqIdx };
80 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
81 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
82 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
84 enum { dimWorld = GridView::dimensionworld };
85 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
86 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
87 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
89 enum { gasCompIdx = FluidSystem::gasCompIdx };
90 enum { oilCompIdx = FluidSystem::oilCompIdx };
91 enum { waterCompIdx = FluidSystem::waterCompIdx };
92 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
94 static constexpr bool waterEnabled = Indices::waterEnabled;
95 static constexpr bool gasEnabled = Indices::gasEnabled;
96 static constexpr bool oilEnabled = Indices::oilEnabled;
97 static constexpr bool compositionSwitchEnabled = compositionSwitchIdx >= 0;
99 static constexpr bool blackoilConserveSurfaceVolume
100 = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
102 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
103 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
104 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
105 static constexpr bool enableFullyImplicitThermal
106 = (getPropValue<TypeTag, Properties::EnergyModuleType>()
107 == EnergyModules::FullyImplicitThermal);
108 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
109 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
110 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
111 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
112 static constexpr bool enableConvectiveMixing
113 = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
114 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
115 static constexpr bool enableSaltPrecipitation
116 = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
117 static constexpr bool enableMICP = Indices::enableMICP;
118 static constexpr bool runAssemblyOnGpu = getPropValue<TypeTag, Properties::RunAssemblyOnGpu>();
120 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
121 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
122 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
126 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
132 using Toolbox = MathToolbox<Evaluation>;
144 ConditionalStorage<enableFullyImplicitThermal, double>
inAlpha;
145 ConditionalStorage<enableFullyImplicitThermal, double>
outAlpha;
153 template <
class LhsEval>
155 const ElementContext& elemCtx,
157 unsigned timeIdx)
const
159 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
160 computeStorage<LhsEval>(storage, intQuants);
163 template <
class LhsEval,
class StorageType,
class IntensiveQuantitiesType = IntensiveQuantities>
165 const IntensiveQuantitiesType& intQuants)
169 const auto& fs = intQuants.fluidState();
172 const FluidSystem& fsys = intQuants.getFluidSystem();
174 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
175 if (!fsys.phaseIsActive(phaseIdx)) {
178 unsigned activeCompIdx
179 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
180 LhsEval surfaceVolume = Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
181 * Toolbox::template decay<LhsEval>(fs.invB(phaseIdx))
182 * Toolbox::template decay<LhsEval>(intQuants.porosity());
184 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
187 if (phaseIdx == oilPhaseIdx && fsys.enableDissolvedGas()) {
188 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
189 storage[conti0EqIdx + activeGasCompIdx]
190 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs())
195 if (phaseIdx == waterPhaseIdx && fsys.enableDissolvedGasInWater()) {
196 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
197 storage[conti0EqIdx + activeGasCompIdx]
198 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw())
203 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedOil()) {
204 unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
205 storage[conti0EqIdx + activeOilCompIdx]
206 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv())
211 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedWater()) {
212 unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
213 storage[conti0EqIdx + activeWaterCompIdx]
214 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw())
222 if constexpr (enableSolvent) {
223 SolventModule::addStorage(storage, intQuants);
227 if constexpr (enableExtbo) {
228 ExtboModule::addStorage(storage, intQuants);
232 if constexpr (enablePolymer) {
233 PolymerModule::addStorage(storage, intQuants);
240 if constexpr (enableFoam) {
241 FoamModule::addStorage(storage, intQuants);
245 if constexpr (enableBrine) {
246 BrineModule::addStorage(storage, intQuants);
250 if constexpr (enableBioeffects) {
251 BioeffectsModule::addStorage(storage, intQuants);
260 template <
class ModuleParamsT,
262 class IntensiveQuantitiesT,
263 class ResidualNBInfoT>
266 const unsigned globalIndexIn,
267 const unsigned globalIndexEx,
268 const IntensiveQuantitiesT& intQuantsIn,
269 const IntensiveQuantitiesT& intQuantsEx,
270 const ResidualNBInfoT& nbInfo,
271 const ModuleParamsT& moduleParams)
273 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
291 computeFlux(RateVector& flux,
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
293 OPM_TIMEBLOCK_LOCAL(
computeFlux, Subsystem::Assembly);
294 assert(timeIdx == 0);
297 RateVector darcy = 0.0;
299 const auto& problem = elemCtx.problem();
300 const auto& stencil = elemCtx.stencil(timeIdx);
301 const auto& scvf = stencil.interiorFace(scvfIdx);
303 unsigned interiorDofIdx = scvf.interiorIndex();
304 unsigned exteriorDofIdx = scvf.exteriorIndex();
305 assert(interiorDofIdx != exteriorDofIdx);
309 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, 0);
310 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, 0);
311 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
312 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
313 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
314 Scalar faceArea = scvf.area();
316 Scalar thpresInToEx = problem.thresholdPressure(globalIndexIn, globalIndexEx);
317 Scalar thpresExToIn = problem.thresholdPressure(globalIndexEx, globalIndexIn);
322 const Scalar g = problem.gravity()[dimWorld - 1];
323 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
324 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
331 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
332 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
335 const Scalar distZ = zIn - zEx;
337 const Scalar inAlpha = problem.thermalHalfTransmissibility(globalIndexIn, globalIndexEx);
338 const Scalar outAlpha = problem.thermalHalfTransmissibility(globalIndexEx, globalIndexIn);
339 const Scalar diffusivity = problem.diffusivity(globalIndexEx, globalIndexIn);
340 const Scalar dispersivity = problem.dispersivity(globalIndexEx, globalIndexIn);
362 problem.moduleParams());
365 template <
class RateVectorT,
366 class IntensiveQuantitiesT,
367 class ResidualNBInfoT,
371 const IntensiveQuantitiesT& intQuantsIn,
372 const IntensiveQuantitiesT& intQuantsEx,
373 const unsigned& globalIndexIn,
374 const unsigned& globalIndexEx,
375 const ResidualNBInfoT& nbInfo,
376 const ModuleParamsT& moduleParams)
378 OPM_TIMEBLOCK_LOCAL(calculateFluxes, Subsystem::Assembly);
379 const Scalar Vin = nbInfo.Vin;
380 const Scalar Vex = nbInfo.Vex;
381 const Scalar distZg = nbInfo.dZg;
382 const Scalar thpresInToEx = nbInfo.thpresInToEx;
383 const Scalar thpresExToIn = nbInfo.thpresExToIn;
384 const Scalar trans = nbInfo.trans;
385 const Scalar faceArea = nbInfo.faceArea;
386 FaceDir::DirEnum facedir = nbInfo.faceDir;
388 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
390 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
391 if (!fsys.phaseIsActive(phaseIdx)) {
399 short interiorDofIdx = 0;
400 short exteriorDofIdx = 1;
401 Evaluation pressureDifference;
402 ExtensiveQuantities::calculatePhasePressureDiff_(upIdx,
419 const IntensiveQuantities& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
420 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
423 Evaluation transMult = (intQuantsIn.rockCompTransMultiplier()
424 + Toolbox::value(intQuantsEx.rockCompTransMultiplier()))
426 if constexpr (enableBioeffects || enableSaltPrecipitation) {
428 *= (intQuantsIn.permFactor() + Toolbox::value(intQuantsEx.permFactor())) / 2;
431 Evaluation darcyFlux;
432 if (globalUpIndex == globalIndexIn) {
433 darcyFlux = pressureDifference * up.mobility(phaseIdx, facedir) * transMult
434 * (-trans / faceArea);
436 darcyFlux = pressureDifference
437 * (Toolbox::value(up.mobility(phaseIdx, facedir)) * transMult
438 * (-trans / faceArea));
441 unsigned activeCompIdx
442 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
444 darcy[conti0EqIdx + activeCompIdx] = darcyFlux.value() * faceArea;
446 unsigned pvtRegionIdx = up.pvtRegionIndex();
448 if (globalUpIndex == globalIndexIn) {
449 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
450 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
451 const auto& surfaceVolumeFlux = invB * darcyFlux;
453 evalPhaseFluxes_<Evaluation>(
454 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
455 if constexpr (enableFullyImplicitThermal) {
456 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
457 flux, phaseIdx, darcyFlux, up.fluidState());
459 if constexpr (enableBioeffects) {
460 BioeffectsModule::template addBioeffectsFluxes_<Evaluation>(
461 flux, phaseIdx, darcyFlux, up);
463 if constexpr (enableBrine) {
464 BrineModule::template addBrineFluxes_<Evaluation, FluidState>(
465 flux, phaseIdx, darcyFlux, up.fluidState());
468 const auto& invB = getInvB_<FluidSystem, FluidState, Scalar>(
469 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
470 const auto& surfaceVolumeFlux = invB * darcyFlux;
471 evalPhaseFluxes_<Scalar>(
472 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
473 if constexpr (enableFullyImplicitThermal) {
474 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
475 flux, phaseIdx, darcyFlux, up.fluidState());
477 if constexpr (enableBioeffects) {
478 BioeffectsModule::template addBioeffectsFluxes_<Scalar>(
479 flux, phaseIdx, darcyFlux, up);
481 if constexpr (enableBrine) {
482 BrineModule::template addBrineFluxes_<Scalar, FluidState>(
483 flux, phaseIdx, darcyFlux, up.fluidState());
491 "Relevant computeFlux() method must be implemented for this module before enabling.");
497 "Relevant computeFlux() method must be implemented for this module before enabling.");
503 "Relevant computeFlux() method must be implemented for this module before enabling.");
507 if constexpr (enableConvectiveMixing) {
508 ConvectiveMixingModule::addConvectiveMixingFlux(
517 moduleParams.convectiveMixingModuleParam);
521 if constexpr (enableFullyImplicitThermal) {
522 const Scalar inAlpha = nbInfo.inAlpha;
523 const Scalar outAlpha = nbInfo.outAlpha;
526 short interiorDofIdx = 0;
527 short exteriorDofIdx = 1;
529 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
535 intQuantsIn.fluidState(),
536 intQuantsEx.fluidState(),
548 "Relevant computeFlux() method must be implemented for this module before enabling.");
553 if constexpr (enableDiffusion) {
554 typename DiffusionModule::ExtensiveQuantities::EvaluationArray
555 effectiveDiffusionCoefficient;
556 DiffusionModule::ExtensiveQuantities::update(
557 effectiveDiffusionCoefficient, intQuantsIn, intQuantsEx);
558 const Scalar diffusivity = nbInfo.diffusivity;
559 const Scalar tmpdiffusivity = diffusivity / faceArea;
560 DiffusionModule::addDiffusiveFlux(
561 flux, intQuantsIn, intQuantsEx, tmpdiffusivity, effectiveDiffusionCoefficient);
566 if constexpr (enableDispersion) {
567 typename DispersionModule::ExtensiveQuantities::ScalarArray normVelocityAvg;
568 DispersionModule::ExtensiveQuantities::update(
569 normVelocityAvg, intQuantsIn, intQuantsEx);
570 const Scalar dispersivity = nbInfo.dispersivity;
571 const Scalar tmpdispersivity = dispersivity / faceArea;
572 DispersionModule::addDispersiveFlux(
573 flux, intQuantsIn, intQuantsEx, tmpdispersivity, normVelocityAvg);
577 if constexpr (enableMICP) {
578 BioeffectsModule::applyScaling(flux);
582 template <
class BoundaryConditionData,
class RateVectorLocal,
class LocalProblem>
584 const LocalProblem& problem,
585 const BoundaryConditionData& bdyInfo,
586 const IntensiveQuantities& insideIntQuants,
587 unsigned globalSpaceIdx)
589#if OPM_IS_INSIDE_HOST_FUNCTION
590 switch (bdyInfo.type) {
598 case BCType::DIRICHLET:
601 case BCType::THERMAL:
605 throw std::logic_error(
"Unknown boundary condition type "
607 +
" in computeBoundaryFlux().");
610 switch (bdyInfo.type) {
614 case BCType::THERMAL:
618 OPM_THROW(std::logic_error,
619 "Boundary condition type " +
std::to_string(
static_cast<int>(bdyInfo.type))
620 +
" is not supported for GPU fluid systems in computeBoundaryFlux().");
625 template <
class BoundaryConditionData>
628 bdyFlux.setMassRate(bdyInfo.massRate, bdyInfo.pvtRegionIdx);
631 template <
class BoundaryConditionData>
634 const BoundaryConditionData& bdyInfo,
635 const IntensiveQuantities& insideIntQuants,
636 unsigned globalSpaceIdx)
639 std::array<short, numPhases> upIdx;
640 std::array<short, numPhases> dnIdx;
641 std::array<Evaluation, numPhases> volumeFlux;
642 std::array<Evaluation, numPhases> pressureDifference;
644 ExtensiveQuantities::calculateBoundaryGradients_(problem,
647 bdyInfo.boundaryFaceIndex,
650 bdyInfo.exFluidState,
660 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
661 if (!FluidSystem::phaseIsActive(phaseIdx)) {
664 const auto& pBoundary = bdyInfo.exFluidState.pressure(phaseIdx);
665 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
666 const unsigned pvtRegionIdx = insideIntQuants.pvtRegionIndex();
669 const auto& darcyFlux = volumeFlux[phaseIdx];
671 if (pBoundary < pInside) {
673 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
674 insideIntQuants.fluidState(), phaseIdx, pvtRegionIdx);
675 Evaluation surfaceVolumeFlux = invB * darcyFlux;
676 evalPhaseFluxes_<Evaluation>(tmp,
678 insideIntQuants.pvtRegionIndex(),
680 insideIntQuants.fluidState());
681 if constexpr (enableFullyImplicitThermal) {
682 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
683 tmp, phaseIdx, darcyFlux, insideIntQuants.fluidState());
685 }
else if (pBoundary > pInside) {
687 using ScalarFluidState =
decltype(bdyInfo.exFluidState);
688 const auto& invB = getInvB_<FluidSystem, ScalarFluidState, Scalar>(
689 bdyInfo.exFluidState, phaseIdx, pvtRegionIdx);
690 Evaluation surfaceVolumeFlux = invB * darcyFlux;
691 evalPhaseFluxes_<Scalar>(tmp,
693 insideIntQuants.pvtRegionIndex(),
695 bdyInfo.exFluidState);
696 if constexpr (enableFullyImplicitThermal) {
697 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
698 tmp, phaseIdx, darcyFlux, bdyInfo.exFluidState);
702 for (
unsigned i = 0; i < tmp.size(); ++i) {
703 bdyFlux[i] += tmp[i];
708 if constexpr (enableFullyImplicitThermal) {
711 Scalar alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
712 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
715 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
720 bdyInfo.exFluidState);
726 "Relevant treatment of boundary conditions must be implemented before enabling.");
729 "Relevant treatment of boundary conditions must be implemented before enabling.");
731 const FluidSystem& fsys = insideIntQuants.getFluidSystem();
737 for (
unsigned i = 0; i < numEq; ++i) {
738 Valgrind::CheckDefined(bdyFlux[i]);
740 Valgrind::CheckDefined(bdyFlux);
744 template <
class ProblemLocal,
class BoundaryConditionData,
class RateVectorLocal>
746 RateVectorLocal& bdyFlux,
747 const BoundaryConditionData& bdyInfo,
748 const IntensiveQuantities& insideIntQuants,
749 [[maybe_unused]]
unsigned globalSpaceIdx)
756 if constexpr (enableFullyImplicitThermal) {
761 if constexpr (runAssemblyOnGpu) {
764 alpha = problem.getAlpha(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
766 alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
767 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
772 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
777 bdyInfo.exFluidState);
782 for (
unsigned i = 0; i < numEq; ++i) {
783 Valgrind::CheckDefined(bdyFlux[i]);
785 Valgrind::CheckDefined(bdyFlux);
790 const Problem& problem,
791 const IntensiveQuantities& insideIntQuants,
792 unsigned globalSpaceIdex,
797 problem.source(source, globalSpaceIdex, timeIdx);
800 if constexpr (enableBioeffects) {
801 BioeffectsModule::addSource(source, problem, insideIntQuants, globalSpaceIdex);
805 if constexpr (enableFullyImplicitThermal) {
806 source[Indices::contiEnergyEqIdx]
807 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
812 const Problem& problem,
813 const IntensiveQuantities& insideIntQuants,
814 unsigned globalSpaceIdex,
818 problem.addToSourceDense(source, globalSpaceIdex, timeIdx);
821 if constexpr (enableBioeffects) {
822 BioeffectsModule::addSource(source, problem, insideIntQuants, globalSpaceIdex);
826 if constexpr (enableFullyImplicitThermal) {
827 source[Indices::contiEnergyEqIdx]
828 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
836 const ElementContext& elemCtx,
838 unsigned timeIdx)
const
842 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
845 if constexpr (enableBioeffects) {
846 BioeffectsModule::addSource(source, elemCtx, dofIdx, timeIdx);
850 if constexpr (enableFullyImplicitThermal) {
851 source[Indices::contiEnergyEqIdx]
852 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
856 template <
class UpEval,
class Flu
idState>
859 unsigned pvtRegionIdx,
860 const ExtensiveQuantities& extQuants,
861 const FluidState& upFs)
863 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
864 const auto& surfaceVolumeFlux = invB * extQuants.volumeFlux(phaseIdx);
865 evalPhaseFluxes_<UpEval>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, upFs);
872 template <
class UpEval,
class Eval,
class Flu
idState,
class RateVectorT = RateVector>
875 unsigned pvtRegionIdx,
876 const Eval& surfaceVolumeFlux,
877 const FluidState& upFs)
879 const FluidSystem& fsys = upFs.fluidSystem();
881 unsigned activeCompIdx
882 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
884 if constexpr (blackoilConserveSurfaceVolume) {
885 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
887 flux[conti0EqIdx + activeCompIdx]
888 += surfaceVolumeFlux * fsys.referenceDensity(phaseIdx, pvtRegionIdx);
891 if (phaseIdx == oilPhaseIdx) {
893 if (fsys.enableDissolvedGas()) {
895 = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
897 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
898 if constexpr (blackoilConserveSurfaceVolume) {
899 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
901 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux
902 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
905 }
else if (phaseIdx == waterPhaseIdx) {
907 if (fsys.enableDissolvedGasInWater()) {
909 = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
911 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
912 if constexpr (blackoilConserveSurfaceVolume) {
913 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
915 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux
916 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
919 }
else if (phaseIdx == gasPhaseIdx) {
921 if (fsys.enableVaporizedOil()) {
923 = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
925 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
926 if constexpr (blackoilConserveSurfaceVolume) {
927 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
929 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux
930 * fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
934 if (fsys.enableVaporizedWater()) {
936 = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
938 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
939 if constexpr (blackoilConserveSurfaceVolume) {
940 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
942 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux
943 * fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
957 template <
class Scalar>
959 unsigned pvtRegionIdx)
978 template <
class ScalarVector,
class FsysType>
980 unsigned pvtRegionIdx,
981 const FsysType& fsys)
983 if constexpr (!blackoilConserveSurfaceVolume) {
988 if constexpr (waterEnabled) {
989 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
990 container[conti0EqIdx + activeWaterCompIdx]
991 *= fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
994 if constexpr (gasEnabled) {
995 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
996 container[conti0EqIdx + activeGasCompIdx]
997 *= fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
1000 if constexpr (oilEnabled) {
1001 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
1002 container[conti0EqIdx + activeOilCompIdx]
1003 *= fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
1011 return dirId < 0 ? FaceDir::DirEnum::Unknown : FaceDir::FromIntersectionIndex(dirId);
Classes required for dynamic convective mixing.
Contains the classes required to extend the black-oil model by energy.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
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.
Provides the auxiliary methods required for consideration of the diffusion equation.
Provides the auxiliary methods required for consideration of the dispersion equation.
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:64
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:159
static OPM_HOST_DEVICE void addHeatFlux(RateVectorT &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:222
Definition: blackoilfoammodules.hh:57
Calculates the local residual of the black oil model.
Definition: blackoillocalresidualtpfa.hh:66
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidualtpfa.hh:835
static void computeBoundaryFluxFree(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:632
static OPM_HOST_DEVICE void computeStorage(StorageType &storage, const IntensiveQuantitiesType &intQuants)
Definition: blackoillocalresidualtpfa.hh:164
static void computeBoundaryFluxRate(RateVector &bdyFlux, const BoundaryConditionData &bdyInfo)
Definition: blackoillocalresidualtpfa.hh:626
static FaceDir::DirEnum faceDirFromDirId(const int dirId)
Definition: blackoillocalresidualtpfa.hh:1009
static OPM_HOST_DEVICE void computeBoundaryThermal(const ProblemLocal &problem, RateVectorLocal &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:745
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:369
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Definition: blackoillocalresidualtpfa.hh:857
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:291
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:873
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:958
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:154
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:264
static OPM_HOST_DEVICE void computeBoundaryFlux(RateVectorLocal &bdyFlux, const LocalProblem &problem, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:583
static void computeSource(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:789
static void computeSourceDense(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:811
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:979
Definition: blackoilpolymermodules.hh:59
Definition: blackoilsolventmodules.hh:64
@ 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:135
FaceDir::DirEnum faceDir
Definition: blackoillocalresidualtpfa.hh:141
double faceArea
Definition: blackoillocalresidualtpfa.hh:137
double thpresInToEx
Definition: blackoillocalresidualtpfa.hh:138
ConditionalStorage< enableFullyImplicitThermal, double > outAlpha
Definition: blackoillocalresidualtpfa.hh:145
ConditionalStorage< enableFullyImplicitThermal, double > inAlpha
Definition: blackoillocalresidualtpfa.hh:144
ConditionalStorage< enableDiffusion, double > diffusivity
Definition: blackoillocalresidualtpfa.hh:146
double thpresExToIn
Definition: blackoillocalresidualtpfa.hh:139
double dZg
Definition: blackoillocalresidualtpfa.hh:140
double Vin
Definition: blackoillocalresidualtpfa.hh:142
ConditionalStorage< enableDispersion, double > dispersivity
Definition: blackoillocalresidualtpfa.hh:147
double trans
Definition: blackoillocalresidualtpfa.hh:136
double Vex
Definition: blackoillocalresidualtpfa.hh:143