version 3.8
flux/staggered/freeflow/fickslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
14#define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
15
16#include <numeric>
17#include <dune/common/fvector.hh>
18
21#include <dumux/common/math.hh>
22
28
29
30namespace Dumux {
31
32// forward declaration
33template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
34class FicksLawImplementation;
35
40template <class TypeTag, ReferenceSystemFormulation referenceSystem>
41class FicksLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
42{
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47 using Extrusion = Extrusion_t<GridGeometry>;
49 using GridView = typename GridGeometry::GridView;
50 using Element = typename GridView::template Codim<0>::Entity;
52 using Indices = typename ModelTraits::Indices;
54
55 static constexpr int numComponents = ModelTraits::numFluidComponents();
56 static constexpr int numPhases = ModelTraits::numFluidPhases();
57
58 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
59
60 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase supported!");
61
62public:
64 // state the discretization method this implementation belongs to
65 static constexpr DiscretizationMethod discMethod{};
66
67 //return the reference system
69 { return referenceSystem; }
70
74
76
83 template<class Problem, class ElementVolumeVariables>
84 static NumEqVector flux(const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace &scvf)
89 {
90 NumEqVector flux(0.0);
91
92 // There is no diffusion over outflow boundaries (grad x == 0).
93 // We assume that if an outflow BC is set for the first transported component, this
94 // also holds for all other components.
95 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
96 return flux;
97
98 const int phaseIdx = 0;
99
100 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
101 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
102 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
103
104 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
105 const Scalar insideDensity = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
106
107 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
108 {
109 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 continue;
111
112 const Scalar massOrMoleFractionInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
113 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
114 const Scalar insideDiffCoeff = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
115
116 if (scvf.boundary())
117 {
118 // When the face lies on a boundary, use the average density
119 const Scalar outsideDensity = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
120 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
121 flux[compIdx] = avgDensity * insideDiffCoeff
122 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
123 }
124 else
125 {
126 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
127 const Scalar outsideDiffCoeff = 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);
131
132 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
133 const Scalar avgDiffCoeff = harmonicMean(insideDiffCoeff, outsideDiffCoeff, insideDistance, outsideDistance);
134
135 flux[compIdx] = avgDensity * avgDiffCoeff
136 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
137 }
138 }
139
140 // Fick's law (for binary systems) states that the net flux of mass within the bulk phase has to be zero:
141 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
142 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
143
144 flux *= Extrusion::area(fvGeometry, scvf);
145
146 return flux;
147 }
148
149private:
150 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
151 {
152 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
153 }
154};
155} // end namespace Dumux
156
157#endif
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
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:84
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:68
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Classes related to flux variables caching.
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:57
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:54
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:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Definition: fluxvariablescaching.hh:55