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>
56template <
class TypeTag>
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 };
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;
82 using Toolbox = MathToolbox<Evaluation>;
89 if (!eclState.getSimulationConfig().rock_config().dispersion()) {
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.");
104 template <
class Context>
106 unsigned spaceIdx,
unsigned timeIdx)
109 if (!context.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
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);
140 template<
class IntensiveQuantities,
class Scalar>
142 const IntensiveQuantities& inIq,
143 const IntensiveQuantities& exIq,
144 const Evaluation& dispersivity,
145 const Scalar& normVelocityAvg)
147 const auto& inFs = inIq.fluidState();
148 const auto& exFs = exIq.fluidState();
149 Evaluation diffR = 0.0;
150 if constexpr(enableBioeffects) {
152 const Evaluation bAvg = (inFs.invB(waterPhaseIdx) + Toolbox::value(exFs.invB(waterPhaseIdx))) / 2;
153 diffR = inIq.microbialConcentration() - Toolbox::value(exIq.microbialConcentration());
154 flux[contiMicrobialEqIdx] +=
156 normVelocityAvg[waterPhaseIdx] *
159 if constexpr(enableMICP) {
160 diffR = inIq.oxygenConcentration() - Toolbox::value(exIq.oxygenConcentration());
161 flux[contiOxygenEqIdx] +=
163 normVelocityAvg[waterPhaseIdx] *
166 diffR = inIq.ureaConcentration() - Toolbox::value(exIq.ureaConcentration());
167 flux[contiUreaEqIdx] +=
169 normVelocityAvg[waterPhaseIdx] *
176 unsigned pvtRegionIndex = inFs.pvtRegionIndex();
177 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
178 if (!FluidSystem::phaseIsActive(phaseIdx)) {
183 if (!FluidSystem::enableDissolvedGasInWater() && FluidSystem::waterPhaseIdx == phaseIdx) {
190 if (FluidSystem::gasPhaseIdx == phaseIdx) {
195 Evaluation bAvg = inFs.invB(phaseIdx);
196 bAvg += Toolbox::value(exFs.invB(phaseIdx));
199 Evaluation convFactor = 1.0;
200 if (FluidSystem::enableDissolvedGas() &&
201 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
202 phaseIdx == FluidSystem::oilPhaseIdx)
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());
208 if (FluidSystem::enableVaporizedOil() &&
209 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
210 phaseIdx == FluidSystem::gasPhaseIdx)
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());
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());
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());
230 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
231 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
232 flux[conti0EqIdx + activeSolventCompIdx] +=
234 normVelocityAvg[phaseIdx] *
240 const unsigned soluteCompIdx = FluidSystem::soluteComponentIndex(phaseIdx);
241 const unsigned activeSoluteCompIdx = FluidSystem::canonicalToActiveCompIdx(soluteCompIdx);
242 flux[conti0EqIdx + activeSoluteCompIdx] +=
244 normVelocityAvg[phaseIdx] *
252 static Scalar toMassFractionGasOil (
unsigned regionIdx)
254 const Scalar rhoO = FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, regionIdx);
255 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
259 static Scalar toMassFractionGasWater (
unsigned regionIdx)
261 const Scalar rhoW = FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, regionIdx);
262 const Scalar rhoG = FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regionIdx);
274template <
class TypeTag>
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;
298 {
return normVelocityCell_[phaseIdx]; }
310 template<
class ElementContext>
311 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
314 if (!elemCtx.simulator().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
317 const auto& problem = elemCtx.simulator().problem();
318 if (problem.model().linearizer().getVelocityInfo().empty()) {
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;
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])]));
341 std::array<Scalar, numPhases> normVelocityCell_{};
350template <
class TypeTag>
357 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
362 const IntensiveQuantities& intQuantsInside,
363 const IntensiveQuantities& intQuantsOutside)
365 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
366 if (!FluidSystem::phaseIsActive(phaseIdx)) {
371 normVelocityAvg[phaseIdx] =
372 0.5 * (intQuantsInside.normVelocityCell(phaseIdx) +
373 intQuantsOutside.normVelocityCell(phaseIdx));
374 Valgrind::CheckDefined(normVelocityAvg[phaseIdx]);
379 template <
class Context,
class Flu
idState>
384 {
throw std::runtime_error(
"Not implemented: Dispersion across boundary not implemented for blackoil"); }
394 {
return dispersivity_; }
404 {
return normVelocityAvg_[phaseIdx]; }
407 {
return normVelocityAvg_; }
410 Scalar dispersivity_;
411 ScalarArray normVelocityAvg_;
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