3.4
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 VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
62 using Element = typename GridView::template Codim<0>::Entity;
63 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
141 static ComponentFluxVector flux(const Problem& problem,
142 const Element& element,
143 const FVElementGeometry& fvGeometry,
144 const VolumeVariables& insideVolVars,
145 const VolumeVariables& outsideVolVars,
146 const SubControlVolumeFace& scvf,
147 const int phaseIdx,
148 const ElementFluxVariablesCache& elemFluxVarsCache)
149 {
150 if constexpr (isMixedDimensional_)
151 if (scvf.numOutsideScv() != 1)
152 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
153
154 // helper lambda to get the outside mole or mass fraction
155 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
156 {
157 return massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
158 };
159
160 // helper lambda to get the averaged density at the scvf
161 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
162 {
163 return 0.5*(rhoInside + rhoOutside);
164 };
165
166 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
167 }
168
169
176 static ComponentFluxVector flux(const Problem& problem,
177 const Element& element,
178 const FVElementGeometry& fvGeometry,
179 const ElementVolumeVariables& elemVolVars,
180 const SubControlVolumeFace& scvf,
181 const int phaseIdx,
182 const ElementFluxVariablesCache& elemFluxVarsCache)
183 {
184 // get inside/outside volume variables
185 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
186 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
187
188 // helper lambda to get the outside mole or mass fraction
189 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
190 {
191 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
192 if constexpr (isMixedDimensional_)
193 {
194 return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
195 : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
196 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
197 }
198 else
199 return massOrMoleFractionOutside;
200 };
201
202 // helper lambda to get the averaged density at the scvf
203 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
204 {
205 if constexpr (isMixedDimensional_)
206 {
207 return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
208 : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
209 }
210 else
211 return 0.5*(rhoInside + rhoOutside);
212 };
213
214 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
215 }
216
218 static Scalar calculateTransmissibility(const Problem& problem,
219 const Element& element,
220 const FVElementGeometry& fvGeometry,
221 const ElementVolumeVariables& elemVolVars,
222 const SubControlVolumeFace& scvf,
223 const int phaseIdx, const int compIdx)
224 {
225
226
227 const auto insideScvIdx = scvf.insideScvIdx();
228 const auto& insideScv = fvGeometry.scv(insideScvIdx);
229 const auto& insideVolVars = elemVolVars[insideScvIdx];
230 const auto getDiffCoeff = [&](const auto& vv)
231 {
232 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
233 if constexpr (FluidSystem::isTracerFluidSystem())
234 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
235 else
236 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
237 };
238
239 const auto insideD = getDiffCoeff(insideVolVars);
240
241 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideD, insideVolVars.extrusionFactor());
242
243 // for the boundary (dirichlet) or at branching points we only need ti
244 Scalar tij;
245 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
246 tij = Extrusion::area(scvf)*ti;
247
248 // otherwise we compute a tpfa harmonic mean
249 else
250 {
251 const auto outsideScvIdx = scvf.outsideScvIdx();
252 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
253 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
254 const auto outsideD = getDiffCoeff(outsideVolVars);
255
256 Scalar tj;
257 if constexpr (dim == dimWorld)
258 // assume the normal vector from outside is anti parallel so we save flipping a vector
259 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideD, outsideVolVars.extrusionFactor());
260 else
261 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()),
262 outsideScv,
263 outsideD,
264 outsideVolVars.extrusionFactor());
265
266 // check if we are dividing by zero!
267 if (ti*tj <= 0.0)
268 tij = 0;
269 else
270 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
271 }
272
273 return tij;
274 }
275
276private:
277 template<class OutsideFractionHelper, class DensityHelper>
278 static ComponentFluxVector flux_(const Element& element,
279 const FVElementGeometry& fvGeometry,
280 const VolumeVariables& insideVolVars,
281 const VolumeVariables& outsideVolVars,
282 const ElementFluxVariablesCache& elemFluxVarsCache,
283 const SubControlVolumeFace& scvf,
284 const int phaseIdx,
285 const OutsideFractionHelper& getOutsideX,
286 const DensityHelper& getRho)
287 {
288 ComponentFluxVector componentFlux(0.0);
289 for (int compIdx = 0; compIdx < numComponents; compIdx++)
290 {
291 using FluidSystem = typename VolumeVariables::FluidSystem;
292 if constexpr (!FluidSystem::isTracerFluidSystem())
293 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
294 continue;
295
296 // diffusion tensors are always solution dependent
297 const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
298
299 // the inside and outside mass/mole fractions fractions
300 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
301 const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
302
303 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
304 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
305
306 const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
307
308 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
309 if constexpr (!FluidSystem::isTracerFluidSystem())
310 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
311 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
312 }
313
314 return componentFlux;
315 }
316
318 static Scalar branchingFacetX_(const Problem& problem,
319 const Element& element,
320 const FVElementGeometry& fvGeometry,
321 const ElementVolumeVariables& elemVolVars,
322 const ElementFluxVariablesCache& elemFluxVarsCache,
323 const SubControlVolumeFace& scvf,
324 const Scalar insideX, const Scalar insideTi,
325 const int phaseIdx, const int compIdx)
326 {
327 Scalar sumTi(insideTi);
328 Scalar sumXTi(insideTi*insideX);
329
330 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
331 {
332 const auto outsideScvIdx = scvf.outsideScvIdx(i);
333 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
334 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
335 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
336
337 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
338 sumTi += outsideTi;
339 sumXTi += outsideTi*massOrMoleFractionOutside;
340 }
341
342 return sumTi > 0 ? sumXTi/sumTi : 0;
343 }
344
346 static Scalar branchingFacetDensity_(const ElementVolumeVariables& elemVolVars,
347 const SubControlVolumeFace& scvf,
348 const int phaseIdx,
349 const Scalar insideRho)
350 {
351 Scalar rho(insideRho);
352 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
353 {
354 const auto outsideScvIdx = scvf.outsideScvIdx(i);
355 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
356 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
357 rho += rhoOutside;
358 }
359 return rho/(scvf.numOutsideScvs()+1);
360 }
361
362 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
363};
364
365} // end namespace Dumux
366
367#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
Definition: propertysystem.hh:150
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:218
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:123
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &insideVolVars, const VolumeVariables &outsideVolVars, 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:141
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:176
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...