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#include <cstddef>
48
49namespace Opm {
50
51 template<class Scalar, template<class> class Storage = Opm::VectorWithDefaultAllocator>
53 {
54 Storage<bool> active_;
55 Storage<Scalar> Xhi_;
56 Storage<Scalar> Psi_;
57 };
58
59#ifdef HAVE_CUDA
60namespace gpuistl
61{
62 template <class Scalar>
64 {
66 view.active_ = gpuistl::make_view(params.active_);
67 view.Xhi_ = gpuistl::make_view(params.Xhi_);
68 view.Psi_ = gpuistl::make_view(params.Psi_);
69 return view;
70 }
71
72 template <class Scalar>
73 ConvectiveMixingModuleParam<Scalar, GpuBuffer> copy_to_gpu(const ConvectiveMixingModuleParam<Scalar, Opm::VectorWithDefaultAllocator>& params)
74 {
75 return ConvectiveMixingModuleParam<Scalar, GpuBuffer>{
76 gpuistl::GpuBuffer(params.active_),
77 gpuistl::GpuBuffer(params.Xhi_),
78 gpuistl::GpuBuffer(params.Psi_)
79 };
80 }
81}
82
83#endif
84
98template <class TypeTag, bool enableConvectiveMixing>
100
105template <class TypeTag>
106class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/false>
107{
112
113public:
114 static void beginEpisode(const EclipseState&,
115 const Schedule&,
116 const int,
118 {}
119
120 template <class Context>
121 static bool active(const Context&)
122 { return false; }
123
124 static void modifyAvgDensity(Evaluation&,
125 const IntensiveQuantities&,
126 const IntensiveQuantities&,
127 const unsigned int,
129 {}
130
131 template <class Context>
132 static void addConvectiveMixingFlux(RateVector&,
133 const Context&,
134 unsigned,
135 unsigned)
136 {}
137
142 static void addConvectiveMixingFlux(RateVector&,
143 const IntensiveQuantities&,
144 const IntensiveQuantities&,
145 const unsigned,
146 const unsigned,
147 const Scalar,
148 const Scalar,
149 const Scalar,
151 {}
152};
153
154template <class TypeTag>
155class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/true>
156{
164 using Toolbox = MathToolbox<Evaluation>;
166
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;
173
174public:
175 static void beginEpisode(const EclipseState& eclState,
176 const Schedule& schedule,
177 const int episodeIdx,
179 {
180 // check that Xhi and Psi didn't change
181 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
182 const auto& control = schedule[episodeIdx].oilvap();
183 if (info.active_.empty()) {
184 info.active_.resize(numRegions);
185 info.Xhi_.resize(numRegions);
186 info.Psi_.resize(numRegions);
187 }
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);
193 }
194 }
195 }
196
197 template <class CMMParam>
198 OPM_HOST_DEVICE static void modifyAvgDensity(Evaluation& rhoAvg,
199 const IntensiveQuantities& intQuantsIn,
200 const IntensiveQuantities& intQuantsEx,
201 const unsigned phaseIdx,
202 const CMMParam& info) {
203
204 if (info.active_.empty()) {
205 return;
206 }
207 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
208 return;
209 }
210
211 FluidSystem fsys = intQuantsIn.getFluidSystem();
212 if (phaseIdx == fsys.gasPhaseIdx) {
213 return;
214 }
215
216 const auto& liquidPhaseIdx =
217 fsys.phaseIsActive(fsys.waterPhaseIdx)
218 ? fsys.waterPhaseIdx
219 : fsys.oilPhaseIdx;
220
221 // Compute avg density based on pure water
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();
225
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));
232
233 const auto& refDensityLiquidIn =
234 fsys.phaseIsActive(waterPhaseIdx)
235 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
236 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
237
238 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
239
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());
243
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});
250
251 const auto& refDensityLiquidEx =
252 fsys.phaseIsActive(waterPhaseIdx)
253 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
254 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
255
256 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
257
258 rhoAvg = (rho_in + rho_ex) / 2;
259 }
260
261 template <class Context>
262 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVector& flux,
263 const Context& elemCtx,
264 unsigned scvfIdx,
265 unsigned timeIdx)
266 {
267 // need for darcy flux calculation
268 const auto& problem = elemCtx.problem();
269 const auto& stencil = elemCtx.stencil(timeIdx);
270 const auto& scvf = stencil.interiorFace(scvfIdx);
271
272 unsigned interiorDofIdx = scvf.interiorIndex();
273 unsigned exteriorDofIdx = scvf.exteriorIndex();
274 assert(interiorDofIdx != exteriorDofIdx);
275
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,
287 intQuantsIn,
288 intQuantsEx,
289 globalIndexIn,
290 globalIndexEx,
291 distZ * g,
292 trans,
293 faceArea,
294 problem.moduleParams().convectiveMixingModuleParam);
295 }
296
301 template <class RateVectorT = RateVector, class CMMParam = ConvectiveMixingModuleParamT>
302 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVectorT& flux,
303 const IntensiveQuantities& intQuantsIn,
304 const IntensiveQuantities& intQuantsEx,
305 const unsigned globalIndexIn,
306 const unsigned globalIndexEx,
307 const Scalar distZg,
308 const Scalar trans,
309 const Scalar faceArea,
310 const CMMParam& info)
311 {
312 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
313
314 if (info.active_.empty()) {
315 return;
316 }
317
318 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
319 return;
320 }
321
322 const auto& liquidPhaseIdx =
323 fsys.phaseIsActive(fsys.waterPhaseIdx)
324 ? fsys.waterPhaseIdx
325 : fsys.oilPhaseIdx;
326
327 // interiour
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();
332
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);
339
340 const auto& densityLiquidIn =
341 fsys.phaseIsActive(fsys.waterPhaseIdx)
342 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
343 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
344
345 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
346 const auto rho_sat_in = bLiquidSatIn *
347 (densityLiquidIn +
348 rssat_in * fsys.referenceDensity(fsys.gasPhaseIdx,
349 intQuantsIn.pvtRegionIndex()));
350
351 // exteriour
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);
362
363 const auto& densityLiquidEx =
364 fsys.phaseIsActive(fsys.waterPhaseIdx)
365 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
366 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
367
368 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
369 const auto rho_sat_ex = bLiquidSatEx *
370 (densityLiquidEx +
371 rssat_ex * fsys.referenceDensity(fsys.gasPhaseIdx,
372 intQuantsEx.pvtRegionIndex()));
373 // rho difference approximation
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;
376
377 // if change in pressure
378 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
379 // find new upstream direction
380 short interiorDofIdx = 0;
381 constexpr short exteriorDofIdx = 1;
382 short upIdx = 0;
383
384 if (pressure_difference_convective_mixing > 0) {
385 upIdx = exteriorDofIdx;
386 }
387
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;
391 const auto& Rsup =
392 fsys.phaseIsActive(fsys.waterPhaseIdx)
393 ? up.fluidState().Rsw()
394 : up.fluidState().Rs();
395
396 const Evaluation& transMult = up.rockCompTransMultiplier();
397 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
398 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
399
400 // We restrict the convective mixing mass flux to rssat * Psi.
401 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
402
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;
408 }
409 else {
410 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
411 }
412
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;
418 }
419 else {
420 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
421 }
422 }
423 }
424 }
425};
426
427template <class TypeTag, bool enableConvectiveMixingV>
429
437template <class TypeTag>
438class BlackOilConvectiveMixingIntensiveQuantities<TypeTag, /*enableConvectiveMixingV=*/true>
439{
443
444public:
450 {
451 const auto liquidPhaseIdx =
452 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
453 ? FluidSystem::waterPhaseIdx
454 : FluidSystem::oilPhaseIdx;
455
456 const Evaluation SoMax = 0.0;
457 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
458 liquidPhaseIdx,
459 asImp_().pvtRegionIndex(),
460 SoMax);
461 }
462
463 OPM_HOST_DEVICE const Evaluation& saturatedDissolutionFactor() const
464 { return saturatedDissolutionFactor_; }
465
466protected:
467 Implementation& asImp_()
468 { return *static_cast<Implementation*>(this); }
469
471};
472
473template <class TypeTag>
475{
476};
477
478}
479
480#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: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