25#ifndef DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
29#include <dune/common/fvector.hh>
45template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
46class FicksLawImplementation;
52template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
57 using FVElementGeometry =
typename GridGeometry::LocalView;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
60 using GridView =
typename GridGeometry::GridView;
61 using Element =
typename GridView::template Codim<0>::Entity;
63 using Indices =
typename ModelTraits::Indices;
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 static constexpr int numPhases = ModelTraits::numFluidPhases();
69 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
71 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase supported!");
78 {
return referenceSystem; }
86 template<
class Problem,
class ElementVolumeVariables>
87 static NumEqVector
flux(
const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace &scvf)
93 NumEqVector flux(0.0);
98 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
101 const int phaseIdx = 0;
103 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
104 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
105 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
107 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
108 const Scalar insideDensity =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
110 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
112 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
115 const Scalar massOrMoleFractionInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
116 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
117 const Scalar insideD = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
121 flux[compIdx] = insideDensity * insideD
122 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
126 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
127 const Scalar outsideD = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
128 * outsideVolVars.extrusionFactor();
129 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
130 const Scalar outsideDensity =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
132 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
133 const Scalar avgD =
harmonicMean(insideD, outsideD, insideDistance, outsideDistance);
135 flux[compIdx] = avgDensity * avgD
136 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
141 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
142 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
150 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
152 if constexpr (Dumux::Deprecated::hasEffDiffCoeff<VolumeVariables>)
153 return volVars.effectiveDiffusionCoefficient(phaseIdx,
154 VolumeVariables::FluidSystem::getMainComponent(phaseIdx), compIdx);
158 return volVars.effectiveDiffusivity(phaseIdx, compIdx);
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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...
Classes related to flux variables caching.
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
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
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
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Definition: fluxvariablescaching.hh:64
static NumEqVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: flux/staggered/freeflow/fickslaw.hh:87
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:77
Declares all properties used in Dumux.