28#ifndef EWOMS_CONVECTIVEMIXING_MODULE_HH
29#define EWOMS_CONVECTIVEMIXING_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/input/eclipse/Schedule/OilVaporizationProperties.hpp>
34#include <opm/input/eclipse/Schedule/Schedule.hpp>
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/common/Valgrind.hpp>
64template <
class TypeTag,
bool enableConvectiveMixing>
71template <
class TypeTag>
80 struct ConvectiveMixingModuleParam
84 static void beginEpisode(
const EclipseState&,
87 ConvectiveMixingModuleParam&)
91 template <
class Context>
96 const IntensiveQuantities&,
97 const IntensiveQuantities&,
99 const ConvectiveMixingModuleParam&)
102 template <
class Context>
114 const IntensiveQuantities&,
115 const IntensiveQuantities&,
121 const ConvectiveMixingModuleParam&)
125template <
class TypeTag>
135 using Toolbox = MathToolbox<Evaluation>;
137 enum { conti0EqIdx = Indices::conti0EqIdx };
138 enum { dimWorld = GridView::dimensionworld };
139 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
140 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
141 static constexpr bool enableEnergy = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
142 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
145 struct ConvectiveMixingModuleParam
153 static void beginEpisode(
const EclipseState& eclState,
154 const Schedule& schedule,
155 const int episodeIdx,
156 ConvectiveMixingModuleParam& info)
159 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
160 const auto& control = schedule[episodeIdx].oilvap();
161 if (info.active_.empty()) {
162 info.active_.resize(numRegions);
163 info.Xhi_.resize(numRegions);
164 info.Psi_.resize(numRegions);
166 for (std::size_t i = 0; i < numRegions; ++i ) {
167 info.active_[i] = control.drsdtConvective(i);
168 if (control.drsdtConvective(i)) {
169 info.Xhi_[i] = control.getMaxDRSDT(i);
170 info.Psi_[i] = control.getPsi(i);
177 const IntensiveQuantities& intQuantsIn,
178 const IntensiveQuantities& intQuantsEx,
179 const unsigned phaseIdx,
180 const ConvectiveMixingModuleParam& info) {
182 if (info.active_.empty()) {
185 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
189 if (phaseIdx == FluidSystem::gasPhaseIdx) {
193 const auto& liquidPhaseIdx =
194 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
195 ? FluidSystem::waterPhaseIdx
196 : FluidSystem::oilPhaseIdx;
199 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
200 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
201 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
203 const auto& bLiquidIn =
204 FluidSystem::phaseIsActive(waterPhaseIdx)
205 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
206 t_in, p_in, Evaluation(0.0), salt_in)
207 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
208 t_in, p_in, Evaluation(0.0));
210 const auto& refDensityLiquidIn =
211 FluidSystem::phaseIsActive(waterPhaseIdx)
212 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
213 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
215 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
217 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
218 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
219 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
221 const auto bLiquidEx =
222 FluidSystem::phaseIsActive(waterPhaseIdx)
223 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
224 t_ex, p_ex, Scalar{0.0}, salt_ex)
225 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
226 t_ex, p_ex, Scalar{0.0});
228 const auto& refDensityLiquidEx =
229 FluidSystem::phaseIsActive(waterPhaseIdx)
230 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
231 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
233 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
235 rhoAvg = (rho_in + rho_ex) / 2;
238 template <
class Context>
240 const Context& elemCtx,
245 const auto& problem = elemCtx.problem();
246 const auto& stencil = elemCtx.stencil(timeIdx);
247 const auto& scvf = stencil.interiorFace(scvfIdx);
249 unsigned interiorDofIdx = scvf.interiorIndex();
250 unsigned exteriorDofIdx = scvf.exteriorIndex();
251 assert(interiorDofIdx != exteriorDofIdx);
253 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
254 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
255 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
256 Scalar faceArea = scvf.area();
257 const Scalar g = problem.gravity()[dimWorld - 1];
258 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
259 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
260 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
261 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
262 const Scalar distZ = zIn - zEx;
263 addConvectiveMixingFlux(flux,
271 problem.moduleParams().convectiveMixingModuleParam);
279 const IntensiveQuantities& intQuantsIn,
280 const IntensiveQuantities& intQuantsEx,
281 const unsigned globalIndexIn,
282 const unsigned globalIndexEx,
285 const Scalar faceArea,
286 const ConvectiveMixingModuleParam& info)
288 if (info.active_.empty()) {
292 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
296 const auto& liquidPhaseIdx =
297 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
298 ? FluidSystem::waterPhaseIdx
299 : FluidSystem::oilPhaseIdx;
302 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
303 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
304 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
305 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
307 const auto bLiquidSatIn =
308 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
309 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
310 t_in, p_in, rssat_in, salt_in)
311 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
312 t_in, p_in, rssat_in);
314 const auto& densityLiquidIn =
315 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
316 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
317 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
319 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
320 const auto rho_sat_in = bLiquidSatIn *
322 rssat_in * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
323 intQuantsIn.pvtRegionIndex()));
326 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
327 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
328 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
329 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
330 const auto bLiquidSatEx =
331 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
332 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
333 t_ex, p_ex, rssat_ex, salt_ex)
334 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
335 t_ex, p_ex, rssat_ex);
337 const auto& densityLiquidEx =
338 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
339 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
340 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
342 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
343 const auto rho_sat_ex = bLiquidSatEx *
345 rssat_ex * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
346 intQuantsEx.pvtRegionIndex()));
348 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
349 const auto pressure_difference_convective_mixing = delta_rho * distZg;
352 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
354 short interiorDofIdx = 0;
355 constexpr short exteriorDofIdx = 1;
358 if (pressure_difference_convective_mixing > 0) {
359 upIdx = exteriorDofIdx;
362 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
363 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
364 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
366 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
367 ? up.fluidState().Rsw()
368 : up.fluidState().Rs();
370 const Evaluation& transMult = up.rockCompTransMultiplier();
371 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
372 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
375 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
377 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
378 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
379 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
380 if (globalUpIndex == globalIndexIn) {
381 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
384 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
387 if constexpr (enableEnergy) {
388 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
389 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, up.pvtRegionIndex());
390 if (globalUpIndex == globalIndexIn) {
391 flux[contiEnergyEqIdx] += convectiveFlux * h;
394 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
401template <
class TypeTag,
bool enableConvectiveMixingV>
411template <
class TypeTag>
425 const auto liquidPhaseIdx =
426 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
427 ? FluidSystem::waterPhaseIdx
428 : FluidSystem::oilPhaseIdx;
430 const Evaluation SoMax = 0.0;
431 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
433 asImp_().pvtRegionIndex(),
438 {
return saturatedDissolutionFactor_; }
442 {
return *
static_cast<Implementation*
>(
this); }
447template <
class TypeTag>
Contains the classes required to extend the black-oil model by energy.
void updateSaturatedDissolutionFactor_()
Compute the intensive quantities needed to handle convective dissolution.
Definition: blackoilconvectivemixingmodule.hh:423
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:441
const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:437
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:444
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
Definition: blackoilconvectivemixingmodule.hh:402
static void addConvectiveMixingFlux(RateVector &, const Context &, unsigned, unsigned)
Definition: blackoilconvectivemixingmodule.hh:103
static bool active(const Context &)
Definition: blackoilconvectivemixingmodule.hh:92
static void modifyAvgDensity(Evaluation &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned int, const ConvectiveMixingModuleParam &)
Definition: blackoilconvectivemixingmodule.hh:95
static void addConvectiveMixingFlux(RateVector &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned, const unsigned, const Scalar, const Scalar, const Scalar, const ConvectiveMixingModuleParam &)
Adds the convective mixing mass flux flux to the flux vector over a flux integration point.
Definition: blackoilconvectivemixingmodule.hh:113
static void addConvectiveMixingFlux(RateVector &flux, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned globalIndexIn, const unsigned globalIndexEx, const Scalar distZg, const Scalar trans, const Scalar faceArea, const ConvectiveMixingModuleParam &info)
Adds the convective mixing mass flux flux to the flux vector over a flux integration point.
Definition: blackoilconvectivemixingmodule.hh:278
static void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:239
static void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const ConvectiveMixingModuleParam &info)
Definition: blackoilconvectivemixingmodule.hh:176
Definition: blackoilconvectivemixingmodule.hh:65
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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::vector< Scalar > Xhi_
Definition: blackoilconvectivemixingmodule.hh:148
std::vector< Scalar > Psi_
Definition: blackoilconvectivemixingmodule.hh:149
std::vector< bool > active_
Definition: blackoilconvectivemixingmodule.hh:147