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 GridGeometry::SubControlVolumeFace;
61 using GridView =
typename GridGeometry::GridView;
62 using Element =
typename GridView::template Codim<0>::Entity;
64 using Indices =
typename ModelTraits::Indices;
67 static constexpr int numComponents = ModelTraits::numFluidComponents();
68 static constexpr int numPhases = ModelTraits::numFluidPhases();
70 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
72 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase supported!");
79 {
return referenceSystem; }
93 template<
class Problem,
class ElementVolumeVariables>
94 static NumEqVector
flux(
const Problem& problem,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const SubControlVolumeFace &scvf)
100 NumEqVector flux(0.0);
105 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
108 const int phaseIdx = 0;
110 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
111 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
112 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
114 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
115 const Scalar insideDensity =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
117 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
119 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
122 const Scalar massOrMoleFractionInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
123 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
124 const Scalar insideD = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
128 flux[compIdx] = insideDensity * insideD
129 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
133 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
134 const Scalar outsideD = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
135 * outsideVolVars.extrusionFactor();
136 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
137 const Scalar outsideDensity =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
139 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
140 const Scalar avgD =
harmonicMean(insideD, outsideD, insideDistance, outsideDistance);
142 flux[compIdx] = avgDensity * avgD
143 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
148 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
149 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
151 flux *= Extrusion::area(scvf);
157 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
159 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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...
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:69
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
Definition: propertysystem.hh:150
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:67
static NumEqVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/staggered/freeflow/fickslaw.hh:94
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:78
Declares all properties used in Dumux.