blackoilconvectivemixingmodule.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_CONVECTIVEMIXING_MODULE_HH
29#define EWOMS_CONVECTIVEMIXING_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/input/eclipse/Schedule/OilVaporizationProperties.hpp>
34#include <opm/input/eclipse/Schedule/Schedule.hpp>
35
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/common/Valgrind.hpp>
38
42
43#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
44
45#include <opm/common/utility/gpuistl_if_available.hpp>
46
47#if HAVE_ECL_INPUT
48#include <cstddef>
49#endif
50
51#include <vector>
52
53namespace Opm {
54
55 template<class Scalar, template<class> class Storage = Opm::VectorWithDefaultAllocator>
57 {
58 Storage<bool> active_;
59 Storage<Scalar> Xhi_;
60 Storage<Scalar> Psi_;
61 };
62
63#ifdef HAVE_CUDA
64namespace gpuistl
65{
66 template <class Scalar>
68 {
70 view.active_ = gpuistl::make_view(params.active_);
71 view.Xhi_ = gpuistl::make_view(params.Xhi_);
72 view.Psi_ = gpuistl::make_view(params.Psi_);
73 return view;
74 }
75
76 template <class Scalar>
77 ConvectiveMixingModuleParam<Scalar, GpuBuffer> copy_to_gpu(const ConvectiveMixingModuleParam<Scalar, Opm::VectorWithDefaultAllocator>& params)
78 {
79 return ConvectiveMixingModuleParam<Scalar, GpuBuffer>{
80 gpuistl::GpuBuffer(params.active_),
81 gpuistl::GpuBuffer(params.Xhi_),
82 gpuistl::GpuBuffer(params.Psi_)
83 };
84 }
85}
86
87#endif
88
102template <class TypeTag, bool enableConvectiveMixing>
104
109template <class TypeTag>
110class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/false>
111{
116
117public:
118
119 #if HAVE_ECL_INPUT
120 static void beginEpisode(const EclipseState&,
121 const Schedule&,
122 const int,
124 {}
125 #endif
126
127 template <class Context>
128 static bool active(const Context&)
129 { return false; }
130
131 static void modifyAvgDensity(Evaluation&,
132 const IntensiveQuantities&,
133 const IntensiveQuantities&,
134 const unsigned int,
136 {}
137
138 template <class Context>
139 static void addConvectiveMixingFlux(RateVector&,
140 const Context&,
141 unsigned,
142 unsigned)
143 {}
144
149 static void addConvectiveMixingFlux(RateVector&,
150 const IntensiveQuantities&,
151 const IntensiveQuantities&,
152 const unsigned,
153 const unsigned,
154 const Scalar,
155 const Scalar,
156 const Scalar,
158 {}
159};
160
161template <class TypeTag>
162class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/true>
163{
171 using Toolbox = MathToolbox<Evaluation>;
173
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;
180
181public:
182 #if HAVE_ECL_INPUT
183 static void beginEpisode(const EclipseState& eclState,
184 const Schedule& schedule,
185 const int episodeIdx,
187 {
188 // check that Xhi and Psi didn't change
189 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
190 const auto& control = schedule[episodeIdx].oilvap();
191 if (info.active_.empty()) {
192 info.active_.resize(numRegions);
193 info.Xhi_.resize(numRegions);
194 info.Psi_.resize(numRegions);
195 }
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);
201 }
202 }
203 }
204 #endif
205
206 template <class CMMParam>
207 OPM_HOST_DEVICE static void modifyAvgDensity(Evaluation& rhoAvg,
208 const IntensiveQuantities& intQuantsIn,
209 const IntensiveQuantities& intQuantsEx,
210 const unsigned phaseIdx,
211 const CMMParam& info) {
212
213 if (info.active_.empty()) {
214 return;
215 }
216 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
217 return;
218 }
219
220 FluidSystem fsys = intQuantsIn.getFluidSystem();
221 if (phaseIdx == fsys.gasPhaseIdx) {
222 return;
223 }
224
225 const auto& liquidPhaseIdx =
226 fsys.phaseIsActive(fsys.waterPhaseIdx)
227 ? fsys.waterPhaseIdx
228 : fsys.oilPhaseIdx;
229
230 // Compute avg density based on pure water
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();
234
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));
241
242 const auto& refDensityLiquidIn =
243 fsys.phaseIsActive(waterPhaseIdx)
244 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
245 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
246
247 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
248
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());
252
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});
259
260 const auto& refDensityLiquidEx =
261 fsys.phaseIsActive(waterPhaseIdx)
262 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
263 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
264
265 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
266
267 rhoAvg = (rho_in + rho_ex) / 2;
268 }
269
270 template <class Context>
271 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVector& flux,
272 const Context& elemCtx,
273 unsigned scvfIdx,
274 unsigned timeIdx)
275 {
276 // need for darcy flux calculation
277 const auto& problem = elemCtx.problem();
278 const auto& stencil = elemCtx.stencil(timeIdx);
279 const auto& scvf = stencil.interiorFace(scvfIdx);
280
281 unsigned interiorDofIdx = scvf.interiorIndex();
282 unsigned exteriorDofIdx = scvf.exteriorIndex();
283 assert(interiorDofIdx != exteriorDofIdx);
284
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,
296 intQuantsIn,
297 intQuantsEx,
298 globalIndexIn,
299 globalIndexEx,
300 distZ * g,
301 trans,
302 faceArea,
303 problem.moduleParams().convectiveMixingModuleParam);
304 }
305
310 template <class RateVectorT = RateVector, class CMMParam = ConvectiveMixingModuleParamT>
311 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVectorT& flux,
312 const IntensiveQuantities& intQuantsIn,
313 const IntensiveQuantities& intQuantsEx,
314 const unsigned globalIndexIn,
315 const unsigned globalIndexEx,
316 const Scalar distZg,
317 const Scalar trans,
318 const Scalar faceArea,
319 const CMMParam& info)
320 {
321 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
322
323 if (info.active_.empty()) {
324 return;
325 }
326
327 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
328 return;
329 }
330
331 const auto& liquidPhaseIdx =
332 fsys.phaseIsActive(fsys.waterPhaseIdx)
333 ? fsys.waterPhaseIdx
334 : fsys.oilPhaseIdx;
335
336 // interiour
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();
341
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);
348
349 const auto& densityLiquidIn =
350 fsys.phaseIsActive(fsys.waterPhaseIdx)
351 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
352 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
353
354 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
355 const auto rho_sat_in = bLiquidSatIn *
356 (densityLiquidIn +
357 rssat_in * fsys.referenceDensity(fsys.gasPhaseIdx,
358 intQuantsIn.pvtRegionIndex()));
359
360 // exteriour
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);
371
372 const auto& densityLiquidEx =
373 fsys.phaseIsActive(fsys.waterPhaseIdx)
374 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
375 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
376
377 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
378 const auto rho_sat_ex = bLiquidSatEx *
379 (densityLiquidEx +
380 rssat_ex * fsys.referenceDensity(fsys.gasPhaseIdx,
381 intQuantsEx.pvtRegionIndex()));
382 // rho difference approximation
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;
385
386 // if change in pressure
387 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
388 // find new upstream direction
389 short interiorDofIdx = 0;
390 constexpr short exteriorDofIdx = 1;
391 short upIdx = 0;
392
393 if (pressure_difference_convective_mixing > 0) {
394 upIdx = exteriorDofIdx;
395 }
396
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;
400 const auto& Rsup =
401 fsys.phaseIsActive(fsys.waterPhaseIdx)
402 ? up.fluidState().Rsw()
403 : up.fluidState().Rs();
404
405 const Evaluation& transMult = up.rockCompTransMultiplier();
406 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
407 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
408
409 // We restrict the convective mixing mass flux to rssat * Psi.
410 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
411
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;
417 }
418 else {
419 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
420 }
421
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;
427 }
428 else {
429 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
430 }
431 }
432 }
433 }
434};
435
436template <class TypeTag, bool enableConvectiveMixingV>
438
446template <class TypeTag>
447class BlackOilConvectiveMixingIntensiveQuantities<TypeTag, /*enableConvectiveMixingV=*/true>
448{
452
453public:
459 {
460 const auto liquidPhaseIdx =
461 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
462 ? FluidSystem::waterPhaseIdx
463 : FluidSystem::oilPhaseIdx;
464
465 const Evaluation SoMax = 0.0;
466 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
467 liquidPhaseIdx,
468 asImp_().pvtRegionIndex(),
469 SoMax);
470 }
471
472 OPM_HOST_DEVICE const Evaluation& saturatedDissolutionFactor() const
473 { return saturatedDissolutionFactor_; }
474
475protected:
476 Implementation& asImp_()
477 { return *static_cast<Implementation*>(this); }
478
480};
481
482template <class TypeTag>
484{
485};
486
487}
488
489#endif
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