28#ifndef EWOMS_BLACK_OIL_LOCAL_RESIDUAL_HH
29#define EWOMS_BLACK_OIL_LOCAL_RESIDUAL_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/fluidstates/BlackOilFluidState.hpp>
49template <
class TypeTag>
61 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
63 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
64 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
66 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
67 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
68 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
70 enum { gasCompIdx = FluidSystem::gasCompIdx };
71 enum { oilCompIdx = FluidSystem::oilCompIdx };
72 enum { waterCompIdx = FluidSystem::waterCompIdx };
73 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
75 static constexpr bool waterEnabled = Indices::waterEnabled;
76 static constexpr bool gasEnabled = Indices::gasEnabled;
77 static constexpr bool oilEnabled = Indices::oilEnabled;
78 static constexpr bool compositionSwitchEnabled = (compositionSwitchIdx >= 0);
80 static constexpr bool blackoilConserveSurfaceVolume =
81 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
82 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
83 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
84 static constexpr bool enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
85 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
86 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
87 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
88 static constexpr bool enableFullyImplicitThermal =
89 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal;
90 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
91 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
93 using Toolbox = MathToolbox<Evaluation>;
95 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
96 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
97 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
100 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
109 template <
class LhsEval>
111 const ElementContext& elemCtx,
113 unsigned timeIdx)
const
116 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
117 const auto& fs = intQuants.fluidState();
121 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
122 if (!FluidSystem::phaseIsActive(phaseIdx)) {
123 if (Indices::numPhases == 3) {
124 const unsigned activeCompIdx =
125 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
127 storage[conti0EqIdx + activeCompIdx] = variable<LhsEval>(0.0, conti0EqIdx + activeCompIdx);
130 storage[conti0EqIdx + activeCompIdx] = 0.0;
136 const unsigned activeCompIdx =
137 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
138 const LhsEval surfaceVolume =
139 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
140 Toolbox::template decay<LhsEval>(fs.invB(phaseIdx)) *
141 Toolbox::template decay<LhsEval>(intQuants.porosity());
143 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
146 if (phaseIdx == oilPhaseIdx && FluidSystem::enableDissolvedGas()) {
147 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
148 storage[conti0EqIdx + activeGasCompIdx] +=
149 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs()) *
154 if (phaseIdx == waterPhaseIdx && FluidSystem::enableDissolvedGasInWater()) {
155 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
156 storage[conti0EqIdx + activeGasCompIdx] +=
157 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw()) *
162 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedOil()) {
163 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
164 storage[conti0EqIdx + activeOilCompIdx] +=
165 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv()) *
170 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedWater()) {
171 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
172 storage[conti0EqIdx + activeWaterCompIdx] +=
173 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw()) *
181 if constexpr (enableSolvent) {
182 SolventModule::addStorage(storage, intQuants);
186 if constexpr (enableExtbo) {
187 ExtboModule::addStorage(storage, intQuants);
191 if constexpr (enablePolymer) {
192 PolymerModule::addStorage(storage, intQuants);
199 if constexpr (enableFoam) {
200 FoamModule::addStorage(storage, intQuants);
204 if constexpr (enableBrine) {
205 BrineModule::addStorage(storage, intQuants);
209 if constexpr (enableBioeffects) {
210 BioeffectsModule::addStorage(storage, intQuants);
218 const ElementContext& elemCtx,
220 unsigned timeIdx)
const
222 assert(timeIdx == 0);
226 const ExtensiveQuantities& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
227 unsigned focusDofIdx = elemCtx.focusDofIndex();
228 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
229 if (!FluidSystem::phaseIsActive(phaseIdx)) {
233 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
234 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
235 const unsigned pvtRegionIdx = up.pvtRegionIndex();
236 if (upIdx == focusDofIdx) {
237 evalPhaseFluxes_<Evaluation>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
240 evalPhaseFluxes_<Scalar>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
245 if constexpr (enableSolvent) {
246 SolventModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
250 if constexpr (enableExtbo) {
251 ExtboModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
255 if constexpr (enablePolymer) {
256 PolymerModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
263 if constexpr (enableFoam) {
264 FoamModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
268 if constexpr (enableBrine) {
269 BrineModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
273 if constexpr (enableBioeffects) {
274 BioeffectsModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
278 if constexpr (enableDiffusion) {
279 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
283 if constexpr (enableConvectiveMixing) {
284 ConvectiveMixingModule::addConvectiveMixingFlux(flux, elemCtx, scvfIdx, timeIdx);
292 const ElementContext& elemCtx,
294 unsigned timeIdx)
const
297 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
300 if constexpr (enableBioeffects) {
301 BioeffectsModule::addSource(source, elemCtx, dofIdx, timeIdx);
305 if constexpr (enableFullyImplicitThermal) {
306 source[Indices::contiEnergyEqIdx] *=
307 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
315 template <
class UpEval,
class Flu
idState>
318 unsigned pvtRegionIdx,
319 const ExtensiveQuantities& extQuants,
320 const FluidState& upFs)
322 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
323 const auto& surfaceVolumeFlux = invB*extQuants.volumeFlux(phaseIdx);
324 const unsigned activeCompIdx =
325 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
327 if constexpr (blackoilConserveSurfaceVolume) {
328 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
331 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux *
332 FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
335 if (phaseIdx == oilPhaseIdx) {
337 if (FluidSystem::enableDissolvedGas()) {
338 const auto& Rs = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
340 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
341 if constexpr (blackoilConserveSurfaceVolume) {
342 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
345 flux[conti0EqIdx + activeGasCompIdx] +=
346 Rs * surfaceVolumeFlux *
347 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
350 }
else if (phaseIdx == waterPhaseIdx) {
352 if (FluidSystem::enableDissolvedGasInWater()) {
353 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
355 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
356 if constexpr (blackoilConserveSurfaceVolume) {
357 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
360 flux[conti0EqIdx + activeGasCompIdx] +=
361 Rsw * surfaceVolumeFlux *
362 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
366 else if (phaseIdx == gasPhaseIdx) {
368 if (FluidSystem::enableVaporizedOil()) {
369 const auto& Rv = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
371 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
372 if constexpr (blackoilConserveSurfaceVolume) {
373 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
376 flux[conti0EqIdx + activeOilCompIdx] +=
377 Rv * surfaceVolumeFlux *
378 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
382 if (FluidSystem::enableVaporizedWater()) {
383 const auto& Rvw = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
385 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
386 if constexpr (blackoilConserveSurfaceVolume) {
387 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
390 flux[conti0EqIdx + activeWaterCompIdx] +=
391 Rvw * surfaceVolumeFlux *
392 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
409 template <
class Scalar>
411 unsigned pvtRegionIdx)
413 if constexpr (!blackoilConserveSurfaceVolume) {
418 if constexpr (waterEnabled) {
419 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
420 container[conti0EqIdx + activeWaterCompIdx] *=
421 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
424 if constexpr (gasEnabled) {
425 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
426 container[conti0EqIdx + activeGasCompIdx] *=
427 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
430 if constexpr (oilEnabled) {
431 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
432 container[conti0EqIdx + activeOilCompIdx] *=
433 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
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.
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 computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:191
Definition: blackoilfoammodules.hh:57
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:51
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: blackoillocalresidual.hh:217
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidual.hh:291
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Helper function to calculate the flux of mass in terms of conservation quantities via specific fluid ...
Definition: blackoillocalresidual.hh:316
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: blackoillocalresidual.hh:410
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: blackoillocalresidual.hh:110
Definition: blackoilpolymermodules.hh:59
Definition: blackoilsolventmodules.hh:64
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