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>
50template <
class TypeTag,
bool enableDiffusion>
56template <
class TypeTag>
72 template <
class Context>
83template <
class TypeTag>
92 enum { numPhases = FluidSystem::numPhases };
93 enum { numComponents = FluidSystem::numComponents };
94 enum { conti0EqIdx = Indices::conti0EqIdx };
96 using Toolbox = Opm::MathToolbox<Evaluation>;
109 template <
class Context>
111 unsigned spaceIdx,
unsigned timeIdx)
113 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
115 const auto& fluidStateI = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx).fluidState();
116 const auto& fluidStateJ = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx).fluidState();
118 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
120 Evaluation rhoMolar = fluidStateI.molarDensity(phaseIdx);
121 rhoMolar += Toolbox::value(fluidStateJ.molarDensity(phaseIdx));
124 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
126 flux[conti0EqIdx + compIdx] +=
128 * extQuants.moleFractionGradientNormal(phaseIdx, compIdx)
129 * extQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
142template <
class TypeTag,
bool enableDiffusion>
148template <
class TypeTag>
162 throw std::logic_error(
"Method tortuosity() does not make sense "
163 "if diffusion is disabled");
172 throw std::logic_error(
"Method diffusionCoefficient() does not "
173 "make sense if diffusion is disabled");
182 throw std::logic_error(
"Method effectiveDiffusionCoefficient() "
183 "does not make sense if diffusion is disabled");
191 template <
class Flu
idState>
193 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
194 const ElementContext&,
203template <
class TypeTag>
211 enum { numPhases = FluidSystem::numPhases };
212 enum { numComponents = FluidSystem::numComponents };
220 {
return diffusionCoefficient_[phaseIdx][compIdx]; }
227 {
return tortuosity_[phaseIdx]; }
234 {
return tortuosity_[phaseIdx] * diffusionCoefficient_[phaseIdx][compIdx]; }
241 template <
class Flu
idState>
243 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
244 const ElementContext& elemCtx,
248 using Toolbox = Opm::MathToolbox<Evaluation>;
250 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
251 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
252 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
260 const Evaluation& base =
263 * intQuants.fluidState().saturation(phaseIdx));
264 tortuosity_[phaseIdx] =
265 1.0 / (intQuants.porosity() * intQuants.porosity())
266 * Toolbox::pow(base, 10.0/3.0);
268 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
269 diffusionCoefficient_[phaseIdx][compIdx] =
270 FluidSystem::diffusionCoefficient(fluidState,
279 std::array<Evaluation, numPhases> tortuosity_;
280 std::array<std::array<Evaluation, numComponents>, numPhases> diffusionCoefficient_;
289template <
class TypeTag,
bool enableDiffusion>
295template <
class TypeTag>
311 template <
class Context,
class Flu
idState>
328 throw std::logic_error(
"The method moleFractionGradient() does not "
329 "make sense if diffusion is disabled.");
342 throw std::logic_error(
"The method effectiveDiffusionCoefficient() "
343 "does not make sense if diffusion is disabled.");
350template <
class TypeTag>
358 enum { dimWorld = GridView::dimensionworld };
359 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
360 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
362 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
363 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
370 void update_(
const ElementContext& elemCtx,
unsigned faceIdx,
unsigned timeIdx)
372 const auto& gradCalc = elemCtx.gradientCalculator();
375 const auto& face = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
376 const auto& normal = face.normal();
377 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
379 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
380 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
382 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
383 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
388 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
391 DimEvalVector moleFractionGradient(0.0);
392 gradCalc.calculateGradient(moleFractionGradient,
395 moleFractionCallback);
397 moleFractionGradientNormal_[phaseIdx][compIdx] = 0.0;
398 for (
unsigned i = 0; i < normal.size(); ++i) {
399 moleFractionGradientNormal_[phaseIdx][compIdx] +=
400 normal[i]*moleFractionGradient[i];
402 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
406 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
407 (intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx)
409 intQuantsOutside.effectiveDiffusionCoefficient(phaseIdx, compIdx))
412 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
417 template <
class Context,
class Flu
idState>
421 const FluidState& fluidState)
423 const auto& stencil = context.stencil(timeIdx);
424 const auto& face = stencil.boundaryFace(bfIdx);
426 const auto& elemCtx = context.elementContext();
427 const unsigned insideScvIdx = face.interiorIndex();
428 const auto& insideScv = stencil.subControlVolume(insideScvIdx);
430 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
431 const auto& fluidStateInside = intQuantsInside.fluidState();
434 DimVector distVec = face.integrationPos();
435 distVec -= context.element().geometry().global(insideScv.localGeometry().center());
437 const Scalar dist = distVec * face.normal();
443 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
444 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
448 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
451 moleFractionGradientNormal_[phaseIdx][compIdx] =
452 (fluidState.moleFraction(phaseIdx, compIdx)
454 fluidStateInside.moleFraction(phaseIdx, compIdx))
456 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
460 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
461 intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx);
463 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
476 {
return moleFractionGradientNormal_[phaseIdx][compIdx]; }
486 {
return effectiveDiffusionCoefficient_[phaseIdx][compIdx]; }
489 std::array<std::array<Evaluation, numComponents>, numPhases> moleFractionGradientNormal_;
490 std::array<std::array<Evaluation, numComponents>, numPhases> effectiveDiffusionCoefficient_;
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:306
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:339
const Evaluation & moleFractionGradientNormal(unsigned, unsigned) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:325
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: diffusionmodule.hh:312
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: diffusionmodule.hh:418
const Evaluation & moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:475
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:485
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:370
Provides the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:290
Scalar diffusionCoefficient(unsigned, unsigned) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:170
Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:180
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:192
Scalar tortuosity(unsigned) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:160
Evaluation tortuosity(unsigned phaseIdx) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:226
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:242
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:233
Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:219
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:143
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:73
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:65
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:110
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:102
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:51
Callback class for a mole fraction of a component in a phase.
Definition: quantitycallbacks.hh:422
void setComponentIndex(unsigned compIdx)
Set the index of the component for which the mole fraction should be returned.
Definition: quantitycallbacks.hh:457
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which a mole fraction should be returned.
Definition: quantitycallbacks.hh:450
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:43
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
This method contains all callback classes for quantities that are required by some extensive quantiti...