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,
class DiscretizationMethod, 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!");
81 {
return referenceSystem; }
95 template<
class Problem,
class ElementVolumeVariables>
96 static NumEqVector
flux(
const Problem& problem,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const SubControlVolumeFace &scvf)
102 NumEqVector flux(0.0);
107 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
110 const int phaseIdx = 0;
112 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
113 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
116 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
117 const Scalar insideDensity =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
119 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
121 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
124 const Scalar massOrMoleFractionInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
125 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
126 const Scalar insideDiffCoeff = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
130 flux[compIdx] = insideDensity * insideDiffCoeff
131 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
135 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
136 const Scalar outsideDiffCoeff = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
137 * outsideVolVars.extrusionFactor();
138 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
139 const Scalar outsideDensity =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
141 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
142 const Scalar avgDiffCoeff =
harmonicMean(insideDiffCoeff, outsideDiffCoeff, insideDistance, outsideDistance);
144 flux[compIdx] = avgDensity * avgDiffCoeff
145 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
150 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
151 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
153 flux *= Extrusion::area(fvGeometry, scvf);
159 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
161 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
Classes related to flux variables caching.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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.
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
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
Definition: method.hh:103
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:44
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:96
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:80
Declares all properties used in Dumux.