3.6-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;
64 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
65 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
67
69 static const int dim = GridView::dimension;
70 static const int dimWorld = GridView::dimensionworld;
71 static const int numPhases = ModelTraits::numFluidPhases();
72 static const int numComponents = ModelTraits::numFluidComponents();
73
74 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
75
77 class TpfaFicksLawCacheFiller
78 {
79 public:
82 template<class FluxVariablesCacheFiller>
83 static void fill(FluxVariablesCache& scvfFluxVarsCache,
84 unsigned int phaseIdx, unsigned int compIdx,
85 const Problem& problem,
86 const Element& element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const SubControlVolumeFace& scvf,
90 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
91 {
92 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
93 }
94 };
95
97 class TpfaFicksLawCache
98 {
99 public:
100 using Filler = TpfaFicksLawCacheFiller;
101
102 void updateDiffusion(const Problem& problem,
103 const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const ElementVolumeVariables& elemVolVars,
106 const SubControlVolumeFace &scvf,
107 const unsigned int phaseIdx,
108 const unsigned int compIdx)
109 {
110 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
111 }
112
113 const Scalar& diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
114 { return tij_[phaseIdx][compIdx]; }
115
116 private:
117 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
118 };
119
120public:
123 static constexpr DiscretizationMethod discMethod{};
124
127 { return referenceSystem; }
128
130 using Cache = TpfaFicksLawCache;
131
133
144 static ComponentFluxVector flux(const Problem& problem,
145 const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const VolumeVariables& insideVolVars,
148 const VolumeVariables& outsideVolVars,
149 const SubControlVolumeFace& scvf,
150 const int phaseIdx,
151 const ElementFluxVariablesCache& elemFluxVarsCache)
152 {
153 if constexpr (isMixedDimensional_)
154 if (scvf.numOutsideScv() != 1)
155 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
156
157 // helper lambda to get the outside mole or mass fraction
158 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
159 {
160 return massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
161 };
162
163 // helper lambda to get the averaged density at the scvf
164 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
165 {
166 return 0.5*(rhoInside + rhoOutside);
167 };
168
169 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
170 }
171
172
179 static ComponentFluxVector flux(const Problem& problem,
180 const Element& element,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& elemVolVars,
183 const SubControlVolumeFace& scvf,
184 const int phaseIdx,
185 const ElementFluxVariablesCache& elemFluxVarsCache)
186 {
187 // get inside/outside volume variables
188 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
189 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
190
191 // helper lambda to get the outside mole or mass fraction
192 const auto getOutsideX = [&](const Scalar xInside, const Scalar tij, const int phaseIdx, const int compIdx)
193 {
194 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
195 if constexpr (isMixedDimensional_)
196 {
197 return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
198 : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
199 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
200 }
201 else
202 return massOrMoleFractionOutside;
203 };
204
205 // helper lambda to get the averaged density at the scvf
206 const auto getRho = [&](const int phaseIdx, const Scalar rhoInside, const Scalar rhoOutside)
207 {
208 if constexpr (isMixedDimensional_)
209 {
210 return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
211 : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
212 }
213 else
214 return 0.5*(rhoInside + rhoOutside);
215 };
216
217 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
218 }
219
221 static Scalar calculateTransmissibility(const Problem& problem,
222 const Element& element,
223 const FVElementGeometry& fvGeometry,
224 const ElementVolumeVariables& elemVolVars,
225 const SubControlVolumeFace& scvf,
226 const int phaseIdx, const int compIdx)
227 {
228
229
230 const auto insideScvIdx = scvf.insideScvIdx();
231 const auto& insideScv = fvGeometry.scv(insideScvIdx);
232 const auto& insideVolVars = elemVolVars[insideScvIdx];
233 const auto getDiffCoeff = [&](const auto& vv)
234 {
235 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
236 if constexpr (FluidSystem::isTracerFluidSystem())
237 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
238 else
239 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
240 };
241
242 const auto insideDiffCoeff = getDiffCoeff(insideVolVars);
243
244 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideDiffCoeff, insideVolVars.extrusionFactor());
245
246 // for the boundary (dirichlet) or at branching points we only need ti
247 Scalar tij;
248 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
249 tij = Extrusion::area(fvGeometry, scvf)*ti;
250
251 // otherwise we compute a tpfa harmonic mean
252 else
253 {
254 const auto outsideScvIdx = scvf.outsideScvIdx();
255 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const auto outsideDiffCoeff = getDiffCoeff(outsideVolVars);
258
259 Scalar tj;
260 if constexpr (dim == dimWorld)
261 // assume the normal vector from outside is anti parallel so we save flipping a vector
262 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, outsideDiffCoeff, outsideVolVars.extrusionFactor());
263 else
264 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()),
265 outsideScv,
266 outsideDiffCoeff,
267 outsideVolVars.extrusionFactor());
268
269 // check if we are dividing by zero!
270 if (ti*tj <= 0.0)
271 tij = 0;
272 else
273 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
274 }
275
276 return tij;
277 }
278
279private:
280 template<class OutsideFractionHelper, class DensityHelper>
281 static ComponentFluxVector flux_(const Element& element,
282 const FVElementGeometry& fvGeometry,
283 const VolumeVariables& insideVolVars,
284 const VolumeVariables& outsideVolVars,
285 const ElementFluxVariablesCache& elemFluxVarsCache,
286 const SubControlVolumeFace& scvf,
287 const int phaseIdx,
288 const OutsideFractionHelper& getOutsideX,
289 const DensityHelper& getRho)
290 {
291 ComponentFluxVector componentFlux(0.0);
292 for (int compIdx = 0; compIdx < numComponents; compIdx++)
293 {
294 using FluidSystem = typename VolumeVariables::FluidSystem;
295 if constexpr (!FluidSystem::isTracerFluidSystem())
296 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
297 continue;
298
299 // diffusion tensors are always solution dependent
300 const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
301
302 // the inside and outside mass/mole fractions fractions
303 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
304 const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
305
306 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
307 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
308
309 const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
310
311 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
312 if constexpr (!FluidSystem::isTracerFluidSystem())
313 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
314 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
315 }
316
317 return componentFlux;
318 }
319
321 static Scalar branchingFacetX_(const Problem& problem,
322 const Element& element,
323 const FVElementGeometry& fvGeometry,
324 const ElementVolumeVariables& elemVolVars,
325 const ElementFluxVariablesCache& elemFluxVarsCache,
326 const SubControlVolumeFace& scvf,
327 const Scalar insideX, const Scalar insideTi,
328 const int phaseIdx, const int compIdx)
329 {
330 Scalar sumTi(insideTi);
331 Scalar sumXTi(insideTi*insideX);
332
333 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
334 {
335 const auto outsideScvIdx = scvf.outsideScvIdx(i);
336 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
337 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
338 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
339
340 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
341 sumTi += outsideTi;
342 sumXTi += outsideTi*massOrMoleFractionOutside;
343 }
344
345 return sumTi > 0 ? sumXTi/sumTi : 0;
346 }
347
349 static Scalar branchingFacetDensity_(const ElementVolumeVariables& elemVolVars,
350 const SubControlVolumeFace& scvf,
351 const int phaseIdx,
352 const Scalar insideRho)
353 {
354 Scalar rho(insideRho);
355 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
356 {
357 const auto outsideScvIdx = scvf.outsideScvIdx(i);
358 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
359 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
360 rho += rhoOutside;
361 }
362 return rho/(scvf.numOutsideScvs()+1);
363 }
364
365 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
366};
367
368} // end namespace Dumux
369
370#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:48
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: method.hh:37
forward declaration of the method-specific implementation
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:179
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:126
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:221
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:130
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:144
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...