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