28#ifndef OPM_DIFFUSION_MODULE_HH
29#define OPM_DIFFUSION_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/common/Valgrind.hpp>
47template <
class TypeTag,
bool enableDiffusion>
53template <
class TypeTag>
71 template <
class Context>
82template <
class TypeTag>
91 enum { numPhases = FluidSystem::numPhases };
92 enum { numComponents = FluidSystem::numComponents };
93 enum { conti0EqIdx = Indices::conti0EqIdx };
95 using Toolbox = Opm::MathToolbox<Evaluation>;
108 template <
class Context>
110 unsigned spaceIdx,
unsigned timeIdx)
112 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
114 const auto& fluidStateI = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx).fluidState();
115 const auto& fluidStateJ = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx).fluidState();
117 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
119 Evaluation rhoMolar = fluidStateI.molarDensity(phaseIdx);
120 rhoMolar += Toolbox::value(fluidStateJ.molarDensity(phaseIdx));
123 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
125 flux[conti0EqIdx + compIdx] +=
127 * extQuants.moleFractionGradientNormal(phaseIdx, compIdx)
128 * extQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
140template <
class TypeTag,
bool enableDiffusion>
146template <
class TypeTag>
160 throw std::logic_error(
"Method tortuosity() does not make sense "
161 "if diffusion is disabled");
170 throw std::logic_error(
"Method diffusionCoefficient() does not "
171 "make sense if diffusion is disabled");
180 throw std::logic_error(
"Method effectiveDiffusionCoefficient() "
181 "does not make sense if diffusion is disabled");
189 template <
class Flu
idState>
191 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
192 const ElementContext&,
201template <
class TypeTag>
209 enum { numPhases = FluidSystem::numPhases };
210 enum { numComponents = FluidSystem::numComponents };
218 {
return diffusionCoefficient_[phaseIdx][compIdx]; }
225 {
return tortuosity_[phaseIdx]; }
232 {
return tortuosity_[phaseIdx] * diffusionCoefficient_[phaseIdx][compIdx]; }
239 template <
class Flu
idState>
241 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
242 const ElementContext& elemCtx,
246 using Toolbox = Opm::MathToolbox<Evaluation>;
248 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
249 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
250 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
257 const Evaluation& base =
260 * intQuants.fluidState().saturation(phaseIdx));
261 tortuosity_[phaseIdx] =
262 1.0 / (intQuants.porosity() * intQuants.porosity())
263 * Toolbox::pow(base, 10.0/3.0);
265 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
266 diffusionCoefficient_[phaseIdx][compIdx] =
267 FluidSystem::diffusionCoefficient(fluidState,
276 Evaluation tortuosity_[numPhases];
277 Evaluation diffusionCoefficient_[numPhases][numComponents];
286template <
class TypeTag,
bool enableDiffusion>
292template <
class TypeTag>
309 template <
class Context,
class Flu
idState>
326 throw std::logic_error(
"The method moleFractionGradient() does not "
327 "make sense if diffusion is disabled.");
340 throw std::logic_error(
"The method effectiveDiffusionCoefficient() "
341 "does not make sense if diffusion is disabled.");
348template <
class TypeTag>
356 enum { dimWorld = GridView::dimensionworld };
357 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
358 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
360 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
361 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
368 void update_(
const ElementContext& elemCtx,
unsigned faceIdx,
unsigned timeIdx)
370 const auto& gradCalc = elemCtx.gradientCalculator();
373 const auto& face = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
374 const auto& normal = face.normal();
375 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
377 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
378 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
380 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
381 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
385 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
388 DimEvalVector moleFractionGradient(0.0);
389 gradCalc.calculateGradient(moleFractionGradient,
392 moleFractionCallback);
394 moleFractionGradientNormal_[phaseIdx][compIdx] = 0.0;
395 for (
unsigned i = 0; i < normal.size(); ++i)
396 moleFractionGradientNormal_[phaseIdx][compIdx] +=
397 normal[i]*moleFractionGradient[i];
398 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
402 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
403 (intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx)
405 intQuantsOutside.effectiveDiffusionCoefficient(phaseIdx, compIdx))
408 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
413 template <
class Context,
class Flu
idState>
417 const FluidState& fluidState)
419 const auto& stencil = context.stencil(timeIdx);
420 const auto& face = stencil.boundaryFace(bfIdx);
422 const auto& elemCtx = context.elementContext();
423 unsigned insideScvIdx = face.interiorIndex();
424 const auto& insideScv = stencil.subControlVolume(insideScvIdx);
426 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
427 const auto& fluidStateInside = intQuantsInside.fluidState();
430 DimVector distVec = face.integrationPos();
431 distVec -= context.element().geometry().global(insideScv.localGeometry().center());
433 Scalar dist = distVec * face.normal();
439 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
440 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
443 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
446 moleFractionGradientNormal_[phaseIdx][compIdx] =
447 (fluidState.moleFraction(phaseIdx, compIdx)
449 fluidStateInside.moleFraction(phaseIdx, compIdx))
451 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
455 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
456 intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx);
458 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
471 {
return moleFractionGradientNormal_[phaseIdx][compIdx]; }
481 {
return effectiveDiffusionCoefficient_[phaseIdx][compIdx]; }
484 Evaluation moleFractionGradientNormal_[numPhases][numComponents];
485 Evaluation effectiveDiffusionCoefficient_[numPhases][numComponents];
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:304
const Evaluation & effectiveDiffusionCoefficient(unsigned, unsigned) const
The effective diffusion coeffcient of a component in a fluid phase at the face's integration point.
Definition: diffusionmodule.hh:337
const Evaluation & moleFractionGradientNormal(unsigned, unsigned) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:323
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: diffusionmodule.hh:310
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: diffusionmodule.hh:414
const Evaluation & moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:470
const Evaluation & effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
The effective diffusion coeffcient of a component in a fluid phase at the face's integration point.
Definition: diffusionmodule.hh:480
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:368
Provides the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:287
Scalar diffusionCoefficient(unsigned, unsigned) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:168
Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:178
void update_(FluidState &, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &, const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:190
Scalar tortuosity(unsigned) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:158
Evaluation tortuosity(unsigned phaseIdx) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:224
void update_(FluidState &fluidState, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > ¶mCache, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:240
Evaluation effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:231
Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:217
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:141
static void addDiffusiveFlux(RateVector &, const Context &, unsigned, unsigned)
Adds the diffusive mass flux flux to the flux vector over a flux integration point.
Definition: diffusionmodule.hh:72
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:64
static void addDiffusiveFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Adds the mass flux due to molecular diffusion to the flux vector over the flux integration point.
Definition: diffusionmodule.hh:109
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:101
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:48
Callback class for a mole fraction of a component in a phase.
Definition: quantitycallbacks.hh:425
void setComponentIndex(unsigned compIdx)
Set the index of the component for which the mole fraction should be returned.
Definition: quantitycallbacks.hh:460
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which a mole fraction should be returned.
Definition: quantitycallbacks.hh:453
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:37
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:235
This method contains all callback classes for quantities that are required by some extensive quantiti...