3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
27
28#include <numeric>
29#include <dune/common/fvector.hh>
30
33#include <dumux/common/math.hh>
34
38
39namespace Dumux {
40
41// forward declaration
42template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
43class FicksLawImplementation;
44
49template <class TypeTag, ReferenceSystemFormulation referenceSystem>
50class FicksLawImplementation<TypeTag, DiscretizationMethod::staggered, referenceSystem>
51{
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;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60 using Indices = typename ModelTraits::Indices;
62
63 static constexpr int numComponents = ModelTraits::numFluidComponents();
64 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
65
66 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase supported!");
67
68public:
69 // state the discretization method this implementation belongs to
71 //return the reference system
73 { return referenceSystem; }
74
78
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)
85 {
86 NumEqVector flux(0.0);
87
88 // There is no diffusion over outflow boundaries (grad x == 0).
89 // We assume that if an outflow BC is set for the first transported component, this
90 // also holds for all other components.
91 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
92 return flux;
93
94 const int phaseIdx = 0;
95
96 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
97 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
98 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
99
100 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
101 const Scalar insideDensity = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
102
103 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
104 {
105 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
106 continue;
107
108 const Scalar massOrMoleFractionInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
109 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
110
111 const Scalar insideD = insideVolVars.effectiveDiffusivity(phaseIdx, compIdx) * insideVolVars.extrusionFactor();
112
113 if (scvf.boundary())
114 {
115 flux[compIdx] = insideDensity * insideD
116 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
117 }
118 else
119 {
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);
124
125 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
126 const Scalar avgD = harmonicMean(insideD, outsideD, insideDistance, outsideDistance);
127
128 flux[compIdx] = avgDensity * avgD
129 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
130 }
131 }
132
133 // Fick's law (for binary systems) states that the net flux of mass within the bulk phase has to be zero:
134 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
135 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
136
137 flux *= scvf.area();
138
139 return flux;
140 }
141};
142} // end namespace
143
144#endif
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.