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>
51 template<
class Scalar,
template<
class>
class Storage = Opm::VectorWithDefaultAllocator>
62 template <
class Scalar>
72 template <
class Scalar>
73 ConvectiveMixingModuleParam<Scalar, GpuBuffer> copy_to_gpu(
const ConvectiveMixingModuleParam<Scalar, Opm::VectorWithDefaultAllocator>& params)
75 return ConvectiveMixingModuleParam<Scalar, GpuBuffer>{
76 gpuistl::GpuBuffer(params.active_),
77 gpuistl::GpuBuffer(params.Xhi_),
78 gpuistl::GpuBuffer(params.Psi_)
98template <
class TypeTag,
bool enableConvectiveMixing>
105template <
class TypeTag>
120 template <
class Context>
125 const IntensiveQuantities&,
126 const IntensiveQuantities&,
131 template <
class Context>
143 const IntensiveQuantities&,
144 const IntensiveQuantities&,
154template <
class TypeTag>
164 using Toolbox = MathToolbox<Evaluation>;
167 enum { conti0EqIdx = Indices::conti0EqIdx };
168 enum { dimWorld = GridView::dimensionworld };
169 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
170 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
171 static constexpr bool enableFullyImplicitThermal = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
172 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
176 const Schedule& schedule,
177 const int episodeIdx,
181 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
182 const auto& control = schedule[episodeIdx].oilvap();
184 info.
active_.resize(numRegions);
185 info.
Xhi_.resize(numRegions);
186 info.
Psi_.resize(numRegions);
188 for (std::size_t i = 0; i < numRegions; ++i ) {
189 info.
active_[i] = control.drsdtConvective(i);
190 if (control.drsdtConvective(i)) {
191 info.
Xhi_[i] = control.getMaxDRSDT(i);
192 info.
Psi_[i] = control.getPsi(i);
197 template <
class CMMParam>
199 const IntensiveQuantities& intQuantsIn,
200 const IntensiveQuantities& intQuantsEx,
201 const unsigned phaseIdx,
202 const CMMParam& info) {
204 if (info.active_.empty()) {
207 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
211 FluidSystem fsys = intQuantsIn.getFluidSystem();
212 if (phaseIdx == fsys.gasPhaseIdx) {
216 const auto& liquidPhaseIdx =
217 fsys.phaseIsActive(fsys.waterPhaseIdx)
222 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
223 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
224 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
226 const auto& bLiquidIn =
227 fsys.phaseIsActive(waterPhaseIdx)
228 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
229 t_in, p_in, Evaluation(0.0), salt_in)
230 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
231 t_in, p_in, Evaluation(0.0));
233 const auto& refDensityLiquidIn =
234 fsys.phaseIsActive(waterPhaseIdx)
235 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
236 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
238 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
240 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
241 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
242 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
244 const auto bLiquidEx =
245 fsys.phaseIsActive(waterPhaseIdx)
246 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
247 t_ex, p_ex, Scalar{0.0}, salt_ex)
248 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
249 t_ex, p_ex, Scalar{0.0});
251 const auto& refDensityLiquidEx =
252 fsys.phaseIsActive(waterPhaseIdx)
253 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
254 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
256 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
258 rhoAvg = (rho_in + rho_ex) / 2;
261 template <
class Context>
263 const Context& elemCtx,
268 const auto& problem = elemCtx.problem();
269 const auto& stencil = elemCtx.stencil(timeIdx);
270 const auto& scvf = stencil.interiorFace(scvfIdx);
272 unsigned interiorDofIdx = scvf.interiorIndex();
273 unsigned exteriorDofIdx = scvf.exteriorIndex();
274 assert(interiorDofIdx != exteriorDofIdx);
276 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
277 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
278 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
279 Scalar faceArea = scvf.area();
280 const Scalar g = problem.gravity()[dimWorld - 1];
281 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
282 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
283 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
284 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
285 const Scalar distZ = zIn - zEx;
286 addConvectiveMixingFlux(flux,
294 problem.moduleParams().convectiveMixingModuleParam);
301 template <
class RateVectorT = RateVector,
class CMMParam = ConvectiveMixingModuleParamT>
303 const IntensiveQuantities& intQuantsIn,
304 const IntensiveQuantities& intQuantsEx,
305 const unsigned globalIndexIn,
306 const unsigned globalIndexEx,
309 const Scalar faceArea,
310 const CMMParam& info)
312 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
314 if (info.active_.empty()) {
318 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
322 const auto& liquidPhaseIdx =
323 fsys.phaseIsActive(fsys.waterPhaseIdx)
328 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
329 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
330 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
331 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
333 const auto bLiquidSatIn =
334 fsys.phaseIsActive(fsys.waterPhaseIdx)
335 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
336 t_in, p_in, rssat_in, salt_in)
337 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
338 t_in, p_in, rssat_in);
340 const auto& densityLiquidIn =
341 fsys.phaseIsActive(fsys.waterPhaseIdx)
342 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
343 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
345 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
346 const auto rho_sat_in = bLiquidSatIn *
348 rssat_in * fsys.referenceDensity(fsys.gasPhaseIdx,
349 intQuantsIn.pvtRegionIndex()));
352 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
353 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
354 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
355 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
356 const auto bLiquidSatEx =
357 fsys.phaseIsActive(fsys.waterPhaseIdx)
358 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
359 t_ex, p_ex, rssat_ex, salt_ex)
360 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
361 t_ex, p_ex, rssat_ex);
363 const auto& densityLiquidEx =
364 fsys.phaseIsActive(fsys.waterPhaseIdx)
365 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
366 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
368 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
369 const auto rho_sat_ex = bLiquidSatEx *
371 rssat_ex * fsys.referenceDensity(fsys.gasPhaseIdx,
372 intQuantsEx.pvtRegionIndex()));
374 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
375 const auto pressure_difference_convective_mixing = delta_rho * distZg;
378 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
380 short interiorDofIdx = 0;
381 constexpr short exteriorDofIdx = 1;
384 if (pressure_difference_convective_mixing > 0) {
385 upIdx = exteriorDofIdx;
388 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
389 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
390 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
392 fsys.phaseIsActive(fsys.waterPhaseIdx)
393 ? up.fluidState().Rsw()
394 : up.fluidState().Rs();
396 const Evaluation& transMult = up.rockCompTransMultiplier();
397 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
398 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
401 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
403 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
404 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
405 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(fsys.gasCompIdx);
406 if (globalUpIndex == globalIndexIn) {
407 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
410 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
413 if constexpr (enableFullyImplicitThermal) {
414 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
415 fsys.referenceDensity(fsys.gasPhaseIdx, up.pvtRegionIndex());
416 if (globalUpIndex == globalIndexIn) {
417 flux[contiEnergyEqIdx] += convectiveFlux * h;
420 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
427template <
class TypeTag,
bool enableConvectiveMixingV>
437template <
class TypeTag>
451 const auto liquidPhaseIdx =
452 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
453 ? FluidSystem::waterPhaseIdx
454 : FluidSystem::oilPhaseIdx;
456 const Evaluation SoMax = 0.0;
457 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
459 asImp_().pvtRegionIndex(),
464 {
return saturatedDissolutionFactor_; }
468 {
return *
static_cast<Implementation*
>(
this); }
473template <
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:449
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:467
OPM_HOST_DEVICE const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:463
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:470
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
Definition: blackoilconvectivemixingmodule.hh:428
static void addConvectiveMixingFlux(RateVector &, const Context &, unsigned, unsigned)
Definition: blackoilconvectivemixingmodule.hh:132
static bool active(const Context &)
Definition: blackoilconvectivemixingmodule.hh:121
static void beginEpisode(const EclipseState &, const Schedule &, const int, ConvectiveMixingModuleParam< Scalar > &)
Definition: blackoilconvectivemixingmodule.hh:114
static void modifyAvgDensity(Evaluation &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned int, const ConvectiveMixingModuleParam< Scalar > &)
Definition: blackoilconvectivemixingmodule.hh:124
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:142
static void beginEpisode(const EclipseState &eclState, const Schedule &schedule, const int episodeIdx, ConvectiveMixingModuleParamT &info)
Definition: blackoilconvectivemixingmodule.hh:175
static OPM_HOST_DEVICE void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const CMMParam &info)
Definition: blackoilconvectivemixingmodule.hh:198
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:302
static OPM_HOST_DEVICE void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:262
Definition: blackoilconvectivemixingmodule.hh:99
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:53
Storage< Scalar > Xhi_
Definition: blackoilconvectivemixingmodule.hh:55
Storage< bool > active_
Definition: blackoilconvectivemixingmodule.hh:54
Storage< Scalar > Psi_
Definition: blackoilconvectivemixingmodule.hh:56