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>
41#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
42#include <opm/common/OpmLog/OpmLog.hpp>
57template <
class TypeTag,
bool enableDispersion>
60template <
class TypeTag,
bool enableDispersion>
66template <
class TypeTag>
75 enum { numPhases = FluidSystem::numPhases };
81 static void initFromState(
const EclipseState&)
89 template <
class Context>
96 template<
class IntensiveQuantities,
class Scalar>
98 const IntensiveQuantities&,
99 const IntensiveQuantities&,
108template <
class TypeTag>
123 enum { numPhases = FluidSystem::numPhases };
124 enum { numComponents = FluidSystem::numComponents };
125 enum { conti0EqIdx = Indices::conti0EqIdx };
126 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
127 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
129 static constexpr unsigned contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
130 static constexpr unsigned contiOxygenEqIdx = Indices::contiOxygenEqIdx;
131 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
132 static constexpr unsigned contiUreaEqIdx = Indices::contiUreaEqIdx;
134 using Toolbox = MathToolbox<Evaluation>;
140 static void initFromState(
const EclipseState& eclState)
142 if (!eclState.getSimulationConfig().rock_config().dispersion()) {
146 if (eclState.getSimulationConfig().hasVAPWAT() || eclState.getSimulationConfig().hasVAPOIL()) {
147 OpmLog::warning(
"Dispersion is activated in combination with VAPWAT/VAPOIL. \n"
148 "Water/oil is still allowed to vaporize, but dispersion in the "
149 "gas phase is ignored.");
158 template <
class Context>
160 unsigned spaceIdx,
unsigned timeIdx)
163 if (!context.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
166 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
167 const auto& inIq = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
168 const auto& exIq = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
169 const auto& dispersivity = extQuants.dispersivity();
170 const auto& normVelocityAvg = extQuants.normVelocityAvg();
171 addDispersiveFlux(flux, inIq, exIq, dispersivity, normVelocityAvg);
194 template<
class IntensiveQuantities,
class Scalar>
196 const IntensiveQuantities& inIq,
197 const IntensiveQuantities& exIq,
198 const Evaluation& dispersivity,
199 const Scalar& normVelocityAvg)
201 const auto& inFs = inIq.fluidState();
202 const auto& exFs = exIq.fluidState();
203 Evaluation diffR = 0.0;
204 if constexpr(enableMICP) {
206 const Evaluation bAvg = (inFs.invB(waterPhaseIdx) + Toolbox::value(exFs.invB(waterPhaseIdx))) / 2;
207 diffR = inIq.microbialConcentration() - Toolbox::value(exIq.microbialConcentration());
208 flux[contiMicrobialEqIdx] +=
210 normVelocityAvg[waterPhaseIdx] *
213 diffR = inIq.oxygenConcentration() - Toolbox::value(exIq.oxygenConcentration());
214 flux[contiOxygenEqIdx] +=
216 normVelocityAvg[waterPhaseIdx] *
219 diffR = inIq.ureaConcentration() - Toolbox::value(exIq.ureaConcentration());
220 flux[contiUreaEqIdx] +=
222 normVelocityAvg[waterPhaseIdx] *
228 unsigned pvtRegionIndex = inFs.pvtRegionIndex();
229 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
230 if (!FluidSystem::phaseIsActive(phaseIdx)) {
235 if (!FluidSystem::enableDissolvedGasInWater() && FluidSystem::waterPhaseIdx == phaseIdx) {
242 if (FluidSystem::gasPhaseIdx == phaseIdx) {
247 Evaluation bAvg = inFs.invB(phaseIdx);
248 bAvg += Toolbox::value(exFs.invB(phaseIdx));
251 Evaluation convFactor = 1.0;
252 if (FluidSystem::enableDissolvedGas() &&
253 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
254 phaseIdx == FluidSystem::oilPhaseIdx)
256 const Evaluation rsAvg = (inFs.Rs() + Toolbox::value(exFs.Rs())) / 2;
257 convFactor = 1.0 / (toMassFractionGasOil(pvtRegionIndex) + rsAvg);
258 diffR = inFs.Rs() - Toolbox::value(exFs.Rs());
260 if (FluidSystem::enableVaporizedOil() &&
261 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
262 phaseIdx == FluidSystem::gasPhaseIdx)
264 const Evaluation rvAvg = (inFs.Rv() + Toolbox::value(exFs.Rv())) / 2;
265 convFactor = toMassFractionGasOil(pvtRegionIndex) /
266 (1.0 + rvAvg * toMassFractionGasOil(pvtRegionIndex));
267 diffR = inFs.Rv() - Toolbox::value(exFs.Rv());
269 if (FluidSystem::enableDissolvedGasInWater() && phaseIdx == FluidSystem::waterPhaseIdx) {
270 const Evaluation rsAvg = (inFs.Rsw() + Toolbox::value(exFs.Rsw())) / 2;
271 convFactor = 1.0 / (toMassFractionGasWater(pvtRegionIndex) + rsAvg);
272 diffR = inFs.Rsw() - Toolbox::value(exFs.Rsw());
274 if (FluidSystem::enableVaporizedWater() && phaseIdx == FluidSystem::gasPhaseIdx) {
275 const Evaluation rvAvg = (inFs.Rvw() + Toolbox::value(exFs.Rvw())) / 2;
276 convFactor = toMassFractionGasWater(pvtRegionIndex) /
277 (1.0 + rvAvg * toMassFractionGasWater(pvtRegionIndex));
278 diffR = inFs.Rvw() - Toolbox::value(exFs.Rvw());
282 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
283 const unsigned activeSolventCompIdx = Indices::canonicalToActiveComponentIndex(solventCompIdx);
284 flux[conti0EqIdx + activeSolventCompIdx] +=
286 normVelocityAvg[phaseIdx] *
292 const unsigned soluteCompIdx = FluidSystem::soluteComponentIndex(phaseIdx);
293 const unsigned activeSoluteCompIdx = Indices::canonicalToActiveComponentIndex(soluteCompIdx);
294 flux[conti0EqIdx + activeSoluteCompIdx] +=
296 normVelocityAvg[phaseIdx] *
304 static Scalar toMassFractionGasOil (
unsigned regionIdx)
306 const Scalar rhoO = FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, regionIdx);
307 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
311 static Scalar toMassFractionGasWater (
unsigned regionIdx)
313 const Scalar rhoW = FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, regionIdx);
314 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
326template <
class TypeTag,
bool enableDispersion>
332template <
class TypeTag>
345 throw std::logic_error(
"Method normVelocityCell() "
346 "does not make sense if dispersion is disabled");
354 template<
class ElementContext>
364template <
class TypeTag>
373 enum { numPhases = FluidSystem::numPhases };
374 enum { numComponents = FluidSystem::numComponents };
375 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
376 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
377 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
378 enum { gasCompIdx = FluidSystem::gasCompIdx };
379 enum { oilCompIdx = FluidSystem::oilCompIdx };
380 enum { waterCompIdx = FluidSystem::waterCompIdx };
381 enum { conti0EqIdx = Indices::conti0EqIdx };
382 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
389 {
return normVelocityCell_[phaseIdx]; }
401 template<
class ElementContext>
402 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
405 if (!elemCtx.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
408 const auto& problem = elemCtx.simulator().problem();
409 if (problem.model().linearizer().getVelocityInfo().empty()) {
412 const std::array<int, 3> phaseIdxs = {gasPhaseIdx, oilPhaseIdx, waterPhaseIdx};
413 const std::array<int, 3> compIdxs = {gasCompIdx, oilCompIdx, waterCompIdx};
414 const auto& velocityInf = problem.model().linearizer().getVelocityInfo();
415 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
416 const auto& velocityInfos = velocityInf[globalDofIdx];
417 for (
unsigned i = 0; i < phaseIdxs.size(); ++i) {
418 normVelocityCell_[i] = 0;
420 for (
const auto& velocityInfo : velocityInfos) {
421 for (
unsigned i = 0; i < phaseIdxs.size(); ++i) {
422 if (FluidSystem::phaseIsActive(phaseIdxs[i])) {
423 normVelocityCell_[phaseIdxs[i]] = max( normVelocityCell_[phaseIdxs[i]],
424 std::abs(velocityInfo.velocity[conti0EqIdx +
425 Indices::canonicalToActiveComponentIndex(compIdxs[i])]));
432 std::array<Scalar, numPhases> normVelocityCell_{};
441template <
class TypeTag,
bool enableDispersion>
442class BlackOilDispersionExtensiveQuantities;
447template <
class TypeTag>
456 enum { numPhases = FluidSystem::numPhases };
468 template <
class Context,
class Flu
idState>
479 const IntensiveQuantities&,
480 const IntensiveQuantities&)
489 throw std::logic_error(
"The method dispersivity() does not "
490 "make sense if dispersion is disabled.");
502 throw std::logic_error(
"The method normVelocityAvg() "
503 "does not make sense if dispersion is disabled.");
510template <
class TypeTag>
518 using Toolbox = MathToolbox<Evaluation>;
521 enum { dimWorld = GridView::dimensionworld };
522 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
523 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
525 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
526 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
531 const IntensiveQuantities& intQuantsInside,
532 const IntensiveQuantities& intQuantsOutside)
534 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
535 if (!FluidSystem::phaseIsActive(phaseIdx)) {
539 if (!FluidSystem::enableDissolvedGasInWater() && FluidSystem::waterPhaseIdx == phaseIdx) {
545 if (FluidSystem::gasPhaseIdx == phaseIdx) {
550 normVelocityAvg[phaseIdx] =
551 0.5 * (intQuantsInside.normVelocityCell(phaseIdx) +
552 intQuantsOutside.normVelocityCell(phaseIdx));
553 Valgrind::CheckDefined(normVelocityAvg[phaseIdx]);
558 template <
class Context,
class Flu
idState>
563 {
throw std::runtime_error(
"Not implemented: Dispersion across boundary not implemented for blackoil"); }
573 {
return dispersivity_; }
583 {
return normVelocityAvg_[phaseIdx]; }
586 {
return normVelocityAvg_; }
589 Scalar dispersivity_;
590 ScalarArray normVelocityAvg_;
Declares the properties required by the black oil model.
Definition: blackoildispersionmodule.hh:449
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: blackoildispersionmodule.hh:469
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the dispersive fluxes.
Definition: blackoildispersionmodule.hh:463
Scalar[numPhases] ScalarArray
Definition: blackoildispersionmodule.hh:476
static void update(ScalarArray &, const IntensiveQuantities &, const IntensiveQuantities &)
Definition: blackoildispersionmodule.hh:478
Scalar normVelocityAvg(unsigned) const
The effective filter velocity coefficient in a fluid phase at the face's integration point.
Definition: blackoildispersionmodule.hh:500
Scalar dispersivity() const
The dispersivity the face.
Definition: blackoildispersionmodule.hh:487
Provides the quantities required to calculate dispersive mass fluxes.
Definition: blackoildispersionmodule.hh:512
Scalar dispersivity() const
The dispersivity of the face.
Definition: blackoildispersionmodule.hh:572
static void update(ScalarArray &normVelocityAvg, const IntensiveQuantities &intQuantsInside, const IntensiveQuantities &intQuantsOutside)
Definition: blackoildispersionmodule.hh:530
std::array< Scalar, numPhases > ScalarArray
Definition: blackoildispersionmodule.hh:529
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: blackoildispersionmodule.hh:559
const auto & normVelocityAvg() const
Definition: blackoildispersionmodule.hh:585
Scalar normVelocityAvg(unsigned phaseIdx) const
The effective velocity coefficient in a fluid phase at the face's integration point.
Definition: blackoildispersionmodule.hh:582
Provides the quantities required to calculate dispersive mass fluxes.
Definition: blackoildispersionmodule.hh:61
Scalar normVelocityCell(unsigned, unsigned) const
Returns the max. norm of the filter velocity of the cell.
Definition: blackoildispersionmodule.hh:343
void update_(ElementContext &, unsigned, unsigned)
Update the quantities required to calculate dispersive fluxes.
Definition: blackoildispersionmodule.hh:355
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:402
Scalar normVelocityCell(unsigned phaseIdx) const
Returns the max. norm of the filter velocity of the cell.
Definition: blackoildispersionmodule.hh:388
Provides the volumetric quantities required for the calculation of dispersive fluxes.
Definition: blackoildispersionmodule.hh:327
static void addDispersiveFlux(RateVector &, const Context &, unsigned, unsigned)
Adds the dispersive flux to the flux vector over a flux integration point.
Definition: blackoildispersionmodule.hh:90
static void addDispersiveFlux(RateVector &, const IntensiveQuantities &, const IntensiveQuantities &, const Evaluation &, const Scalar &)
Definition: blackoildispersionmodule.hh:97
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:195
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:159
Provides the auxiliary methods required for consideration of the dispersion equation.
Definition: blackoildispersionmodule.hh:58
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: blackoilboundaryratevector.hh:39
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