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>
63template <
class TypeTag,
bool enableConvectiveMixing>
70template <
class TypeTag>
81 enum { conti0EqIdx = Indices::conti0EqIdx };
82 enum { dimWorld = GridView::dimensionworld };
85 struct ConvectiveMixingModuleParam
89 static void beginEpisode(
const EclipseState&,
92 ConvectiveMixingModuleParam&)
96 template <
class Context>
101 const IntensiveQuantities&,
102 const IntensiveQuantities&,
104 const ConvectiveMixingModuleParam&)
107 template <
class Context>
119 const IntensiveQuantities&,
120 const IntensiveQuantities&,
126 const ConvectiveMixingModuleParam&)
130template <
class TypeTag>
140 using Toolbox = MathToolbox<Evaluation>;
142 enum { conti0EqIdx = Indices::conti0EqIdx };
143 enum { dimWorld = GridView::dimensionworld };
144 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
145 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
146 static constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
147 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
150 struct ConvectiveMixingModuleParam
158 static void beginEpisode(
const EclipseState& eclState,
159 const Schedule& schedule,
160 const int episodeIdx,
161 ConvectiveMixingModuleParam& info)
164 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
165 const auto& control = schedule[episodeIdx].oilvap();
166 if (info.active_.empty()) {
167 info.active_.resize(numRegions);
168 info.Xhi_.resize(numRegions);
169 info.Psi_.resize(numRegions);
171 for (std::size_t i = 0; i < numRegions; ++i ) {
172 info.active_[i] = control.drsdtConvective(i);
173 if (control.drsdtConvective(i)) {
174 info.Xhi_[i] = control.getMaxDRSDT(i);
175 info.Psi_[i] = control.getPsi(i);
182 const IntensiveQuantities& intQuantsIn,
183 const IntensiveQuantities& intQuantsEx,
184 const unsigned phaseIdx,
185 const ConvectiveMixingModuleParam& info) {
187 if (info.active_.empty()) {
190 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
194 if (phaseIdx == FluidSystem::gasPhaseIdx) {
198 const auto& liquidPhaseIdx =
199 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
200 ? FluidSystem::waterPhaseIdx
201 : FluidSystem::oilPhaseIdx;
204 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
205 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
206 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
208 const auto& bLiquidIn =
209 FluidSystem::phaseIsActive(waterPhaseIdx)
210 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
211 t_in, p_in, Evaluation(0.0), salt_in)
212 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
213 t_in, p_in, Evaluation(0.0));
215 const auto& refDensityLiquidIn =
216 FluidSystem::phaseIsActive(waterPhaseIdx)
217 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
218 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
220 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
222 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
223 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
224 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
226 const auto bLiquidEx =
227 FluidSystem::phaseIsActive(waterPhaseIdx)
228 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
229 t_ex, p_ex, Scalar{0.0}, salt_ex)
230 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
231 t_ex, p_ex, Scalar{0.0});
233 const auto& refDensityLiquidEx =
234 FluidSystem::phaseIsActive(waterPhaseIdx)
235 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
236 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
238 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
240 rhoAvg = (rho_in + rho_ex) / 2;
243 template <
class Context>
245 const Context& elemCtx,
250 const auto& problem = elemCtx.problem();
251 const auto& stencil = elemCtx.stencil(timeIdx);
252 const auto& scvf = stencil.interiorFace(scvfIdx);
254 unsigned interiorDofIdx = scvf.interiorIndex();
255 unsigned exteriorDofIdx = scvf.exteriorIndex();
256 assert(interiorDofIdx != exteriorDofIdx);
258 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
259 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
260 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
261 Scalar faceArea = scvf.area();
262 const Scalar g = problem.gravity()[dimWorld - 1];
263 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
264 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
265 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
266 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
267 const Scalar distZ = zIn - zEx;
268 addConvectiveMixingFlux(flux,
276 problem.moduleParams().convectiveMixingModuleParam);
284 const IntensiveQuantities& intQuantsIn,
285 const IntensiveQuantities& intQuantsEx,
286 const unsigned globalIndexIn,
287 const unsigned globalIndexEx,
290 const Scalar faceArea,
291 const ConvectiveMixingModuleParam& info)
293 if (info.active_.empty()) {
297 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
301 const auto& liquidPhaseIdx =
302 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
303 ? FluidSystem::waterPhaseIdx
304 : FluidSystem::oilPhaseIdx;
307 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
308 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
309 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
310 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
312 const auto bLiquidSatIn =
313 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
314 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
315 t_in, p_in, rssat_in, salt_in)
316 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
317 t_in, p_in, rssat_in);
319 const auto& densityLiquidIn =
320 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
321 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
322 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
324 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
325 const auto rho_sat_in = bLiquidSatIn *
327 rssat_in * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
328 intQuantsIn.pvtRegionIndex()));
331 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
332 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
333 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
334 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
335 const auto bLiquidSatEx =
336 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
337 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
338 t_ex, p_ex, rssat_ex, salt_ex)
339 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
340 t_ex, p_ex, rssat_ex);
342 const auto& densityLiquidEx =
343 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
344 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
345 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
347 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
348 const auto rho_sat_ex = bLiquidSatEx *
350 rssat_ex * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
351 intQuantsEx.pvtRegionIndex()));
353 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
354 const auto pressure_difference_convective_mixing = delta_rho * distZg;
357 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
359 short interiorDofIdx = 0;
360 constexpr short exteriorDofIdx = 1;
363 if (pressure_difference_convective_mixing > 0) {
364 upIdx = exteriorDofIdx;
367 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
368 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
369 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
371 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
372 ? up.fluidState().Rsw()
373 : up.fluidState().Rs();
375 const Evaluation& transMult = up.rockCompTransMultiplier();
376 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
377 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
380 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
382 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
383 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
384 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
385 if (globalUpIndex == globalIndexIn) {
386 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
389 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
392 if constexpr (enableEnergy) {
393 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
394 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, up.pvtRegionIndex());
395 if (globalUpIndex == globalIndexIn) {
396 flux[contiEnergyEqIdx] += convectiveFlux * h;
399 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
406template <
class TypeTag,
bool enableConvectiveMixingV>
416template <
class TypeTag>
430 const auto liquidPhaseIdx =
431 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
432 ? FluidSystem::waterPhaseIdx
433 : FluidSystem::oilPhaseIdx;
435 const Evaluation SoMax = 0.0;
436 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
438 asImp_().pvtRegionIndex(),
443 {
return saturatedDissolutionFactor_; }
447 {
return *
static_cast<Implementation*
>(
this); }
452template <
class TypeTag>
void updateSaturatedDissolutionFactor_()
Compute the intensive quantities needed to handle convective dissolution.
Definition: blackoilconvectivemixingmodule.hh:428
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:446
const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:442
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:449
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
Definition: blackoilconvectivemixingmodule.hh:407
static void addConvectiveMixingFlux(RateVector &, const Context &, unsigned, unsigned)
Definition: blackoilconvectivemixingmodule.hh:108
static bool active(const Context &)
Definition: blackoilconvectivemixingmodule.hh:97
static void modifyAvgDensity(Evaluation &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned int, const ConvectiveMixingModuleParam &)
Definition: blackoilconvectivemixingmodule.hh:100
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:118
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:283
static void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:244
static void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const ConvectiveMixingModuleParam &info)
Definition: blackoilconvectivemixingmodule.hh:181
Definition: blackoilconvectivemixingmodule.hh:64
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: blackoilboundaryratevector.hh:39
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:153
std::vector< Scalar > Psi_
Definition: blackoilconvectivemixingmodule.hh:154
std::vector< bool > active_
Definition: blackoilconvectivemixingmodule.hh:152