version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
14
15#include <dune/common/fvector.hh>
16
19
24
26
27namespace Dumux {
28
29// forward declaration
30template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
31class FicksLawImplementation;
32
37template <class TypeTag, ReferenceSystemFormulation referenceSystem>
38class FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
39{
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
46 using Extrusion = Extrusion_t<GridGeometry>;
47 using GridView = typename GridGeometry::GridView;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
50 using Element = typename GridView::template Codim<0>::Entity;
52 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
53 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
55
57 static const int dim = GridView::dimension;
58 static const int dimWorld = GridView::dimensionworld;
59 static const int numPhases = ModelTraits::numFluidPhases();
60 static const int numComponents = ModelTraits::numFluidComponents();
61
62 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
63
65 class TpfaFicksLawCacheFiller
66 {
67 public:
70 template<class FluxVariablesCacheFiller>
71 static void fill(FluxVariablesCache& scvfFluxVarsCache,
72 unsigned int phaseIdx, unsigned int compIdx,
73 const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolumeFace& scvf,
78 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
79 {
80 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
81 }
82 };
83
85 class TpfaFicksLawCache
86 {
87 public:
88 using Filler = TpfaFicksLawCacheFiller;
89
90 void updateDiffusion(const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const SubControlVolumeFace &scvf,
95 const unsigned int phaseIdx,
96 const unsigned int compIdx)
97 {
98 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
99 }
100
101 const Scalar& diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
102 { return tij_[phaseIdx][compIdx]; }
103
104 private:
105 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
106 };
107
108public:
111 static constexpr DiscretizationMethod discMethod{};
112
115 { return referenceSystem; }
116
118 using Cache = TpfaFicksLawCache;
119
121
132 static ComponentFluxVector flux(const Problem& problem,
133 const Element& element,
134 const FVElementGeometry& fvGeometry,
135 const VolumeVariables& insideVolVars,
136 const VolumeVariables& outsideVolVars,
137 const SubControlVolumeFace& scvf,
138 const int phaseIdx,
139 const ElementFluxVariablesCache& elemFluxVarsCache)
140 {
141 if constexpr (isMixedDimensional_)
142 if (scvf.numOutsideScv() != 1)
143 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
144
145 // helper lambda to get the outside mole or mass fraction
146 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
147 {
148 return massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
149 };
150
151 // helper lambda to get the averaged density at the scvf
152 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
153 {
154 return 0.5*(rhoInside + rhoOutside);
155 };
156
157 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
158 }
159
160
167 static ComponentFluxVector flux(const Problem& problem,
168 const Element& element,
169 const FVElementGeometry& fvGeometry,
170 const ElementVolumeVariables& elemVolVars,
171 const SubControlVolumeFace& scvf,
172 const int phaseIdx,
173 const ElementFluxVariablesCache& elemFluxVarsCache)
174 {
175 // get inside/outside volume variables
176 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
177 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
178
179 // helper lambda to get the outside mole or mass fraction
180 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
181 {
182 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
183 if constexpr (isMixedDimensional_)
184 {
185 return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
186 : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
187 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
188 }
189 else
190 return massOrMoleFractionOutside;
191 };
192
193 // helper lambda to get the averaged density at the scvf
194 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
195 {
196 if constexpr (isMixedDimensional_)
197 {
198 return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
199 : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
200 }
201 else
202 return 0.5*(rhoInside + rhoOutside);
203 };
204
205 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
206 }
207
209 static Scalar calculateTransmissibility(const Problem& problem,
210 const Element& element,
211 const FVElementGeometry& fvGeometry,
212 const ElementVolumeVariables& elemVolVars,
213 const SubControlVolumeFace& scvf,
214 const int phaseIdx, const int compIdx)
215 {
216
217
218 const auto insideScvIdx = scvf.insideScvIdx();
219 const auto& insideScv = fvGeometry.scv(insideScvIdx);
220 const auto& insideVolVars = elemVolVars[insideScvIdx];
221 const auto getDiffCoeff = [&](const auto& vv)
222 {
223 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
224 if constexpr (FluidSystem::isTracerFluidSystem())
225 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
226 else
227 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
228 };
229
230 const auto insideDiffCoeff = getDiffCoeff(insideVolVars);
231
232 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideDiffCoeff, insideVolVars.extrusionFactor());
233
234 // for the boundary (dirichlet) or at branching points we only need ti
235 Scalar tij;
236 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
237 tij = Extrusion::area(fvGeometry, scvf)*ti;
238
239 // otherwise we compute a tpfa harmonic mean
240 else
241 {
242 const auto outsideScvIdx = scvf.outsideScvIdx();
243 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
244 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
245 const auto outsideDiffCoeff = getDiffCoeff(outsideVolVars);
246
247 Scalar tj;
248 if constexpr (dim == dimWorld)
249 // assume the normal vector from outside is anti parallel so we save flipping a vector
250 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, outsideDiffCoeff, outsideVolVars.extrusionFactor());
251 else
252 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()),
253 outsideScv,
254 outsideDiffCoeff,
255 outsideVolVars.extrusionFactor());
256
257 // check if we are dividing by zero!
258 if (ti*tj <= 0.0)
259 tij = 0;
260 else
261 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
262 }
263
264 return tij;
265 }
266
267private:
268 template<class OutsideFractionHelper, class DensityHelper>
269 static ComponentFluxVector flux_(const Element& element,
270 const FVElementGeometry& fvGeometry,
271 const VolumeVariables& insideVolVars,
272 const VolumeVariables& outsideVolVars,
273 const ElementFluxVariablesCache& elemFluxVarsCache,
274 const SubControlVolumeFace& scvf,
275 const int phaseIdx,
276 const OutsideFractionHelper& getOutsideX,
277 const DensityHelper& getRho)
278 {
279 ComponentFluxVector componentFlux(0.0);
280 for (int compIdx = 0; compIdx < numComponents; compIdx++)
281 {
282 using FluidSystem = typename VolumeVariables::FluidSystem;
283 if constexpr (!FluidSystem::isTracerFluidSystem())
284 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
285 continue;
286
287 // diffusion tensors are always solution dependent
288 const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
289
290 // the inside and outside mass/mole fractions fractions
291 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
292 const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
293
294 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
295 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
296
297 const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
298
299 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
300 if constexpr (!FluidSystem::isTracerFluidSystem())
301 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
302 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
303 }
304
305 return componentFlux;
306 }
307
309 static Scalar branchingFacetX_(const Problem& problem,
310 const Element& element,
311 const FVElementGeometry& fvGeometry,
312 const ElementVolumeVariables& elemVolVars,
313 const ElementFluxVariablesCache& elemFluxVarsCache,
314 const SubControlVolumeFace& scvf,
315 const Scalar insideX, const Scalar insideTi,
316 const int phaseIdx, const int compIdx)
317 {
318 Scalar sumTi(insideTi);
319 Scalar sumXTi(insideTi*insideX);
320
321 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
322 {
323 const auto outsideScvIdx = scvf.outsideScvIdx(i);
324 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
325 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
326 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
327
328 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
329 sumTi += outsideTi;
330 sumXTi += outsideTi*massOrMoleFractionOutside;
331 }
332
333 return sumTi > 0 ? sumXTi/sumTi : 0;
334 }
335
337 static Scalar branchingFacetDensity_(const ElementVolumeVariables& elemVolVars,
338 const SubControlVolumeFace& scvf,
339 const int phaseIdx,
340 const Scalar insideRho)
341 {
342 Scalar rho(insideRho);
343 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
344 {
345 const auto outsideScvIdx = scvf.outsideScvIdx(i);
346 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
347 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
348 rho += rhoOutside;
349 }
350 return rho/(scvf.numOutsideScvs()+1);
351 }
352
353 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
354};
355
356} // end namespace Dumux
357
358#endif
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:39
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:167
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:114
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:209
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:118
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:132
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...
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::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:36
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
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: method.hh:25
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...