blackoildispersionmodule.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_DISPERSION_MODULE_HH
29#define EWOMS_DISPERSION_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/material/common/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
35
40
41#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
42#include <opm/common/OpmLog/OpmLog.hpp>
43
44#include <array>
45#include <cmath>
46#include <stdexcept>
47
48namespace Opm {
49
56template <class TypeTag>
57class BlackOilDispersionModule<TypeTag, /*enableDispersion=*/true>
58{
69
70 enum { numPhases = FluidSystem::numPhases };
71 enum { numComponents = FluidSystem::numComponents };
72 enum { conti0EqIdx = Indices::conti0EqIdx };
73 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
74 static constexpr bool enableDispersion = true;
75 enum { enableMICP = Indices::enableMICP };
76
77 static constexpr unsigned contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
78 static constexpr unsigned contiOxygenEqIdx = Indices::contiOxygenEqIdx;
79 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
80 static constexpr unsigned contiUreaEqIdx = Indices::contiUreaEqIdx;
81
82 using Toolbox = MathToolbox<Evaluation>;
83
84public:
86
87 static void initFromState(const EclipseState& eclState)
88 {
89 if (!eclState.getSimulationConfig().rock_config().dispersion()) {
90 return;
91 }
92
93 if (eclState.getSimulationConfig().hasVAPWAT() || eclState.getSimulationConfig().hasVAPOIL()) {
94 OpmLog::warning("Dispersion is activated in combination with VAPWAT/VAPOIL. \n"
95 "Water/oil is still allowed to vaporize, but dispersion in the "
96 "gas phase is ignored.");
97 }
98 }
99
104 template <class Context>
105 static void addDispersiveFlux(RateVector& flux, const Context& context,
106 unsigned spaceIdx, unsigned timeIdx)
107 {
108 // Only work if dispersion is enabled by DISPERC in the deck
109 if (!context.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
110 return;
111 }
112 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
113 const auto& inIq = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
114 const auto& exIq = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
115 const auto& dispersivity = extQuants.dispersivity();
116 const auto& normVelocityAvg = extQuants.normVelocityAvg();
117 addDispersiveFlux(flux, inIq, exIq, dispersivity, normVelocityAvg);
118 }
119
140 template<class IntensiveQuantities, class Scalar>
141 static void addDispersiveFlux(RateVector& flux,
142 const IntensiveQuantities& inIq,
143 const IntensiveQuantities& exIq,
144 const Evaluation& dispersivity,
145 const Scalar& normVelocityAvg)
146 {
147 const auto& inFs = inIq.fluidState();
148 const auto& exFs = exIq.fluidState();
149 Evaluation diffR = 0.0;
150 if constexpr(enableBioeffects) {
151 // The dispersion coefficients are given for mass concentrations
152 const Evaluation bAvg = (inFs.invB(waterPhaseIdx) + Toolbox::value(exFs.invB(waterPhaseIdx))) / 2;
153 diffR = inIq.microbialConcentration() - Toolbox::value(exIq.microbialConcentration());
154 flux[contiMicrobialEqIdx] +=
155 bAvg *
156 normVelocityAvg[waterPhaseIdx] *
157 dispersivity *
158 diffR;
159 if constexpr(enableMICP) {
160 diffR = inIq.oxygenConcentration() - Toolbox::value(exIq.oxygenConcentration());
161 flux[contiOxygenEqIdx] +=
162 bAvg *
163 normVelocityAvg[waterPhaseIdx] *
164 dispersivity *
165 diffR;
166 diffR = inIq.ureaConcentration() - Toolbox::value(exIq.ureaConcentration());
167 flux[contiUreaEqIdx] +=
168 bAvg *
169 normVelocityAvg[waterPhaseIdx] *
170 dispersivity *
171 diffR;
172 return;
173 }
174 }
175
176 unsigned pvtRegionIndex = inFs.pvtRegionIndex();
177 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
178 if (!FluidSystem::phaseIsActive(phaseIdx)) {
179 continue;
180 }
181
182 // no dispersion in water for blackoil models unless water can contain dissolved gas
183 if (!FluidSystem::enableDissolvedGasInWater() && FluidSystem::waterPhaseIdx == phaseIdx) {
184 continue;
185 }
186
187 // adding dispersion in the gas phase leads to
188 // convergence issues and unphysical results.
189 // we disable dispersion in the gas phase for now
190 if (FluidSystem::gasPhaseIdx == phaseIdx) {
191 continue;
192 }
193
194 // arithmetic mean of the phase's b factor
195 Evaluation bAvg = inFs.invB(phaseIdx);
196 bAvg += Toolbox::value(exFs.invB(phaseIdx));
197 bAvg /= 2;
198
199 Evaluation convFactor = 1.0;
200 if (FluidSystem::enableDissolvedGas() &&
201 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
202 phaseIdx == FluidSystem::oilPhaseIdx)
203 {
204 const Evaluation rsAvg = (inFs.Rs() + Toolbox::value(exFs.Rs())) / 2;
205 convFactor = 1.0 / (toMassFractionGasOil(pvtRegionIndex) + rsAvg);
206 diffR = inFs.Rs() - Toolbox::value(exFs.Rs());
207 }
208 if (FluidSystem::enableVaporizedOil() &&
209 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
210 phaseIdx == FluidSystem::gasPhaseIdx)
211 {
212 const Evaluation rvAvg = (inFs.Rv() + Toolbox::value(exFs.Rv())) / 2;
213 convFactor = toMassFractionGasOil(pvtRegionIndex) /
214 (1.0 + rvAvg * toMassFractionGasOil(pvtRegionIndex));
215 diffR = inFs.Rv() - Toolbox::value(exFs.Rv());
216 }
217 if (FluidSystem::enableDissolvedGasInWater() && phaseIdx == FluidSystem::waterPhaseIdx) {
218 const Evaluation rsAvg = (inFs.Rsw() + Toolbox::value(exFs.Rsw())) / 2;
219 convFactor = 1.0 / (toMassFractionGasWater(pvtRegionIndex) + rsAvg);
220 diffR = inFs.Rsw() - Toolbox::value(exFs.Rsw());
221 }
222 if (FluidSystem::enableVaporizedWater() && phaseIdx == FluidSystem::gasPhaseIdx) {
223 const Evaluation rvAvg = (inFs.Rvw() + Toolbox::value(exFs.Rvw())) / 2;
224 convFactor = toMassFractionGasWater(pvtRegionIndex) /
225 (1.0 + rvAvg * toMassFractionGasWater(pvtRegionIndex));
226 diffR = inFs.Rvw() - Toolbox::value(exFs.Rvw());
227 }
228
229 // mass flux of solvent component
230 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
231 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
232 flux[conti0EqIdx + activeSolventCompIdx] +=
233 -bAvg *
234 normVelocityAvg[phaseIdx] *
235 convFactor *
236 dispersivity *
237 diffR;
238
239 // mass flux of solute component
240 const unsigned soluteCompIdx = FluidSystem::soluteComponentIndex(phaseIdx);
241 const unsigned activeSoluteCompIdx = FluidSystem::canonicalToActiveCompIdx(soluteCompIdx);
242 flux[conti0EqIdx + activeSoluteCompIdx] +=
243 bAvg *
244 normVelocityAvg[phaseIdx] *
245 convFactor *
246 dispersivity *
247 diffR;
248 }
249 }
250
251private:
252 static Scalar toMassFractionGasOil (unsigned regionIdx)
253 {
254 const Scalar rhoO = FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, regionIdx);
255 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
256 return rhoO / rhoG;
257 }
258
259 static Scalar toMassFractionGasWater (unsigned regionIdx)
260 {
261 const Scalar rhoW = FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, regionIdx);
262 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
263 return rhoW / rhoG;
264 }
265};
266
274template <class TypeTag>
275class BlackOilDispersionIntensiveQuantities<TypeTag, /*enableDispersion=*/true>
276{
281
282 enum { numPhases = FluidSystem::numPhases };
283 enum { numComponents = FluidSystem::numComponents };
284 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
285 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
286 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
287 enum { gasCompIdx = FluidSystem::gasCompIdx };
288 enum { oilCompIdx = FluidSystem::oilCompIdx };
289 enum { waterCompIdx = FluidSystem::waterCompIdx };
290 enum { conti0EqIdx = Indices::conti0EqIdx };
291 static constexpr bool enableDispersion = true;
292
293public:
297 Scalar normVelocityCell(unsigned phaseIdx) const
298 { return normVelocityCell_[phaseIdx]; }
299
300protected:
310 template<class ElementContext>
311 void update_(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
312 {
313 // Only work if dispersion is enabled by DISPERC in the deck
314 if (!elemCtx.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
315 return;
316 }
317 const auto& problem = elemCtx.simulator().problem();
318 if (problem.model().linearizer().getVelocityInfo().empty()) {
319 return;
320 }
321 const std::array<int, 3> phaseIdxs = {gasPhaseIdx, oilPhaseIdx, waterPhaseIdx};
322 const std::array<int, 3> compIdxs = {gasCompIdx, oilCompIdx, waterCompIdx};
323 const auto& velocityInf = problem.model().linearizer().getVelocityInfo();
324 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
325 const auto& velocityInfos = velocityInf[globalDofIdx];
326 for (unsigned i = 0; i < phaseIdxs.size(); ++i) {
327 normVelocityCell_[i] = 0;
328 }
329 for (const auto& velocityInfo : velocityInfos) {
330 for (unsigned i = 0; i < phaseIdxs.size(); ++i) {
331 if (FluidSystem::phaseIsActive(phaseIdxs[i])) {
332 normVelocityCell_[phaseIdxs[i]] = std::max(normVelocityCell_[phaseIdxs[i]],
333 std::abs(velocityInfo.velocity[conti0EqIdx +
334 FluidSystem::canonicalToActiveCompIdx(compIdxs[i])]));
335 }
336 }
337 }
338 }
339
340private:
341 std::array<Scalar, numPhases> normVelocityCell_{};
342};
343
350template <class TypeTag>
351class BlackOilDispersionExtensiveQuantities<TypeTag, /*enableDispersion=*/true>
352{
356
357 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
358
359public:
360 using ScalarArray = std::array<Scalar, numPhases>;
361 static void update(ScalarArray& normVelocityAvg,
362 const IntensiveQuantities& intQuantsInside,
363 const IntensiveQuantities& intQuantsOutside)
364 {
365 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
366 if (!FluidSystem::phaseIsActive(phaseIdx)) {
367 continue;
368 }
369 // use the arithmetic average for the effective
370 // velocity coefficients at the face's integration point.
371 normVelocityAvg[phaseIdx] =
372 0.5 * (intQuantsInside.normVelocityCell(phaseIdx) +
373 intQuantsOutside.normVelocityCell(phaseIdx));
374 Valgrind::CheckDefined(normVelocityAvg[phaseIdx]);
375 }
376 }
377
378protected:
379 template <class Context, class FluidState>
380 void updateBoundary_(const Context&,
381 unsigned,
382 unsigned,
383 const FluidState&)
384 { throw std::runtime_error("Not implemented: Dispersion across boundary not implemented for blackoil"); }
385
386public:
393 Scalar dispersivity() const
394 { return dispersivity_; }
395
403 Scalar normVelocityAvg(unsigned phaseIdx) const
404 { return normVelocityAvg_[phaseIdx]; }
405
406 const auto& normVelocityAvg() const
407 { return normVelocityAvg_; }
408
409private:
410 Scalar dispersivity_;
411 ScalarArray normVelocityAvg_;
412};
413
414} // namespace Opm
415
416#endif
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
Definition: blackoildispersionmodule.hh:352
Scalar dispersivity() const
The dispersivity of the face.
Definition: blackoildispersionmodule.hh:393
static void update(ScalarArray &normVelocityAvg, const IntensiveQuantities &intQuantsInside, const IntensiveQuantities &intQuantsOutside)
Definition: blackoildispersionmodule.hh:361
std::array< Scalar, numPhases > ScalarArray
Definition: blackoildispersionmodule.hh:360
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: blackoildispersionmodule.hh:380
const auto & normVelocityAvg() const
Definition: blackoildispersionmodule.hh:406
Scalar normVelocityAvg(unsigned phaseIdx) const
The effective velocity coefficient in a fluid phase at the face's integration point.
Definition: blackoildispersionmodule.hh:403
Provides the quantities required to calculate dispersive mass fluxes.
void update_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the quantities required to calculate dispersive mass fluxes. This considers the linear disperi...
Definition: blackoildispersionmodule.hh:311
Scalar normVelocityCell(unsigned phaseIdx) const
Returns the max. norm of the filter velocity of the cell.
Definition: blackoildispersionmodule.hh:297
Provides the volumetric quantities required for the calculation of dispersive fluxes.
static void initFromState(const EclipseState &eclState)
Definition: blackoildispersionmodule.hh:87
static void addDispersiveFlux(RateVector &flux, const IntensiveQuantities &inIq, const IntensiveQuantities &exIq, const Evaluation &dispersivity, const Scalar &normVelocityAvg)
Adds the mass flux due to dispersion to the flux vector over the integration point....
Definition: blackoildispersionmodule.hh:141
static void addDispersiveFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Adds the mass flux due to dispersion to the flux vector over the flux integration point.
Definition: blackoildispersionmodule.hh:105
Provides the auxiliary methods required for consideration of the dispersion equation.
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