13#ifndef DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
31template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
38template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
46 using FVElementGeometry =
typename GridGeometry::LocalView;
47 using SubControlVolume =
typename GridGeometry::SubControlVolume;
48 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
52 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
53 using FluxVarCache =
typename GridFluxVariablesCache::FluxVariablesCache;
56 using Element =
typename GridView::template Codim<0>::Entity;
59 enum { dim = GridView::dimension} ;
60 enum { dimWorld = GridView::dimensionworld} ;
63 numPhases = ModelTraits::numFluidPhases(),
64 numComponents = ModelTraits::numFluidComponents()
66 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
67 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
74 {
return referenceSystem; }
82 static ComponentFluxVector
flux(
const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
88 const ElementFluxVariablesCache& elemFluxVarsCache)
91 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
92 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
95 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
96 const auto& shapeValues = fluxVarsCache.shapeValues();
100 for (
auto&& scv : scvs(fvGeometry))
101 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
103 ComponentFluxVector componentFlux(0.0);
104 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
106 if constexpr (!FluidSystem::isTracerFluidSystem())
107 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
113 const auto massOrMoleFrac = [&](
const SubControlVolume& scv){
return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
114 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, diffCoeff, rho);
117 if constexpr (!FluidSystem::isTracerFluidSystem())
118 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
119 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
121 return componentFlux;
125 static std::array<std::vector<Scalar>, numComponents>
127 const Element& element,
128 const FVElementGeometry& fvGeometry,
129 const ElementVolumeVariables& elemVolVars,
130 const SubControlVolumeFace& scvf,
131 const FluxVarCache& fluxVarCache,
135 const auto& shapeValues = fluxVarCache.shapeValues();
136 for (
auto&& scv : scvs(fvGeometry))
137 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
139 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
140 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
142 std::array<std::vector<Scalar>, numComponents> ti;
143 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
145 if constexpr (!FluidSystem::isTracerFluidSystem())
146 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
149 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
151 ti[compIdx].resize(fvGeometry.numScv());
152 for (
auto&& scv : scvs(fvGeometry))
153 ti[compIdx][scv.indexInElement()] = -rho*
vtmv(scvf.unitOuterNormal(), diffCoeff, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(fvGeometry, scvf);
160 static Scalar averageDiffusionCoefficient_(
const int phaseIdx,
const int compIdx,
161 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV,
162 const Problem& problem,
163 const SubControlVolumeFace& scvf)
166 auto [insideDiffCoeff, outsideDiffCoeff] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
169 insideDiffCoeff *= insideVV.extrusionFactor();
170 outsideDiffCoeff *= outsideVV.extrusionFactor();
173 return faceTensorAverage(insideDiffCoeff, outsideDiffCoeff, scvf.unitOuterNormal());
176 static std::pair<Scalar, Scalar>
177 diffusionCoefficientsAtInterface_([[maybe_unused]]
const int phaseIdx,
const int compIdx,
178 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV)
180 if constexpr (!FluidSystem::isTracerFluidSystem())
182 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
183 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
184 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
185 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
189 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
190 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
191 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
195 template<
class EvaluateVariable,
class Tensor>
196 static Scalar discreteFlux_(
const FVElementGeometry& fvGeometry,
197 const SubControlVolumeFace& scvf,
198 const FluxVarCache& fluxVarsCache,
200 const Tensor& D,
const Scalar preFactor)
202 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
203 for (
auto&& scv : scvs(fvGeometry))
205 return -1.0*preFactor*
vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(fvGeometry, scvf);
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
static std::array< std::vector< Scalar >, numComponents > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, const int phaseIdx)
Definition: flux/box/fickslaw.hh:126
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:73
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/box/fickslaw.hh:82
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
A free function to average a Tensor at an interface.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:880
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:54
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29