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#if HAVE_ECL_INPUT
44#include <cstddef>
45#endif
46
47#include <vector>
48
49namespace Opm {
50
64template <class TypeTag, bool enableConvectiveMixing>
66
71template <class TypeTag>
72class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/false>
73{
78
79public:
80 struct ConvectiveMixingModuleParam
81 {};
82
83 #if HAVE_ECL_INPUT
84 static void beginEpisode(const EclipseState&,
85 const Schedule&,
86 const int,
87 ConvectiveMixingModuleParam&)
88 {}
89 #endif
90
91 template <class Context>
92 static bool active(const Context&)
93 { return false; }
94
95 static void modifyAvgDensity(Evaluation&,
96 const IntensiveQuantities&,
97 const IntensiveQuantities&,
98 const unsigned int,
99 const ConvectiveMixingModuleParam&)
100 {}
101
102 template <class Context>
103 static void addConvectiveMixingFlux(RateVector&,
104 const Context&,
105 unsigned,
106 unsigned)
107 {}
108
113 static void addConvectiveMixingFlux(RateVector&,
114 const IntensiveQuantities&,
115 const IntensiveQuantities&,
116 const unsigned,
117 const unsigned,
118 const Scalar,
119 const Scalar,
120 const Scalar,
121 const ConvectiveMixingModuleParam&)
122 {}
123};
124
125template <class TypeTag>
126class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/true>
127{
135 using Toolbox = MathToolbox<Evaluation>;
136
137 enum { conti0EqIdx = Indices::conti0EqIdx };
138 enum { dimWorld = GridView::dimensionworld };
139 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
140 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
141 static constexpr bool enableEnergy = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
142 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
143
144public:
145 struct ConvectiveMixingModuleParam
146 {
147 std::vector<bool> active_;
148 std::vector<Scalar> Xhi_;
149 std::vector<Scalar> Psi_;
150 };
151
152 #if HAVE_ECL_INPUT
153 static void beginEpisode(const EclipseState& eclState,
154 const Schedule& schedule,
155 const int episodeIdx,
156 ConvectiveMixingModuleParam& info)
157 {
158 // check that Xhi and Psi didn't change
159 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
160 const auto& control = schedule[episodeIdx].oilvap();
161 if (info.active_.empty()) {
162 info.active_.resize(numRegions);
163 info.Xhi_.resize(numRegions);
164 info.Psi_.resize(numRegions);
165 }
166 for (std::size_t i = 0; i < numRegions; ++i ) {
167 info.active_[i] = control.drsdtConvective(i);
168 if (control.drsdtConvective(i)) {
169 info.Xhi_[i] = control.getMaxDRSDT(i);
170 info.Psi_[i] = control.getPsi(i);
171 }
172 }
173 }
174 #endif
175
176 static void modifyAvgDensity(Evaluation& rhoAvg,
177 const IntensiveQuantities& intQuantsIn,
178 const IntensiveQuantities& intQuantsEx,
179 const unsigned phaseIdx,
180 const ConvectiveMixingModuleParam& info) {
181
182 if (info.active_.empty()) {
183 return;
184 }
185 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
186 return;
187 }
188
189 if (phaseIdx == FluidSystem::gasPhaseIdx) {
190 return;
191 }
192
193 const auto& liquidPhaseIdx =
194 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
195 ? FluidSystem::waterPhaseIdx
196 : FluidSystem::oilPhaseIdx;
197
198 // Compute avg density based on pure water
199 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
200 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
201 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
202
203 const auto& bLiquidIn =
204 FluidSystem::phaseIsActive(waterPhaseIdx)
205 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
206 t_in, p_in, Evaluation(0.0), salt_in)
207 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
208 t_in, p_in, Evaluation(0.0));
209
210 const auto& refDensityLiquidIn =
211 FluidSystem::phaseIsActive(waterPhaseIdx)
212 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
213 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
214
215 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
216
217 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
218 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
219 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
220
221 const auto bLiquidEx =
222 FluidSystem::phaseIsActive(waterPhaseIdx)
223 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
224 t_ex, p_ex, Scalar{0.0}, salt_ex)
225 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
226 t_ex, p_ex, Scalar{0.0});
227
228 const auto& refDensityLiquidEx =
229 FluidSystem::phaseIsActive(waterPhaseIdx)
230 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
231 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
232
233 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
234
235 rhoAvg = (rho_in + rho_ex) / 2;
236 }
237
238 template <class Context>
239 static void addConvectiveMixingFlux(RateVector& flux,
240 const Context& elemCtx,
241 unsigned scvfIdx,
242 unsigned timeIdx)
243 {
244 // need for darcy flux calculation
245 const auto& problem = elemCtx.problem();
246 const auto& stencil = elemCtx.stencil(timeIdx);
247 const auto& scvf = stencil.interiorFace(scvfIdx);
248
249 unsigned interiorDofIdx = scvf.interiorIndex();
250 unsigned exteriorDofIdx = scvf.exteriorIndex();
251 assert(interiorDofIdx != exteriorDofIdx);
252
253 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
254 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
255 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
256 Scalar faceArea = scvf.area();
257 const Scalar g = problem.gravity()[dimWorld - 1];
258 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
259 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
260 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
261 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
262 const Scalar distZ = zIn - zEx;
263 addConvectiveMixingFlux(flux,
264 intQuantsIn,
265 intQuantsEx,
266 globalIndexIn,
267 globalIndexEx,
268 distZ * g,
269 trans,
270 faceArea,
271 problem.moduleParams().convectiveMixingModuleParam);
272 }
273
278 static void addConvectiveMixingFlux(RateVector& flux,
279 const IntensiveQuantities& intQuantsIn,
280 const IntensiveQuantities& intQuantsEx,
281 const unsigned globalIndexIn,
282 const unsigned globalIndexEx,
283 const Scalar distZg,
284 const Scalar trans,
285 const Scalar faceArea,
286 const ConvectiveMixingModuleParam& info)
287 {
288 if (info.active_.empty()) {
289 return;
290 }
291
292 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
293 return;
294 }
295
296 const auto& liquidPhaseIdx =
297 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
298 ? FluidSystem::waterPhaseIdx
299 : FluidSystem::oilPhaseIdx;
300
301 // interiour
302 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
303 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
304 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
305 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
306
307 const auto bLiquidSatIn =
308 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
309 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
310 t_in, p_in, rssat_in, salt_in)
311 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
312 t_in, p_in, rssat_in);
313
314 const auto& densityLiquidIn =
315 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
316 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
317 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
318
319 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
320 const auto rho_sat_in = bLiquidSatIn *
321 (densityLiquidIn +
322 rssat_in * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
323 intQuantsIn.pvtRegionIndex()));
324
325 // exteriour
326 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
327 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
328 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
329 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
330 const auto bLiquidSatEx =
331 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
332 ? FluidSystem::waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
333 t_ex, p_ex, rssat_ex, salt_ex)
334 : FluidSystem::oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
335 t_ex, p_ex, rssat_ex);
336
337 const auto& densityLiquidEx =
338 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
339 ? FluidSystem::waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
340 : FluidSystem::oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
341
342 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
343 const auto rho_sat_ex = bLiquidSatEx *
344 (densityLiquidEx +
345 rssat_ex * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx,
346 intQuantsEx.pvtRegionIndex()));
347 // rho difference approximation
348 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
349 const auto pressure_difference_convective_mixing = delta_rho * distZg;
350
351 // if change in pressure
352 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
353 // find new upstream direction
354 short interiorDofIdx = 0;
355 constexpr short exteriorDofIdx = 1;
356 short upIdx = 0;
357
358 if (pressure_difference_convective_mixing > 0) {
359 upIdx = exteriorDofIdx;
360 }
361
362 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
363 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
364 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
365 const auto& Rsup =
366 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
367 ? up.fluidState().Rsw()
368 : up.fluidState().Rs();
369
370 const Evaluation& transMult = up.rockCompTransMultiplier();
371 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
372 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
373
374 // We restrict the convective mixing mass flux to rssat * Psi.
375 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
376
377 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
378 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
379 unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
380 if (globalUpIndex == globalIndexIn) {
381 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
382 }
383 else {
384 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
385 }
386
387 if constexpr (enableEnergy) {
388 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
389 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, up.pvtRegionIndex());
390 if (globalUpIndex == globalIndexIn) {
391 flux[contiEnergyEqIdx] += convectiveFlux * h;
392 }
393 else {
394 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
395 }
396 }
397 }
398 }
399};
400
401template <class TypeTag, bool enableConvectiveMixingV>
403
411template <class TypeTag>
412class BlackOilConvectiveMixingIntensiveQuantities<TypeTag, /*enableConvectiveMixingV=*/true>
413{
417
418public:
424 {
425 const auto liquidPhaseIdx =
426 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
427 ? FluidSystem::waterPhaseIdx
428 : FluidSystem::oilPhaseIdx;
429
430 const Evaluation SoMax = 0.0;
431 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
432 liquidPhaseIdx,
433 asImp_().pvtRegionIndex(),
434 SoMax);
435 }
436
437 const Evaluation& saturatedDissolutionFactor() const
438 { return saturatedDissolutionFactor_; }
439
440protected:
441 Implementation& asImp_()
442 { return *static_cast<Implementation*>(this); }
443
445};
446
447template <class TypeTag>
449{
450};
451
452}
453
454#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:423
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:441
const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:437
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:444
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
Definition: blackoilconvectivemixingmodule.hh:402
static void addConvectiveMixingFlux(RateVector &, const Context &, unsigned, unsigned)
Definition: blackoilconvectivemixingmodule.hh:103
static bool active(const Context &)
Definition: blackoilconvectivemixingmodule.hh:92
static void modifyAvgDensity(Evaluation &, const IntensiveQuantities &, const IntensiveQuantities &, const unsigned int, const ConvectiveMixingModuleParam &)
Definition: blackoilconvectivemixingmodule.hh:95
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:113
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:278
static void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:239
static void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const ConvectiveMixingModuleParam &info)
Definition: blackoilconvectivemixingmodule.hh:176
Definition: blackoilconvectivemixingmodule.hh:65
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: blackoilbioeffectsmodules.hh:43
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:148
std::vector< Scalar > Psi_
Definition: blackoilconvectivemixingmodule.hh:149
std::vector< bool > active_
Definition: blackoilconvectivemixingmodule.hh:147