3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/cctpfa/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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
26
27#include <dune/common/fvector.hh>
28
31
36
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::cctpfa, referenceSystem>
51{
56 using FVElementGeometry = typename GridGeometry::LocalView;
57 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
58 using Extrusion = Extrusion_t<GridGeometry>;
59 using GridView = typename GridGeometry::GridView;
60 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61 using Element = typename GridView::template Codim<0>::Entity;
62 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
66
68 static const int dim = GridView::dimension;
69 static const int dimWorld = GridView::dimensionworld;
70 static const int numPhases = ModelTraits::numFluidPhases();
71 static const int numComponents = ModelTraits::numFluidComponents();
72
73 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
74
76 class TpfaFicksLawCacheFiller
77 {
78 public:
81 template<class FluxVariablesCacheFiller>
82 static void fill(FluxVariablesCache& scvfFluxVarsCache,
83 unsigned int phaseIdx, unsigned int compIdx,
84 const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace& scvf,
89 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
90 {
91 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
92 }
93 };
94
96 class TpfaFicksLawCache
97 {
98 public:
99 using Filler = TpfaFicksLawCacheFiller;
100
101 void updateDiffusion(const Problem& problem,
102 const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const ElementVolumeVariables& elemVolVars,
105 const SubControlVolumeFace &scvf,
106 const unsigned int phaseIdx,
107 const unsigned int compIdx)
108 {
109 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
110 }
111
112 const Scalar& diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
113 { return tij_[phaseIdx][compIdx]; }
114
115 private:
116 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
117 };
118
119public:
124 { return referenceSystem; }
125
127 using Cache = TpfaFicksLawCache;
128
130
137 static ComponentFluxVector flux(const Problem& problem,
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf,
142 const int phaseIdx,
143 const ElementFluxVariablesCache& elemFluxVarsCache)
144 {
145 ComponentFluxVector componentFlux(0.0);
146 for (int compIdx = 0; compIdx < numComponents; compIdx++)
147 {
148 if constexpr (!FluidSystem::isTracerFluidSystem())
149 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
150 continue;
151
152 // diffusion tensors are always solution dependent
153 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
154
155 // get inside/outside volume variables
156 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
157 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
158
159 // the inside and outside mass/mole fractions fractions
160 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
161 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
162 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
163 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
164 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
165
166 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
167 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
168
169 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
170 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
171
172 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
173 if constexpr (!FluidSystem::isTracerFluidSystem())
174 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
175 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
176 }
177
178 return componentFlux;
179 }
180
182 static Scalar calculateTransmissibility(const Problem& problem,
183 const Element& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars,
186 const SubControlVolumeFace& scvf,
187 const int phaseIdx, const int compIdx)
188 {
189
190
191 const auto insideScvIdx = scvf.insideScvIdx();
192 const auto& insideScv = fvGeometry.scv(insideScvIdx);
193 const auto& insideVolVars = elemVolVars[insideScvIdx];
194 const auto getDiffCoeff = [&](const auto& vv)
195 {
196 if constexpr (FluidSystem::isTracerFluidSystem())
197 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
198 else
199 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
200 };
201
202 const auto insideD = getDiffCoeff(insideVolVars);
203
204 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideD, insideVolVars.extrusionFactor());
205
206 // for the boundary (dirichlet) or at branching points we only need ti
207 Scalar tij;
208 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
209 tij = Extrusion::area(scvf)*ti;
210
211 // otherwise we compute a tpfa harmonic mean
212 else
213 {
214 const auto outsideScvIdx = scvf.outsideScvIdx();
215 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
216 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
217 const auto outsideD = getDiffCoeff(outsideVolVars);
218
219 Scalar tj;
220 if (dim == dimWorld)
221 // assume the normal vector from outside is anti parallel so we save flipping a vector
222 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideD, outsideVolVars.extrusionFactor());
223 else
224 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()),
225 outsideScv,
226 outsideD,
227 outsideVolVars.extrusionFactor());
228
229 // check if we are dividing by zero!
230 if (ti*tj <= 0.0)
231 tij = 0;
232 else
233 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
234 }
235
236 return tij;
237 }
238
239private:
241 static Scalar branchingFacetX(const Problem& problem,
242 const Element& element,
243 const FVElementGeometry& fvGeometry,
244 const ElementVolumeVariables& elemVolVars,
245 const ElementFluxVariablesCache& elemFluxVarsCache,
246 const SubControlVolumeFace& scvf,
247 const Scalar insideX, const Scalar insideTi,
248 const int phaseIdx, const int compIdx)
249 {
250 Scalar sumTi(insideTi);
251 Scalar sumXTi(insideTi*insideX);
252
253 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
254 {
255 const auto outsideScvIdx = scvf.outsideScvIdx(i);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
258 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
259
260 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
261 sumTi += outsideTi;
262 sumXTi += outsideTi*massOrMoleFractionOutside;
263 }
264
265 return sumTi > 0 ? sumXTi/sumTi : 0;
266 }
267
269 static Scalar branchingFacetDensity(const ElementVolumeVariables& elemVolVars,
270 const SubControlVolumeFace& scvf,
271 const int phaseIdx,
272 const Scalar insideRho)
273 {
274 Scalar rho(insideRho);
275 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
276 {
277 const auto outsideScvIdx = scvf.outsideScvIdx(i);
278 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
279 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
280 rho += rhoOutside;
281 }
282 return rho/(scvf.numOutsideScvs()+1);
283 }
284};
285
286} // end namespace Dumux
287
288#endif
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...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:47
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
Definition: adapt.hh:29
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:51
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const int compIdx)
compute diffusive transmissibilities
Definition: flux/cctpfa/fickslaw.hh:182
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:123
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:127
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/cctpfa/fickslaw.hh:137
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...