28#ifndef EWOMS_DISPERSION_MODULE_HH
29#define EWOMS_DISPERSION_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/common/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
40#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
41#include <opm/common/OpmLog/OpmLog.hpp>
55template <
class TypeTag,
bool enableDispersion>
58template <
class TypeTag,
bool enableDispersion>
64template <
class TypeTag>
83 template <
class Context>
90 template<
class IntensiveQuantities,
class Scalar>
92 const IntensiveQuantities&,
93 const IntensiveQuantities&,
102template <
class TypeTag>
116 enum { numPhases = FluidSystem::numPhases };
117 enum { numComponents = FluidSystem::numComponents };
118 enum { conti0EqIdx = Indices::conti0EqIdx };
119 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
120 enum { enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>() };
121 enum { enableMICP = Indices::enableMICP };
123 static constexpr unsigned contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
124 static constexpr unsigned contiOxygenEqIdx = Indices::contiOxygenEqIdx;
125 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
126 static constexpr unsigned contiUreaEqIdx = Indices::contiUreaEqIdx;
128 using Toolbox = MathToolbox<Evaluation>;
135 if (!eclState.getSimulationConfig().rock_config().dispersion()) {
139 if (eclState.getSimulationConfig().hasVAPWAT() || eclState.getSimulationConfig().hasVAPOIL()) {
140 OpmLog::warning(
"Dispersion is activated in combination with VAPWAT/VAPOIL. \n"
141 "Water/oil is still allowed to vaporize, but dispersion in the "
142 "gas phase is ignored.");
150 template <
class Context>
152 unsigned spaceIdx,
unsigned timeIdx)
155 if (!context.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
158 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
159 const auto& inIq = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
160 const auto& exIq = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
161 const auto& dispersivity = extQuants.dispersivity();
162 const auto& normVelocityAvg = extQuants.normVelocityAvg();
163 addDispersiveFlux(flux, inIq, exIq, dispersivity, normVelocityAvg);
186 template<
class IntensiveQuantities,
class Scalar>
188 const IntensiveQuantities& inIq,
189 const IntensiveQuantities& exIq,
190 const Evaluation& dispersivity,
191 const Scalar& normVelocityAvg)
193 const auto& inFs = inIq.fluidState();
194 const auto& exFs = exIq.fluidState();
195 Evaluation diffR = 0.0;
196 if constexpr(enableBioeffects) {
198 const Evaluation bAvg = (inFs.invB(waterPhaseIdx) + Toolbox::value(exFs.invB(waterPhaseIdx))) / 2;
199 diffR = inIq.microbialConcentration() - Toolbox::value(exIq.microbialConcentration());
200 flux[contiMicrobialEqIdx] +=
202 normVelocityAvg[waterPhaseIdx] *
205 if constexpr(enableMICP) {
206 diffR = inIq.oxygenConcentration() - Toolbox::value(exIq.oxygenConcentration());
207 flux[contiOxygenEqIdx] +=
209 normVelocityAvg[waterPhaseIdx] *
212 diffR = inIq.ureaConcentration() - Toolbox::value(exIq.ureaConcentration());
213 flux[contiUreaEqIdx] +=
215 normVelocityAvg[waterPhaseIdx] *
222 unsigned pvtRegionIndex = inFs.pvtRegionIndex();
223 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
224 if (!FluidSystem::phaseIsActive(phaseIdx)) {
229 if (!FluidSystem::enableDissolvedGasInWater() && FluidSystem::waterPhaseIdx == phaseIdx) {
236 if (FluidSystem::gasPhaseIdx == phaseIdx) {
241 Evaluation bAvg = inFs.invB(phaseIdx);
242 bAvg += Toolbox::value(exFs.invB(phaseIdx));
245 Evaluation convFactor = 1.0;
246 if (FluidSystem::enableDissolvedGas() &&
247 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
248 phaseIdx == FluidSystem::oilPhaseIdx)
250 const Evaluation rsAvg = (inFs.Rs() + Toolbox::value(exFs.Rs())) / 2;
251 convFactor = 1.0 / (toMassFractionGasOil(pvtRegionIndex) + rsAvg);
252 diffR = inFs.Rs() - Toolbox::value(exFs.Rs());
254 if (FluidSystem::enableVaporizedOil() &&
255 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
256 phaseIdx == FluidSystem::gasPhaseIdx)
258 const Evaluation rvAvg = (inFs.Rv() + Toolbox::value(exFs.Rv())) / 2;
259 convFactor = toMassFractionGasOil(pvtRegionIndex) /
260 (1.0 + rvAvg * toMassFractionGasOil(pvtRegionIndex));
261 diffR = inFs.Rv() - Toolbox::value(exFs.Rv());
263 if (FluidSystem::enableDissolvedGasInWater() && phaseIdx == FluidSystem::waterPhaseIdx) {
264 const Evaluation rsAvg = (inFs.Rsw() + Toolbox::value(exFs.Rsw())) / 2;
265 convFactor = 1.0 / (toMassFractionGasWater(pvtRegionIndex) + rsAvg);
266 diffR = inFs.Rsw() - Toolbox::value(exFs.Rsw());
268 if (FluidSystem::enableVaporizedWater() && phaseIdx == FluidSystem::gasPhaseIdx) {
269 const Evaluation rvAvg = (inFs.Rvw() + Toolbox::value(exFs.Rvw())) / 2;
270 convFactor = toMassFractionGasWater(pvtRegionIndex) /
271 (1.0 + rvAvg * toMassFractionGasWater(pvtRegionIndex));
272 diffR = inFs.Rvw() - Toolbox::value(exFs.Rvw());
276 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
277 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
278 flux[conti0EqIdx + activeSolventCompIdx] +=
280 normVelocityAvg[phaseIdx] *
286 const unsigned soluteCompIdx = FluidSystem::soluteComponentIndex(phaseIdx);
287 const unsigned activeSoluteCompIdx = FluidSystem::canonicalToActiveCompIdx(soluteCompIdx);
288 flux[conti0EqIdx + activeSoluteCompIdx] +=
290 normVelocityAvg[phaseIdx] *
298 static Scalar toMassFractionGasOil (
unsigned regionIdx)
300 const Scalar rhoO = FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, regionIdx);
301 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
305 static Scalar toMassFractionGasWater (
unsigned regionIdx)
307 const Scalar rhoW = FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, regionIdx);
308 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
320template <
class TypeTag,
bool enableDispersion>
326template <
class TypeTag>
338 throw std::logic_error(
"Method normVelocityCell() "
339 "does not make sense if dispersion is disabled");
347 template<
class ElementContext>
357template <
class TypeTag>
365 enum { numPhases = FluidSystem::numPhases };
366 enum { numComponents = FluidSystem::numComponents };
367 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
368 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
369 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
370 enum { gasCompIdx = FluidSystem::gasCompIdx };
371 enum { oilCompIdx = FluidSystem::oilCompIdx };
372 enum { waterCompIdx = FluidSystem::waterCompIdx };
373 enum { conti0EqIdx = Indices::conti0EqIdx };
374 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
381 {
return normVelocityCell_[phaseIdx]; }
393 template<
class ElementContext>
394 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
397 if (!elemCtx.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
400 const auto& problem = elemCtx.simulator().problem();
401 if (problem.model().linearizer().getVelocityInfo().empty()) {
404 const std::array<int, 3> phaseIdxs = {gasPhaseIdx, oilPhaseIdx, waterPhaseIdx};
405 const std::array<int, 3> compIdxs = {gasCompIdx, oilCompIdx, waterCompIdx};
406 const auto& velocityInf = problem.model().linearizer().getVelocityInfo();
407 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
408 const auto& velocityInfos = velocityInf[globalDofIdx];
409 for (
unsigned i = 0; i < phaseIdxs.size(); ++i) {
410 normVelocityCell_[i] = 0;
412 for (
const auto& velocityInfo : velocityInfos) {
413 for (
unsigned i = 0; i < phaseIdxs.size(); ++i) {
414 if (FluidSystem::phaseIsActive(phaseIdxs[i])) {
415 normVelocityCell_[phaseIdxs[i]] = std::max(normVelocityCell_[phaseIdxs[i]],
416 std::abs(velocityInfo.velocity[conti0EqIdx +
417 FluidSystem::canonicalToActiveCompIdx(compIdxs[i])]));
424 std::array<Scalar, numPhases> normVelocityCell_{};
433template <
class TypeTag,
bool enableDispersion>
434class BlackOilDispersionExtensiveQuantities;
439template <
class TypeTag>
447 enum { numPhases = FluidSystem::numPhases };
459 template <
class Context,
class Flu
idState>
470 const IntensiveQuantities&,
471 const IntensiveQuantities&)
480 throw std::logic_error(
"The method dispersivity() does not "
481 "make sense if dispersion is disabled.");
493 throw std::logic_error(
"The method normVelocityAvg() "
494 "does not make sense if dispersion is disabled.");
501template <
class TypeTag>
508 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
513 const IntensiveQuantities& intQuantsInside,
514 const IntensiveQuantities& intQuantsOutside)
516 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
517 if (!FluidSystem::phaseIsActive(phaseIdx)) {
522 normVelocityAvg[phaseIdx] =
523 0.5 * (intQuantsInside.normVelocityCell(phaseIdx) +
524 intQuantsOutside.normVelocityCell(phaseIdx));
525 Valgrind::CheckDefined(normVelocityAvg[phaseIdx]);
530 template <
class Context,
class Flu
idState>
535 {
throw std::runtime_error(
"Not implemented: Dispersion across boundary not implemented for blackoil"); }
545 {
return dispersivity_; }
555 {
return normVelocityAvg_[phaseIdx]; }
558 {
return normVelocityAvg_; }
561 Scalar dispersivity_;
562 ScalarArray normVelocityAvg_;
Declares the properties required by the black oil model.
Definition: blackoildispersionmodule.hh:441
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: blackoildispersionmodule.hh:460
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the dispersive fluxes.
Definition: blackoildispersionmodule.hh:454
Scalar[numPhases] ScalarArray
Definition: blackoildispersionmodule.hh:467
static void update(ScalarArray &, const IntensiveQuantities &, const IntensiveQuantities &)
Definition: blackoildispersionmodule.hh:469
Scalar normVelocityAvg(unsigned) const
The effective filter velocity coefficient in a fluid phase at the face's integration point.
Definition: blackoildispersionmodule.hh:491
Scalar dispersivity() const
The dispersivity the face.
Definition: blackoildispersionmodule.hh:478
Provides the quantities required to calculate dispersive mass fluxes.
Definition: blackoildispersionmodule.hh:503
Scalar dispersivity() const
The dispersivity of the face.
Definition: blackoildispersionmodule.hh:544
static void update(ScalarArray &normVelocityAvg, const IntensiveQuantities &intQuantsInside, const IntensiveQuantities &intQuantsOutside)
Definition: blackoildispersionmodule.hh:512
std::array< Scalar, numPhases > ScalarArray
Definition: blackoildispersionmodule.hh:511
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: blackoildispersionmodule.hh:531
const auto & normVelocityAvg() const
Definition: blackoildispersionmodule.hh:557
Scalar normVelocityAvg(unsigned phaseIdx) const
The effective velocity coefficient in a fluid phase at the face's integration point.
Definition: blackoildispersionmodule.hh:554
Provides the quantities required to calculate dispersive mass fluxes.
Definition: blackoildispersionmodule.hh:59
Scalar normVelocityCell(unsigned, unsigned) const
Returns the max. norm of the filter velocity of the cell.
Definition: blackoildispersionmodule.hh:336
void update_(ElementContext &, unsigned, unsigned)
Update the quantities required to calculate dispersive fluxes.
Definition: blackoildispersionmodule.hh:348
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:394
Scalar normVelocityCell(unsigned phaseIdx) const
Returns the max. norm of the filter velocity of the cell.
Definition: blackoildispersionmodule.hh:380
Provides the volumetric quantities required for the calculation of dispersive fluxes.
Definition: blackoildispersionmodule.hh:321
static void addDispersiveFlux(RateVector &, const Context &, unsigned, unsigned)
Adds the dispersive flux to the flux vector over a flux integration point.
Definition: blackoildispersionmodule.hh:84
static void addDispersiveFlux(RateVector &, const IntensiveQuantities &, const IntensiveQuantities &, const Evaluation &, const Scalar &)
Definition: blackoildispersionmodule.hh:91
static void initFromState(const EclipseState &)
Definition: blackoildispersionmodule.hh:76
static void initFromState(const EclipseState &eclState)
Definition: blackoildispersionmodule.hh:133
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:187
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:151
Provides the auxiliary methods required for consideration of the dispersion equation.
Definition: blackoildispersionmodule.hh:56
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