25#ifndef DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
42template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
49template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
57 using FVElementGeometry =
typename GridGeometry::LocalView;
58 using SubControlVolume =
typename GridGeometry::SubControlVolume;
59 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
66 using Element =
typename GridView::template Codim<0>::Entity;
69 enum { dim = GridView::dimension} ;
70 enum { dimWorld = GridView::dimensionworld} ;
73 numPhases = ModelTraits::numFluidPhases(),
74 numComponents = ModelTraits::numFluidComponents()
76 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
77 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
84 {
return referenceSystem; }
92 static ComponentFluxVector
flux(
const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace& scvf,
98 const ElementFluxVariablesCache& elemFluxVarsCache)
101 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
102 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
106 const auto& shapeValues = fluxVarsCache.shapeValues();
110 for (
auto&& scv : scvs(fvGeometry))
111 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
113 ComponentFluxVector componentFlux(0.0);
114 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
116 if constexpr (!FluidSystem::isTracerFluidSystem())
117 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
120 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
123 const auto massOrMoleFrac = [&](
const SubControlVolume& scv){
return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
124 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, D, rho);
127 if constexpr (!FluidSystem::isTracerFluidSystem())
128 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
129 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
131 return componentFlux;
135 static std::array<std::vector<Scalar>, numComponents>
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const FluxVarCache& fluxVarCache,
145 const auto& shapeValues = fluxVarCache.shapeValues();
146 for (
auto&& scv : scvs(fvGeometry))
147 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
152 std::array<std::vector<Scalar>, numComponents> ti;
153 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
155 if constexpr (!FluidSystem::isTracerFluidSystem())
156 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
159 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
161 ti[compIdx].resize(fvGeometry.numScv());
162 for (
auto&& scv : scvs(fvGeometry))
163 ti[compIdx][scv.indexInElement()] = -rho*
vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(scvf);
170 static Scalar averageDiffusionCoefficient_(
const int phaseIdx,
const int compIdx,
171 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV,
172 const Problem& problem,
173 const SubControlVolumeFace& scvf)
176 auto [insideD, outsideD] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
179 insideD *= insideVV.extrusionFactor();
180 outsideD *= outsideVV.extrusionFactor();
183 return problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
186 static std::pair<Scalar, Scalar>
187 diffusionCoefficientsAtInterface_([[maybe_unused]]
const int phaseIdx,
const int compIdx,
188 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV)
190 if constexpr (!FluidSystem::isTracerFluidSystem())
192 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
193 const auto insideD = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
194 const auto outsideD = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
195 return { std::move(insideD), std::move(outsideD) };
199 const auto insideD = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
200 const auto outsideD = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
201 return { std::move(insideD), std::move(outsideD) };
205 template<
class EvaluateVariable,
class Tensor>
206 static Scalar discreteFlux_(
const FVElementGeometry& fvGeometry,
207 const SubControlVolumeFace& scvf,
208 const FluxVarCache& fluxVarsCache,
210 const Tensor& D,
const Scalar preFactor)
212 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
213 for (
auto&& scv : scvs(fvGeometry))
215 return -1.0*preFactor*
vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(scvf);
Define some often used mathematical functions.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:66
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:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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:849
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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:149
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:83
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:92
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:136
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Declares all properties used in Dumux.