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>
43template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
50template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolume =
typename GridGeometry::SubControlVolume;
60 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
64 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
65 using FluxVarCache =
typename GridFluxVariablesCache::FluxVariablesCache;
68 using Element =
typename GridView::template Codim<0>::Entity;
71 enum { dim = GridView::dimension} ;
72 enum { dimWorld = GridView::dimensionworld} ;
75 numPhases = ModelTraits::numFluidPhases(),
76 numComponents = ModelTraits::numFluidComponents()
78 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
79 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
86 {
return referenceSystem; }
94 static ComponentFluxVector
flux(
const Problem& problem,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const SubControlVolumeFace& scvf,
100 const ElementFluxVariablesCache& elemFluxVarsCache)
103 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
107 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
108 const auto& shapeValues = fluxVarsCache.shapeValues();
112 for (
auto&& scv : scvs(fvGeometry))
113 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
115 ComponentFluxVector componentFlux(0.0);
116 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
118 if constexpr (!FluidSystem::isTracerFluidSystem())
119 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
122 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
125 const auto massOrMoleFrac = [&](
const SubControlVolume& scv){
return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
126 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, diffCoeff, rho);
129 if constexpr (!FluidSystem::isTracerFluidSystem())
130 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
131 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
133 return componentFlux;
137 static std::array<std::vector<Scalar>, numComponents>
139 const Element& element,
140 const FVElementGeometry& fvGeometry,
141 const ElementVolumeVariables& elemVolVars,
142 const SubControlVolumeFace& scvf,
143 const FluxVarCache& fluxVarCache,
147 const auto& shapeValues = fluxVarCache.shapeValues();
148 for (
auto&& scv : scvs(fvGeometry))
149 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
151 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
152 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
154 std::array<std::vector<Scalar>, numComponents> ti;
155 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
157 if constexpr (!FluidSystem::isTracerFluidSystem())
158 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
161 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
163 ti[compIdx].resize(fvGeometry.numScv());
164 for (
auto&& scv : scvs(fvGeometry))
165 ti[compIdx][scv.indexInElement()] = -rho*
vtmv(scvf.unitOuterNormal(), diffCoeff, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(fvGeometry, scvf);
172 static Scalar averageDiffusionCoefficient_(
const int phaseIdx,
const int compIdx,
173 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV,
174 const Problem& problem,
175 const SubControlVolumeFace& scvf)
178 auto [insideDiffCoeff, outsideDiffCoeff] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
181 insideDiffCoeff *= insideVV.extrusionFactor();
182 outsideDiffCoeff *= outsideVV.extrusionFactor();
185 return faceTensorAverage(insideDiffCoeff, outsideDiffCoeff, scvf.unitOuterNormal());
188 static std::pair<Scalar, Scalar>
189 diffusionCoefficientsAtInterface_([[maybe_unused]]
const int phaseIdx,
const int compIdx,
190 const VolumeVariables& insideVV,
const VolumeVariables& outsideVV)
192 if constexpr (!FluidSystem::isTracerFluidSystem())
194 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
195 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
196 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
197 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
201 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
202 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
203 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
207 template<
class EvaluateVariable,
class Tensor>
208 static Scalar discreteFlux_(
const FVElementGeometry& fvGeometry,
209 const SubControlVolumeFace& scvf,
210 const FluxVarCache& fluxVarsCache,
212 const Tensor& D,
const Scalar preFactor)
214 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
215 for (
auto&& scv : scvs(fvGeometry))
217 return -1.0*preFactor*
vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(fvGeometry, scvf);
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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...
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:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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:41
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:44
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:138
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:85
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:94
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.