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/common/utility/gpuistl_if_available.hpp>
34#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
35
36#include <opm/input/eclipse/Schedule/OilVaporizationProperties.hpp>
37#include <opm/input/eclipse/Schedule/Schedule.hpp>
38
39#include <opm/material/common/MathToolbox.hpp>
40#include <opm/material/common/Valgrind.hpp>
41#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
42
46
48
50
51#include <cstddef>
52
53namespace Opm {
54
67template <class TypeTag>
68class BlackOilConvectiveMixingModule<TypeTag, /*enableConvectiveMixing=*/true>
69{
77 using Toolbox = MathToolbox<Evaluation>;
79
80 static constexpr int conti0EqIdx = Indices::conti0EqIdx;
81 static constexpr int dimWorld = GridView::dimensionworld;
82 static constexpr bool enableFullyImplicitThermal =
83 (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
84 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
85
86public:
87 static void beginEpisode(const EclipseState& eclState,
88 const Schedule& schedule,
89 const int episodeIdx,
91 {
92 // check that Xhi and Psi didn't change
93 std::size_t numRegions = eclState.runspec().tabdims().getNumPVTTables();
94 const auto& control = schedule[episodeIdx].oilvap();
95 if (info.active_.empty()) {
96 info.active_.resize(numRegions);
97 info.Xhi_.resize(numRegions);
98 info.Psi_.resize(numRegions);
99 }
100 for (std::size_t i = 0; i < numRegions; ++i ) {
101 info.active_[i] = control.drsdtConvective(i);
102 if (control.drsdtConvective(i)) {
103 info.Xhi_[i] = control.getMaxDRSDT(i);
104 info.Psi_[i] = control.getPsi(i);
105 }
106 }
107 }
108
109 template <class CMMParam>
110 OPM_HOST_DEVICE static void modifyAvgDensity(Evaluation& rhoAvg,
111 const IntensiveQuantities& intQuantsIn,
112 const IntensiveQuantities& intQuantsEx,
113 const unsigned phaseIdx,
114 const CMMParam& info) {
115 // Local constexpr aliases avoid nvcc's failure to resolve class-scope
116 // static constexpr members in device code.
117 constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
118 constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
119
120 if (info.active_.empty()) {
121 return;
122 }
123 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
124 return;
125 }
126
127 FluidSystem fsys = intQuantsIn.getFluidSystem();
128 if (phaseIdx == fsys.gasPhaseIdx) {
129 return;
130 }
131
132 const auto& liquidPhaseIdx =
133 fsys.phaseIsActive(waterPhaseIdx)
134 ? waterPhaseIdx
135 : oilPhaseIdx;
136
137 // Compute avg density based on pure water
138 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
139 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
140 const auto& salt_in = intQuantsIn.fluidState().saltConcentration();
141
142 const auto& bLiquidIn =
143 fsys.phaseIsActive(waterPhaseIdx)
144 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
145 t_in, p_in, Evaluation(0.0), salt_in)
146 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
147 t_in, p_in, Evaluation(0.0));
148
149 const auto& refDensityLiquidIn =
150 fsys.phaseIsActive(waterPhaseIdx)
151 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
152 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
153
154 const auto& rho_in = bLiquidIn * refDensityLiquidIn;
155
156 const auto t_ex = Toolbox::value(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
157 const auto p_ex = Toolbox::value(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
158 const auto salt_ex = Toolbox::value(intQuantsEx.fluidState().saltConcentration());
159
160 const auto bLiquidEx =
161 fsys.phaseIsActive(waterPhaseIdx)
162 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
163 t_ex, p_ex, Scalar{0.0}, salt_ex)
164 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
165 t_ex, p_ex, Scalar{0.0});
166
167 const auto& refDensityLiquidEx =
168 fsys.phaseIsActive(waterPhaseIdx)
169 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
170 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
171
172 const auto rho_ex = bLiquidEx * refDensityLiquidEx;
173
174 rhoAvg = (rho_in + rho_ex) / 2;
175 }
176
177 template <class Context>
178 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVector& flux,
179 const Context& elemCtx,
180 unsigned scvfIdx,
181 unsigned timeIdx)
182 {
183 // need for darcy flux calculation
184 const auto& problem = elemCtx.problem();
185 const auto& stencil = elemCtx.stencil(timeIdx);
186 const auto& scvf = stencil.interiorFace(scvfIdx);
187
188 unsigned interiorDofIdx = scvf.interiorIndex();
189 unsigned exteriorDofIdx = scvf.exteriorIndex();
190 assert(interiorDofIdx != exteriorDofIdx);
191
192 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
193 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
194 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
195 Scalar faceArea = scvf.area();
196 const Scalar g = problem.gravity()[dimWorld - 1];
197 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
198 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
199 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
200 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
201 const Scalar distZ = zIn - zEx;
202 addConvectiveMixingFlux(flux,
203 intQuantsIn,
204 intQuantsEx,
205 globalIndexIn,
206 globalIndexEx,
207 distZ * g,
208 trans,
209 faceArea,
210 problem.moduleParams().convectiveMixingModuleParam);
211 }
212
217 template <class RateVectorT = RateVector, class CMMParam = ConvectiveMixingModuleParamT>
218 OPM_HOST_DEVICE static void addConvectiveMixingFlux(RateVectorT& flux,
219 const IntensiveQuantities& intQuantsIn,
220 const IntensiveQuantities& intQuantsEx,
221 const unsigned globalIndexIn,
222 const unsigned globalIndexEx,
223 const Scalar distZg,
224 const Scalar trans,
225 const Scalar faceArea,
226 const CMMParam& info)
227 {
228 // Local constexpr aliases avoid nvcc's failure to resolve class-scope
229 // static constexpr members in device code.
230 constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
231 constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
232
233 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
234
235 if (info.active_.empty()) {
236 return;
237 }
238
239 if (!info.active_[intQuantsIn.pvtRegionIndex()] || !info.active_[intQuantsEx.pvtRegionIndex()]) {
240 return;
241 }
242
243 const auto& liquidPhaseIdx =
244 fsys.phaseIsActive(waterPhaseIdx)
245 ? waterPhaseIdx
246 : oilPhaseIdx;
247
248 // interiour
249 const auto& t_in = intQuantsIn.fluidState().temperature(liquidPhaseIdx);
250 const auto& p_in = intQuantsIn.fluidState().pressure(liquidPhaseIdx);
251 const auto& rssat_in = intQuantsIn.saturatedDissolutionFactor();
252 const auto& salt_in = intQuantsIn.fluidState().saltSaturation();
253
254 const auto bLiquidSatIn =
255 fsys.phaseIsActive(waterPhaseIdx)
256 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
257 t_in, p_in, rssat_in, salt_in)
258 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsIn.pvtRegionIndex(),
259 t_in, p_in, rssat_in);
260
261 const auto& densityLiquidIn =
262 fsys.phaseIsActive(waterPhaseIdx)
263 ? fsys.waterPvt().waterReferenceDensity(intQuantsIn.pvtRegionIndex())
264 : fsys.oilPvt().oilReferenceDensity(intQuantsIn.pvtRegionIndex());
265
266 const auto rho_in = Opm::getValue(intQuantsIn.fluidState().invB(liquidPhaseIdx)) * densityLiquidIn;
267 const auto rho_sat_in = bLiquidSatIn *
268 (densityLiquidIn +
269 rssat_in * fsys.referenceDensity(fsys.gasPhaseIdx,
270 intQuantsIn.pvtRegionIndex()));
271
272 // exteriour
273 const auto t_ex = Opm::getValue(intQuantsEx.fluidState().temperature(liquidPhaseIdx));
274 const auto p_ex = Opm::getValue(intQuantsEx.fluidState().pressure(liquidPhaseIdx));
275 const auto rssat_ex = Opm::getValue(intQuantsEx.saturatedDissolutionFactor());
276 const auto salt_ex = Opm::getValue(intQuantsEx.fluidState().saltSaturation());
277 const auto bLiquidSatEx =
278 fsys.phaseIsActive(waterPhaseIdx)
279 ? fsys.waterPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
280 t_ex, p_ex, rssat_ex, salt_ex)
281 : fsys.oilPvt().inverseFormationVolumeFactor(intQuantsEx.pvtRegionIndex(),
282 t_ex, p_ex, rssat_ex);
283
284 const auto& densityLiquidEx =
285 fsys.phaseIsActive(waterPhaseIdx)
286 ? fsys.waterPvt().waterReferenceDensity(intQuantsEx.pvtRegionIndex())
287 : fsys.oilPvt().oilReferenceDensity(intQuantsEx.pvtRegionIndex());
288
289 const auto rho_ex = Opm::getValue(intQuantsEx.fluidState().invB(liquidPhaseIdx)) * densityLiquidEx;
290 const auto rho_sat_ex = bLiquidSatEx *
291 (densityLiquidEx +
292 rssat_ex * fsys.referenceDensity(fsys.gasPhaseIdx,
293 intQuantsEx.pvtRegionIndex()));
294 // rho difference approximation
295 const auto delta_rho = (rho_sat_ex + rho_sat_in - rho_in - rho_ex) / 2;
296 const auto pressure_difference_convective_mixing = delta_rho * distZg;
297
298 // if change in pressure
299 if (Opm::abs(pressure_difference_convective_mixing) > 1e-12) {
300 // find new upstream direction
301 short interiorDofIdx = 0;
302 constexpr short exteriorDofIdx = 1;
303 short upIdx = 0;
304
305 if (pressure_difference_convective_mixing > 0) {
306 upIdx = exteriorDofIdx;
307 }
308
309 const auto& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
310 const auto& rssat_up = (upIdx == interiorDofIdx) ? rssat_in : rssat_ex;
311 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
312 const auto& Rsup =
313 fsys.phaseIsActive(waterPhaseIdx)
314 ? up.fluidState().Rsw()
315 : up.fluidState().Rs();
316
317 const Evaluation& transMult = up.rockCompTransMultiplier();
318 const auto& invB = up.fluidState().invB(liquidPhaseIdx);
319 const auto& visc = up.fluidState().viscosity(liquidPhaseIdx);
320
321 // We restrict the convective mixing mass flux to rssat * Psi.
322 const Evaluation RsupRestricted = Opm::min(Rsup, rssat_up*info.Psi_[up.pvtRegionIndex()]);
323
324 const auto convectiveFlux = -trans * transMult * info.Xhi_[up.pvtRegionIndex()] * invB *
325 pressure_difference_convective_mixing * RsupRestricted / (visc * faceArea);
326 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(fsys.gasCompIdx);
327 if (globalUpIndex == globalIndexIn) {
328 flux[conti0EqIdx + activeGasCompIdx] += convectiveFlux;
329 }
330 else {
331 flux[conti0EqIdx + activeGasCompIdx] += Opm::getValue(convectiveFlux);
332 }
333
334 if constexpr (enableFullyImplicitThermal) {
335 const auto& h = up.fluidState().enthalpy(liquidPhaseIdx) *
336 fsys.referenceDensity(fsys.gasPhaseIdx, up.pvtRegionIndex());
337 if (globalUpIndex == globalIndexIn) {
338 flux[contiEnergyEqIdx] += convectiveFlux * h;
339 }
340 else {
341 flux[contiEnergyEqIdx] += Opm::getValue(h) * Opm::getValue(convectiveFlux);
342 }
343 }
344 }
345 }
346};
347
355template <class TypeTag>
356class BlackOilConvectiveMixingIntensiveQuantities<TypeTag, /*enableConvectiveMixingV=*/true>
357{
361
362public:
368 {
369 const auto liquidPhaseIdx =
370 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
371 ? FluidSystem::waterPhaseIdx
372 : FluidSystem::oilPhaseIdx;
373
374 const Evaluation SoMax = 0.0;
375 saturatedDissolutionFactor_ = FluidSystem::saturatedDissolutionFactor(asImp_().fluidState(),
376 liquidPhaseIdx,
377 asImp_().pvtRegionIndex(),
378 SoMax);
379 }
380
381 OPM_HOST_DEVICE const Evaluation& saturatedDissolutionFactor() const
382 { return saturatedDissolutionFactor_; }
383
384protected:
385 Implementation& asImp_()
386 { return *static_cast<Implementation*>(this); }
387
389};
390
391}
392
393#endif
Classes required for dynamic convective mixing.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
void updateSaturatedDissolutionFactor_()
Compute the intensive quantities needed to handle convective dissolution.
Definition: blackoilconvectivemixingmodule.hh:367
Implementation & asImp_()
Definition: blackoilconvectivemixingmodule.hh:385
OPM_HOST_DEVICE const Evaluation & saturatedDissolutionFactor() const
Definition: blackoilconvectivemixingmodule.hh:381
Evaluation saturatedDissolutionFactor_
Definition: blackoilconvectivemixingmodule.hh:388
Provides the volumetric quantities required for the equations needed by the convective mixing (DRSDTC...
static void beginEpisode(const EclipseState &eclState, const Schedule &schedule, const int episodeIdx, ConvectiveMixingModuleParamT &info)
Definition: blackoilconvectivemixingmodule.hh:87
static OPM_HOST_DEVICE void modifyAvgDensity(Evaluation &rhoAvg, const IntensiveQuantities &intQuantsIn, const IntensiveQuantities &intQuantsEx, const unsigned phaseIdx, const CMMParam &info)
Definition: blackoilconvectivemixingmodule.hh:110
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:218
static OPM_HOST_DEVICE void addConvectiveMixingFlux(RateVector &flux, const Context &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilconvectivemixingmodule.hh:178
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: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: blackoilconvectivemixingmoduleparam.hpp:47
Storage< Scalar > Xhi_
Definition: blackoilconvectivemixingmoduleparam.hpp:49
Storage< bool > active_
Definition: blackoilconvectivemixingmoduleparam.hpp:48
Storage< Scalar > Psi_
Definition: blackoilconvectivemixingmoduleparam.hpp:50