3.2-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
39
42
43namespace Dumux {
44
46template<class TypeTag,
47 ReferenceSystemFormulation referenceSystem,
48 bool isNetwork>
50
57template<class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
59 CCTpfaFacetCouplingFicksLawImpl< TypeTag, referenceSystem,
62
67template<class TypeTag, ReferenceSystemFormulation referenceSystem>
68class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/false>
69: public FicksLawImplementation<TypeTag, DiscretizationMethod::cctpfa, referenceSystem>
70{
73
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolume = typename GridGeometry::SubControlVolume;
77 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
78
79 using GridView = typename GridGeometry::GridView;
80 using Element = typename GridView::template Codim<0>::Entity;
81
85
86 static const int numPhases = ModelTraits::numFluidPhases();
87 static const int numComponents = ModelTraits::numFluidComponents();
88
90 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
91
95 class FacetCouplingFicksLawCache
96 {
97 public:
99 using Filler = typename ParentType::Cache::Filler;
100
104 static constexpr int insideTijIdx = 0;
105 static constexpr int outsideTijIdx = 1;
106 static constexpr int facetTijIdx = 2;
107
109 using DiffusionTransmissibilityContainer = std::array<Scalar, 3>;
110
112 template< class Problem, class ElementVolumeVariables >
113 void updateDiffusion(const Problem& problem,
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& elemVolVars,
117 const SubControlVolumeFace &scvf,
118 unsigned int phaseIdx,
119 unsigned int compIdx)
120 {
121 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
122 }
123
128 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
129 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
130
132 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
133 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
134
136 Scalar diffusionTijOutside(unsigned int phaseIdx, unsigned int compIdx) const
137 { return tij_[phaseIdx][compIdx][outsideTijIdx]; }
138
140 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
141 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
142
143 private:
144 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
145 };
146
147public:
149 using Cache = FacetCouplingFicksLawCache;
150
153
156 { return referenceSystem; }
157
159 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
160 static ComponentFluxVector flux(const Problem& problem,
161 const Element& element,
162 const FVElementGeometry& fvGeometry,
163 const ElementVolumeVariables& elemVolVars,
164 const SubControlVolumeFace& scvf,
165 int phaseIdx,
166 const ElementFluxVarsCache& elemFluxVarsCache)
167 {
168 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
169 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
170
171 ComponentFluxVector componentFlux(0.0);
172 for (int compIdx = 0; compIdx < numComponents; compIdx++)
173 {
174 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
175 continue;
176
177 // get inside/outside volume variables
178 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
179 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
180 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
181 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
182
183 // the inside and outside mass/mole fractions fractions
184 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
185 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
186 const Scalar xOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
187
188 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
189 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
190 const Scalar rho = 0.5*(rhoInside + rhoFacet);
191
192 componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
193 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
194
195 if (!scvf.boundary())
196 componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
197 componentFlux[compIdx] *= rho;
198
199 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
200 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
201 }
202
203 return componentFlux;
204 }
205
206 // The flux variables cache has to be bound to an element prior to flux calculations
207 // During the binding, the transmissibility will be computed and stored using the method below.
208 template< class Problem, class ElementVolumeVariables >
209 static typename Cache::DiffusionTransmissibilityContainer
210 calculateTransmissibility(const Problem& problem,
211 const Element& element,
212 const FVElementGeometry& fvGeometry,
213 const ElementVolumeVariables& elemVolVars,
214 const SubControlVolumeFace& scvf,
215 unsigned int phaseIdx, unsigned int compIdx)
216 {
217 typename Cache::DiffusionTransmissibilityContainer tij;
218 if (!problem.couplingManager().isCoupled(element, scvf))
219 {
221 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
222 return tij;
223 }
224
226 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
227
228 const auto insideScvIdx = scvf.insideScvIdx();
229 const auto& insideScv = fvGeometry.scv(insideScvIdx);
230 const auto& insideVolVars = elemVolVars[insideScvIdx];
231 const auto wIn = scvf.area()*computeTpfaTransmissibility(scvf,
232 insideScv,
233 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
234 insideVolVars.extrusionFactor());
235
236 // proceed depending on the interior BC types used
237 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
238
239 // neumann-coupling
240 if (iBcTypes.hasOnlyNeumann())
241 {
242 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
243 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
244 /facetVolVars.extrusionFactor()
245 *vtmv(scvf.unitOuterNormal(),
246 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
247 scvf.unitOuterNormal());
248
249 // The fluxes across this face and the outside face can be expressed in matrix form:
250 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
251 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
252 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
253 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
254 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
255 // that allow the description of the fluxes as functions of the cell and Dirichlet mass/mole fractions only.
256 if (!scvf.boundary())
257 {
258 const auto outsideScvIdx = scvf.outsideScvIdx();
259 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
260 const auto wOut = -1.0*scvf.area()*computeTpfaTransmissibility(scvf,
261 fvGeometry.scv(outsideScvIdx),
262 outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
263 outsideVolVars.extrusionFactor());
264
265 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
266 {
267 // optimized implementation: factorization obtained using sympy
268 // see CCTpfaFacetCouplingDarcysLaw for more details
269 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
270 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
271 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
272 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
273 }
274 else
275 {
276 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
277 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
278 tij[Cache::outsideTijIdx] = 0.0;
279 }
280 }
281 else
282 {
283 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
284 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
285 tij[Cache::outsideTijIdx] = 0.0;
286 }
287 }
288 else if (iBcTypes.hasOnlyDirichlet())
289 {
290 tij[Cache::insideTijIdx] = wIn;
291 tij[Cache::outsideTijIdx] = 0.0;
292 tij[Cache::facetTijIdx] = -wIn;
293 }
294 else
295 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
296
297 return tij;
298 }
299};
300
305template<class TypeTag, ReferenceSystemFormulation referenceSystem>
306class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/true>
307: public FicksLawImplementation<TypeTag, DiscretizationMethod::cctpfa, referenceSystem>
308{
311
313 using FVElementGeometry = typename GridGeometry::LocalView;
314 using SubControlVolume = typename GridGeometry::SubControlVolume;
315 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
316
317 using GridView = typename GridGeometry::GridView;
318 using Element = typename GridView::template Codim<0>::Entity;
319
323
324 static const int numPhases = ModelTraits::numFluidPhases();
325 static const int numComponents = ModelTraits::numFluidComponents();
326
328 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
329
333 class FacetCouplingFicksLawCache
334 {
335 public:
337 using Filler = typename ParentType::Cache::Filler;
338
342 static constexpr int insideTijIdx = 0;
343 static constexpr int facetTijIdx = 1;
344
346 using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
347
349 template< class Problem, class ElementVolumeVariables >
350 void updateDiffusion(const Problem& problem,
351 const Element& element,
352 const FVElementGeometry& fvGeometry,
353 const ElementVolumeVariables& elemVolVars,
354 const SubControlVolumeFace &scvf,
355 unsigned int phaseIdx,
356 unsigned int compIdx)
357 {
358 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
359 }
360
365 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
366 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
367
369 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
370 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
371
373 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
374 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
375
376 private:
377 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
378 };
379
380public:
382 using Cache = FacetCouplingFicksLawCache;
383
386
389 { return referenceSystem; }
390
392 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
393 static ComponentFluxVector flux(const Problem& problem,
394 const Element& element,
395 const FVElementGeometry& fvGeometry,
396 const ElementVolumeVariables& elemVolVars,
397 const SubControlVolumeFace& scvf,
398 int phaseIdx,
399 const ElementFluxVarsCache& elemFluxVarsCache)
400 {
401 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
402 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
403
404 ComponentFluxVector componentFlux(0.0);
405 for (int compIdx = 0; compIdx < numComponents; compIdx++)
406 {
407 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
408 continue;
409
410 // get inside/outside volume variables
411 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
412 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
413 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
414
415 // the inside and outside mass/mole fractions fractions
416 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
417 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
418
419 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
420 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
421 const Scalar rho = 0.5*(rhoInside + rhoFacet);
422
423 componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
424 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
425
426 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
427 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
428 }
429
430 return componentFlux;
431 }
432
433 // The flux variables cache has to be bound to an element prior to flux calculations
434 // During the binding, the transmissibility will be computed and stored using the method below.
435 template< class Problem, class ElementVolumeVariables >
436 static typename Cache::DiffusionTransmissibilityContainer
437 calculateTransmissibility(const Problem& problem,
438 const Element& element,
439 const FVElementGeometry& fvGeometry,
440 const ElementVolumeVariables& elemVolVars,
441 const SubControlVolumeFace& scvf,
442 unsigned int phaseIdx, unsigned int compIdx)
443 {
444 typename Cache::DiffusionTransmissibilityContainer tij;
445 if (!problem.couplingManager().isCoupled(element, scvf))
446 {
448 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
449 return tij;
450 }
451
453 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
454
455 // On surface grids only xi = 1.0 can be used, as the coupling condition
456 // for xi != 1.0 does not generalize for surface grids where the normal
457 // vectors of the inside/outside elements have different orientations.
458 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
459 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
460
461 const auto insideScvIdx = scvf.insideScvIdx();
462 const auto& insideScv = fvGeometry.scv(insideScvIdx);
463 const auto& insideVolVars = elemVolVars[insideScvIdx];
464 const auto wIn = scvf.area()*computeTpfaTransmissibility(scvf,
465 insideScv,
466 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
467 insideVolVars.extrusionFactor());
468
469 // proceed depending on the interior BC types used
470 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
471
472 // neumann-coupling
473 if (iBcTypes.hasOnlyNeumann())
474 {
475 // Here we use the square root of the facet extrusion factor
476 // as an approximate average distance from scvf ip to facet center
477 using std::sqrt;
478 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
479 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
480 /sqrt(facetVolVars.extrusionFactor())
481 *vtmv(scvf.unitOuterNormal(),
482 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
483 scvf.unitOuterNormal());
484
485 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
486 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
487 }
488 else if (iBcTypes.hasOnlyDirichlet())
489 {
490 tij[Cache::insideTijIdx] = wIn;
491 tij[Cache::facetTijIdx] = -wIn;
492 }
493 else
494 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
495
496 return tij;
497 }
498};
499
500} // end namespace Dumux
501
502#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
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:840
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:125
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:49
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:70
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:160
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:210
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:155
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:308
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:388
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:437
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:393
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.