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>
43#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
45#include <opm/common/utility/gpuistl_if_available.hpp>
55 template<
class Scalar,
template<
class>
class Storage = Opm::VectorWithDefaultAllocator>
66 template <
class Scalar>
76 template <
class Scalar>
77 ConvectiveMixingModuleParam<Scalar, GpuBuffer> copy_to_gpu(
const ConvectiveMixingModuleParam<Scalar, Opm::VectorWithDefaultAllocator>& params)
79 return ConvectiveMixingModuleParam<Scalar, GpuBuffer>{
80 gpuistl::GpuBuffer(params.active_),
81 gpuistl::GpuBuffer(params.Xhi_),
82 gpuistl::GpuBuffer(params.Psi_)
102template <
class TypeTag,
bool enableConvectiveMixing>
109template <
class TypeTag>
120 static void beginEpisode(
const EclipseState&,
127 template <
class Context>
132 const IntensiveQuantities&,
133 const IntensiveQuantities&,
138 template <
class Context>
150 const IntensiveQuantities&,
151 const IntensiveQuantities&,
161template <
class TypeTag>
171 using Toolbox = MathToolbox<Evaluation>;
174 enum { conti0EqIdx = Indices::conti0EqIdx };
175 enum { dimWorld = GridView::dimensionworld };
176 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
177 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
178 static constexpr bool enableFullyImplicitThermal = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
179 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
183 static void beginEpisode(
const EclipseState& eclState,
184 const Schedule& schedule,
185 const int episodeIdx,
189 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
190 const auto& control = schedule[episodeIdx].oilvap();
192 info.
active_.resize(numRegions);
193 info.
Xhi_.resize(numRegions);
194 info.
Psi_.resize(numRegions);
196 for (std::size_t i = 0; i < numRegions; ++i ) {
197 info.
active_[i] = control.drsdtConvective(i);
198 if (control.drsdtConvective(i)) {
199 info.
Xhi_[i] = control.getMaxDRSDT(i);
200 info.
Psi_[i] = control.getPsi(i);
206 template <
class CMMParam>
208 const IntensiveQuantities& intQuantsIn,
209 const IntensiveQuantities& intQuantsEx,
210 const unsigned phaseIdx,
211 const CMMParam& info) {
213 if (info.active_.empty()) {
216 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
220 FluidSystem fsys = intQuantsIn.getFluidSystem();
221 if (phaseIdx == fsys.gasPhaseIdx) {
225 const auto& liquidPhaseIdx =
226 fsys.phaseIsActive(fsys.waterPhaseIdx)
231 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
232 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
233 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
235 const auto& bLiquidIn =
236 fsys.phaseIsActive(waterPhaseIdx)
237 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
238 t_in, p_in, Evaluation(0.0), salt_in)
239 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
240 t_in, p_in, Evaluation(0.0));
242 const auto& refDensityLiquidIn =
243 fsys.phaseIsActive(waterPhaseIdx)
244 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
245 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
247 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
249 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
250 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
251 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
253 const auto bLiquidEx =
254 fsys.phaseIsActive(waterPhaseIdx)
255 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
256 t_ex, p_ex, Scalar{0.0}, salt_ex)
257 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
258 t_ex, p_ex, Scalar{0.0});
260 const auto& refDensityLiquidEx =
261 fsys.phaseIsActive(waterPhaseIdx)
262 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
263 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
265 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
267 rhoAvg = (rho_in + rho_ex) / 2;
270 template <
class Context>
272 const Context& elemCtx,
277 const auto& problem = elemCtx.problem();
278 const auto& stencil = elemCtx.stencil(timeIdx);
279 const auto& scvf = stencil.interiorFace(scvfIdx);
281 unsigned interiorDofIdx = scvf.interiorIndex();
282 unsigned exteriorDofIdx = scvf.exteriorIndex();
283 assert(interiorDofIdx != exteriorDofIdx);
285 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
286 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
287 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
288 Scalar faceArea = scvf.area();
289 const Scalar g = problem.gravity()[dimWorld - 1];
290 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
291 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
292 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
293 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
294 const Scalar distZ = zIn - zEx;
295 addConvectiveMixingFlux(flux,
303 problem.moduleParams().convectiveMixingModuleParam);
310 template <
class RateVectorT = RateVector,
class CMMParam = ConvectiveMixingModuleParamT>
312 const IntensiveQuantities& intQuantsIn,
313 const IntensiveQuantities& intQuantsEx,
314 const unsigned globalIndexIn,
315 const unsigned globalIndexEx,
318 const Scalar faceArea,
319 const CMMParam& info)
321 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
323 if (info.active_.empty()) {
327 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
331 const auto& liquidPhaseIdx =
332 fsys.phaseIsActive(fsys.waterPhaseIdx)
337 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
338 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
339 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
340 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
342 const auto bLiquidSatIn =
343 fsys.phaseIsActive(fsys.waterPhaseIdx)
344 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
345 t_in, p_in, rssat_in, salt_in)
346 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
347 t_in, p_in, rssat_in);
349 const auto& densityLiquidIn =
350 fsys.phaseIsActive(fsys.waterPhaseIdx)
351 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
352 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
354 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
355 const auto rho_sat_in = bLiquidSatIn *
357 rssat_in * fsys.referenceDensity(fsys.gasPhaseIdx,
358 intQuantsIn.pvtRegionIndex()));
361 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
362 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
363 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
364 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
365 const auto bLiquidSatEx =
366 fsys.phaseIsActive(fsys.waterPhaseIdx)
367 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
368 t_ex, p_ex, rssat_ex, salt_ex)
369 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
370 t_ex, p_ex, rssat_ex);
372 const auto& densityLiquidEx =
373 fsys.phaseIsActive(fsys.waterPhaseIdx)
374 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
375 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
377 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
378 const auto rho_sat_ex = bLiquidSatEx *
380 rssat_ex * fsys.referenceDensity(fsys.gasPhaseIdx,
381 intQuantsEx.pvtRegionIndex()));
383 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
384 const auto pressure_difference_convective_mixing = delta_rho * distZg;
387 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
389 short interiorDofIdx = 0;
390 constexpr short exteriorDofIdx = 1;
393 if (pressure_difference_convective_mixing > 0) {
394 upIdx = exteriorDofIdx;
397 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
398 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
399 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
401 fsys.phaseIsActive(fsys.waterPhaseIdx)
402 ? up.fluidState().Rsw()
403 : up.fluidState().Rs();
405 const Evaluation& transMult = up.rockCompTransMultiplier();
406 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
407 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
410 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
412 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
413 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
414 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(fsys.gasCompIdx);
415 if (globalUpIndex == globalIndexIn) {
416 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
419 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
422 if constexpr (enableFullyImplicitThermal) {
423 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
424 fsys.referenceDensity(fsys.gasPhaseIdx, up.pvtRegionIndex());
425 if (globalUpIndex == globalIndexIn) {
426 flux[contiEnergyEqIdx] += convectiveFlux * h;
429 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
436template <
class TypeTag,
bool enableConvectiveMixingV>
446template <
class TypeTag>
460 const auto liquidPhaseIdx =
461 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
462 ? FluidSystem::waterPhaseIdx
463 : FluidSystem::oilPhaseIdx;
465 const Evaluation SoMax = 0.0;
466 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
468 asImp_().pvtRegionIndex(),
473 {
return saturatedDissolutionFactor_; }
477 {
return *
static_cast<Implementation*
>(
this); }
482template <
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:458
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:476
OPM_HOST_DEVICE const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:472
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:479
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
Definition: blackoilconvectivemixingmodule.hh:437
static void addConvectiveMixingFlux(RateVector &, const Context &, unsigned, unsigned)
Definition: blackoilconvectivemixingmodule.hh:139
static bool active(const Context &)
Definition: blackoilconvectivemixingmodule.hh:128
static void modifyAvgDensity(Evaluation &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned int, const ConvectiveMixingModuleParam< Scalar > &)
Definition: blackoilconvectivemixingmodule.hh:131
static void addConvectiveMixingFlux(RateVector &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned, const unsigned, const Scalar, const Scalar, const Scalar, const ConvectiveMixingModuleParam< Scalar > &)
Adds the convective mixing mass flux flux to the flux vector over a flux integration point.
Definition: blackoilconvectivemixingmodule.hh:149
static OPM_HOST_DEVICE void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const CMMParam &info)
Definition: blackoilconvectivemixingmodule.hh:207
static OPM_HOST_DEVICE void addConvectiveMixingFlux(RateVectorT &flux, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned globalIndexIn, const unsigned globalIndexEx, const Scalar distZg, const Scalar trans, const Scalar faceArea, const CMMParam &info)
Adds the convective mixing mass flux flux to the flux vector over a flux integration point.
Definition: blackoilconvectivemixingmodule.hh:311
static OPM_HOST_DEVICE void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:271
Definition: blackoilconvectivemixingmodule.hh:103
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.
PointerView< T > make_view(const std::shared_ptr< T > &ptr)
Definition: gpu_smart_pointer.hpp:318
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
Definition: blackoilconvectivemixingmodule.hh:57
Storage< Scalar > Xhi_
Definition: blackoilconvectivemixingmodule.hh:59
Storage< bool > active_
Definition: blackoilconvectivemixingmodule.hh:58
Storage< Scalar > Psi_
Definition: blackoilconvectivemixingmodule.hh:60