3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/cellcentered/tpfa/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_FACET_COUPLING_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH
26
27#include <array>
28#include <cmath>
29
30#include <dune/common/float_cmp.hh>
31#include <dune/common/fvector.hh>
32
33#include <dumux/common/math.hh>
36
40
43
44namespace Dumux {
45
47template<class TypeTag,
48 ReferenceSystemFormulation referenceSystem,
49 bool isNetwork>
51
58template<class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
60 CCTpfaFacetCouplingFicksLawImpl< TypeTag, referenceSystem,
63
68template<class TypeTag, ReferenceSystemFormulation referenceSystem>
69class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/false>
70: public FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
71{
74
76 using FVElementGeometry = typename GridGeometry::LocalView;
77 using SubControlVolume = typename GridGeometry::SubControlVolume;
78 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
79 using Extrusion = Extrusion_t<GridGeometry>;
80
81 using GridView = typename GridGeometry::GridView;
82 using Element = typename GridView::template Codim<0>::Entity;
83
87
88 static const int numPhases = ModelTraits::numFluidPhases();
89 static const int numComponents = ModelTraits::numFluidComponents();
90
92 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
93
97 class FacetCouplingFicksLawCache
98 {
99 public:
101 using Filler = typename ParentType::Cache::Filler;
102
106 static constexpr int insideTijIdx = 0;
107 static constexpr int outsideTijIdx = 1;
108 static constexpr int facetTijIdx = 2;
109
111 using DiffusionTransmissibilityContainer = std::array<Scalar, 3>;
112
114 template< class Problem, class ElementVolumeVariables >
115 void updateDiffusion(const Problem& problem,
116 const Element& element,
117 const FVElementGeometry& fvGeometry,
118 const ElementVolumeVariables& elemVolVars,
119 const SubControlVolumeFace &scvf,
120 unsigned int phaseIdx,
121 unsigned int compIdx)
122 {
123 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
124 }
125
130 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
131 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
132
134 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
135 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
136
138 Scalar diffusionTijOutside(unsigned int phaseIdx, unsigned int compIdx) const
139 { return tij_[phaseIdx][compIdx][outsideTijIdx]; }
140
142 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
143 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
144
145 private:
146 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
147 };
148
149public:
151 using Cache = FacetCouplingFicksLawCache;
152
155 static constexpr DiscretizationMethod discMethod{};
156
159 { return referenceSystem; }
160
162 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
163 static ComponentFluxVector flux(const Problem& problem,
164 const Element& element,
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const SubControlVolumeFace& scvf,
168 int phaseIdx,
169 const ElementFluxVarsCache& elemFluxVarsCache)
170 {
171 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
172 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
173
174 ComponentFluxVector componentFlux(0.0);
175 for (int compIdx = 0; compIdx < numComponents; compIdx++)
176 {
177 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
178 continue;
179
180 // get inside/outside volume variables
181 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
182 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
183 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
184 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
185
186 // the inside and outside mass/mole fractions fractions
187 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
188 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
189 const Scalar xOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
190
191 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
192 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
193 const Scalar rho = 0.5*(rhoInside + rhoFacet);
194
195 componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
196 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
197
198 if (!scvf.boundary())
199 componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
200 componentFlux[compIdx] *= rho;
201
202 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
203 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
204 }
205
206 return componentFlux;
207 }
208
209 // The flux variables cache has to be bound to an element prior to flux calculations
210 // During the binding, the transmissibility will be computed and stored using the method below.
211 template< class Problem, class ElementVolumeVariables >
212 static typename Cache::DiffusionTransmissibilityContainer
213 calculateTransmissibility(const Problem& problem,
214 const Element& element,
215 const FVElementGeometry& fvGeometry,
216 const ElementVolumeVariables& elemVolVars,
217 const SubControlVolumeFace& scvf,
218 unsigned int phaseIdx, unsigned int compIdx)
219 {
220 typename Cache::DiffusionTransmissibilityContainer tij;
221 if (!problem.couplingManager().isCoupled(element, scvf))
222 {
224 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
225 return tij;
226 }
227
229 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
230
231 const auto insideScvIdx = scvf.insideScvIdx();
232 const auto& insideScv = fvGeometry.scv(insideScvIdx);
233 const auto& insideVolVars = elemVolVars[insideScvIdx];
234 const auto wIn = Extrusion::area(fvGeometry, scvf)
235 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
236 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
237 insideVolVars.extrusionFactor());
238
239 // proceed depending on the interior BC types used
240 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
241
242 // neumann-coupling
243 if (iBcTypes.hasOnlyNeumann())
244 {
245 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
246 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
247 /facetVolVars.extrusionFactor()
248 *vtmv(scvf.unitOuterNormal(),
249 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
250 scvf.unitOuterNormal());
251
252 // The fluxes across this face and the outside face can be expressed in matrix form:
253 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
254 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
255 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
256 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
257 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
258 // that allow the description of the fluxes as functions of the cell and Dirichlet mass/mole fractions only.
259 if (!scvf.boundary())
260 {
261 const auto outsideScvIdx = scvf.outsideScvIdx();
262 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
263 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
264 *computeTpfaTransmissibility(fvGeometry, scvf, fvGeometry.scv(outsideScvIdx),
265 outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
266 outsideVolVars.extrusionFactor());
267
268 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
269 {
270 // optimized implementation: factorization obtained using sympy
271 // see CCTpfaFacetCouplingDarcysLaw for more details
272 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
273 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
274 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
275 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
276 }
277 else
278 {
279 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
280 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
281 tij[Cache::outsideTijIdx] = 0.0;
282 }
283 }
284 else
285 {
286 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
287 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
288 tij[Cache::outsideTijIdx] = 0.0;
289 }
290 }
291 else if (iBcTypes.hasOnlyDirichlet())
292 {
293 tij[Cache::insideTijIdx] = wIn;
294 tij[Cache::outsideTijIdx] = 0.0;
295 tij[Cache::facetTijIdx] = -wIn;
296 }
297 else
298 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
299
300 return tij;
301 }
302};
303
308template<class TypeTag, ReferenceSystemFormulation referenceSystem>
309class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/true>
310: public FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
311{
314
316 using FVElementGeometry = typename GridGeometry::LocalView;
317 using SubControlVolume = typename GridGeometry::SubControlVolume;
318 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
319 using Extrusion = Extrusion_t<GridGeometry>;
320
321 using GridView = typename GridGeometry::GridView;
322 using Element = typename GridView::template Codim<0>::Entity;
323
327
328 static const int numPhases = ModelTraits::numFluidPhases();
329 static const int numComponents = ModelTraits::numFluidComponents();
330
332 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
333
337 class FacetCouplingFicksLawCache
338 {
339 public:
341 using Filler = typename ParentType::Cache::Filler;
342
346 static constexpr int insideTijIdx = 0;
347 static constexpr int facetTijIdx = 1;
348
350 using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
351
353 template< class Problem, class ElementVolumeVariables >
354 void updateDiffusion(const Problem& problem,
355 const Element& element,
356 const FVElementGeometry& fvGeometry,
357 const ElementVolumeVariables& elemVolVars,
358 const SubControlVolumeFace &scvf,
359 unsigned int phaseIdx,
360 unsigned int compIdx)
361 {
362 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
363 }
364
369 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
370 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
371
373 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
374 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
375
377 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
378 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
379
380 private:
381 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
382 };
383
384public:
386 using Cache = FacetCouplingFicksLawCache;
387
390 static constexpr DiscretizationMethod discMethod{};
391
394 { return referenceSystem; }
395
397 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
398 static ComponentFluxVector flux(const Problem& problem,
399 const Element& element,
400 const FVElementGeometry& fvGeometry,
401 const ElementVolumeVariables& elemVolVars,
402 const SubControlVolumeFace& scvf,
403 int phaseIdx,
404 const ElementFluxVarsCache& elemFluxVarsCache)
405 {
406 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
407 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
408
409 ComponentFluxVector componentFlux(0.0);
410 for (int compIdx = 0; compIdx < numComponents; compIdx++)
411 {
412 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
413 continue;
414
415 // get inside/outside volume variables
416 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
417 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
418 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
419
420 // the inside and outside mass/mole fractions fractions
421 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
422 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
423
424 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
425 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
426 const Scalar rho = 0.5*(rhoInside + rhoFacet);
427
428 componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
429 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
430
431 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
432 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
433 }
434
435 return componentFlux;
436 }
437
438 // The flux variables cache has to be bound to an element prior to flux calculations
439 // During the binding, the transmissibility will be computed and stored using the method below.
440 template< class Problem, class ElementVolumeVariables >
441 static typename Cache::DiffusionTransmissibilityContainer
442 calculateTransmissibility(const Problem& problem,
443 const Element& element,
444 const FVElementGeometry& fvGeometry,
445 const ElementVolumeVariables& elemVolVars,
446 const SubControlVolumeFace& scvf,
447 unsigned int phaseIdx, unsigned int compIdx)
448 {
449 typename Cache::DiffusionTransmissibilityContainer tij;
450 if (!problem.couplingManager().isCoupled(element, scvf))
451 {
453 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
454 return tij;
455 }
456
458 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
459
460 // On surface grids only xi = 1.0 can be used, as the coupling condition
461 // for xi != 1.0 does not generalize for surface grids where the normal
462 // vectors of the inside/outside elements have different orientations.
463 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
464 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
465
466 const auto insideScvIdx = scvf.insideScvIdx();
467 const auto& insideScv = fvGeometry.scv(insideScvIdx);
468 const auto& insideVolVars = elemVolVars[insideScvIdx];
469 const auto wIn = Extrusion::area(fvGeometry, scvf)
470 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
471 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
472 insideVolVars.extrusionFactor());
473
474 // proceed depending on the interior BC types used
475 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
476
477 // neumann-coupling
478 if (iBcTypes.hasOnlyNeumann())
479 {
480 // Here we use the square root of the facet extrusion factor
481 // as an approximate average distance from scvf ip to facet center
482 using std::sqrt;
483 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
484 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
485 /sqrt(facetVolVars.extrusionFactor())
486 *vtmv(scvf.unitOuterNormal(),
487 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
488 scvf.unitOuterNormal());
489
490 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
491 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
492 }
493 else if (iBcTypes.hasOnlyDirichlet())
494 {
495 tij[Cache::insideTijIdx] = wIn;
496 tij[Cache::facetTijIdx] = -wIn;
497 }
498 else
499 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
500
501 return tij;
502 }
503};
504
505} // end namespace Dumux
506
507#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
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...
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
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
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
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:130
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:50
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:71
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the diffusive fluxes.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:163
static Cache::DiffusionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, unsigned int phaseIdx, unsigned int compIdx)
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:213
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:158
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:311
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:393
static Cache::DiffusionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, unsigned int phaseIdx, unsigned int compIdx)
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:442
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the diffusive fluxes.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:398
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Fick's law for cell-centered finite volume schemes with two-point flux approximation.