25#ifndef DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
29#include <dune/common/fvector.hh>
42template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
43class FicksLawImplementation;
49template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
54 using FVElementGeometry =
typename GridGeometry::LocalView;
55 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using GridView =
typename GridGeometry::GridView;
57 using Element =
typename GridView::template Codim<0>::Entity;
63 static constexpr int numComponents = ModelTraits::numFluidComponents();
64 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
66 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase supported!");
73 {
return referenceSystem; }
79 template<
class Problem>
80 static NumEqVector
flux(
const Problem& problem,
81 const Element& element,
82 const FVElementGeometry& fvGeometry,
83 const ElementVolumeVariables& elemVolVars,
84 const SubControlVolumeFace &scvf)
86 NumEqVector flux(0.0);
91 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
94 const int phaseIdx = 0;
96 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
97 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
98 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
100 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
101 const Scalar insideDensity =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
103 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
105 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
108 const Scalar massOrMoleFractionInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
109 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
111 const Scalar insideD = insideVolVars.effectiveDiffusivity(phaseIdx, compIdx) * insideVolVars.extrusionFactor();
115 flux[compIdx] = insideDensity * insideD
116 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
120 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
121 const Scalar outsideD = outsideVolVars.effectiveDiffusivity(phaseIdx, compIdx) * outsideVolVars.extrusionFactor();
122 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
123 const Scalar outsideDensity =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
125 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
126 const Scalar avgD =
harmonicMean(insideD, outsideD, insideDistance, outsideDistance);
128 flux[compIdx] = avgDensity * avgD
129 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
134 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
135 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Classes related to flux variables caching.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
Definition: fluxvariablescaching.hh:64
static NumEqVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: staggered/freeflow/fickslaw.hh:80
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/fickslaw.hh:72
Declares all properties used in Dumux.