3.3.0
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, DiscretizationMethod::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
158 { return referenceSystem; }
159
161 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
162 static ComponentFluxVector flux(const Problem& problem,
163 const Element& element,
164 const FVElementGeometry& fvGeometry,
165 const ElementVolumeVariables& elemVolVars,
166 const SubControlVolumeFace& scvf,
167 int phaseIdx,
168 const ElementFluxVarsCache& elemFluxVarsCache)
169 {
170 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
171 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
172
173 ComponentFluxVector componentFlux(0.0);
174 for (int compIdx = 0; compIdx < numComponents; compIdx++)
175 {
176 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
177 continue;
178
179 // get inside/outside volume variables
180 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
181 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
182 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
183 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
184
185 // the inside and outside mass/mole fractions fractions
186 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
187 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
188 const Scalar xOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
189
190 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
191 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
192 const Scalar rho = 0.5*(rhoInside + rhoFacet);
193
194 componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
195 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
196
197 if (!scvf.boundary())
198 componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
199 componentFlux[compIdx] *= rho;
200
201 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
202 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
203 }
204
205 return componentFlux;
206 }
207
208 // The flux variables cache has to be bound to an element prior to flux calculations
209 // During the binding, the transmissibility will be computed and stored using the method below.
210 template< class Problem, class ElementVolumeVariables >
211 static typename Cache::DiffusionTransmissibilityContainer
212 calculateTransmissibility(const Problem& problem,
213 const Element& element,
214 const FVElementGeometry& fvGeometry,
215 const ElementVolumeVariables& elemVolVars,
216 const SubControlVolumeFace& scvf,
217 unsigned int phaseIdx, unsigned int compIdx)
218 {
219 typename Cache::DiffusionTransmissibilityContainer tij;
220 if (!problem.couplingManager().isCoupled(element, scvf))
221 {
223 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
224 return tij;
225 }
226
228 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
229
230 const auto insideScvIdx = scvf.insideScvIdx();
231 const auto& insideScv = fvGeometry.scv(insideScvIdx);
232 const auto& insideVolVars = elemVolVars[insideScvIdx];
233 const auto wIn = Extrusion::area(scvf)
234 *computeTpfaTransmissibility(scvf, insideScv,
235 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
236 insideVolVars.extrusionFactor());
237
238 // proceed depending on the interior BC types used
239 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
240
241 // neumann-coupling
242 if (iBcTypes.hasOnlyNeumann())
243 {
244 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
245 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
246 /facetVolVars.extrusionFactor()
247 *vtmv(scvf.unitOuterNormal(),
248 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
249 scvf.unitOuterNormal());
250
251 // The fluxes across this face and the outside face can be expressed in matrix form:
252 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
253 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
254 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
255 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
256 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
257 // that allow the description of the fluxes as functions of the cell and Dirichlet mass/mole fractions only.
258 if (!scvf.boundary())
259 {
260 const auto outsideScvIdx = scvf.outsideScvIdx();
261 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
262 const auto wOut = -1.0*Extrusion::area(scvf)
263 *computeTpfaTransmissibility(scvf, fvGeometry.scv(outsideScvIdx),
264 outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
265 outsideVolVars.extrusionFactor());
266
267 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
268 {
269 // optimized implementation: factorization obtained using sympy
270 // see CCTpfaFacetCouplingDarcysLaw for more details
271 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
272 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
273 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
274 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
275 }
276 else
277 {
278 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
279 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
280 tij[Cache::outsideTijIdx] = 0.0;
281 }
282 }
283 else
284 {
285 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
286 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
287 tij[Cache::outsideTijIdx] = 0.0;
288 }
289 }
290 else if (iBcTypes.hasOnlyDirichlet())
291 {
292 tij[Cache::insideTijIdx] = wIn;
293 tij[Cache::outsideTijIdx] = 0.0;
294 tij[Cache::facetTijIdx] = -wIn;
295 }
296 else
297 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
298
299 return tij;
300 }
301};
302
307template<class TypeTag, ReferenceSystemFormulation referenceSystem>
308class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/true>
309: public FicksLawImplementation<TypeTag, DiscretizationMethod::cctpfa, referenceSystem>
310{
313
315 using FVElementGeometry = typename GridGeometry::LocalView;
316 using SubControlVolume = typename GridGeometry::SubControlVolume;
317 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
318 using Extrusion = Extrusion_t<GridGeometry>;
319
320 using GridView = typename GridGeometry::GridView;
321 using Element = typename GridView::template Codim<0>::Entity;
322
326
327 static const int numPhases = ModelTraits::numFluidPhases();
328 static const int numComponents = ModelTraits::numFluidComponents();
329
331 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
332
336 class FacetCouplingFicksLawCache
337 {
338 public:
340 using Filler = typename ParentType::Cache::Filler;
341
345 static constexpr int insideTijIdx = 0;
346 static constexpr int facetTijIdx = 1;
347
349 using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
350
352 template< class Problem, class ElementVolumeVariables >
353 void updateDiffusion(const Problem& problem,
354 const Element& element,
355 const FVElementGeometry& fvGeometry,
356 const ElementVolumeVariables& elemVolVars,
357 const SubControlVolumeFace &scvf,
358 unsigned int phaseIdx,
359 unsigned int compIdx)
360 {
361 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
362 }
363
368 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
369 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
370
372 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
373 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
374
376 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
377 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
378
379 private:
380 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
381 };
382
383public:
385 using Cache = FacetCouplingFicksLawCache;
386
389
392 { return referenceSystem; }
393
395 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
396 static ComponentFluxVector flux(const Problem& problem,
397 const Element& element,
398 const FVElementGeometry& fvGeometry,
399 const ElementVolumeVariables& elemVolVars,
400 const SubControlVolumeFace& scvf,
401 int phaseIdx,
402 const ElementFluxVarsCache& elemFluxVarsCache)
403 {
404 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
405 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
406
407 ComponentFluxVector componentFlux(0.0);
408 for (int compIdx = 0; compIdx < numComponents; compIdx++)
409 {
410 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
411 continue;
412
413 // get inside/outside volume variables
414 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
415 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
416 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
417
418 // the inside and outside mass/mole fractions fractions
419 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
420 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
421
422 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
423 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
424 const Scalar rho = 0.5*(rhoInside + rhoFacet);
425
426 componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
427 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
428
429 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
430 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
431 }
432
433 return componentFlux;
434 }
435
436 // The flux variables cache has to be bound to an element prior to flux calculations
437 // During the binding, the transmissibility will be computed and stored using the method below.
438 template< class Problem, class ElementVolumeVariables >
439 static typename Cache::DiffusionTransmissibilityContainer
440 calculateTransmissibility(const Problem& problem,
441 const Element& element,
442 const FVElementGeometry& fvGeometry,
443 const ElementVolumeVariables& elemVolVars,
444 const SubControlVolumeFace& scvf,
445 unsigned int phaseIdx, unsigned int compIdx)
446 {
447 typename Cache::DiffusionTransmissibilityContainer tij;
448 if (!problem.couplingManager().isCoupled(element, scvf))
449 {
451 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
452 return tij;
453 }
454
456 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
457
458 // On surface grids only xi = 1.0 can be used, as the coupling condition
459 // for xi != 1.0 does not generalize for surface grids where the normal
460 // vectors of the inside/outside elements have different orientations.
461 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
462 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
463
464 const auto insideScvIdx = scvf.insideScvIdx();
465 const auto& insideScv = fvGeometry.scv(insideScvIdx);
466 const auto& insideVolVars = elemVolVars[insideScvIdx];
467 const auto wIn = Extrusion::area(scvf)
468 *computeTpfaTransmissibility(scvf, insideScv,
469 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
470 insideVolVars.extrusionFactor());
471
472 // proceed depending on the interior BC types used
473 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
474
475 // neumann-coupling
476 if (iBcTypes.hasOnlyNeumann())
477 {
478 // Here we use the square root of the facet extrusion factor
479 // as an approximate average distance from scvf ip to facet center
480 using std::sqrt;
481 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
482 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
483 /sqrt(facetVolVars.extrusionFactor())
484 *vtmv(scvf.unitOuterNormal(),
485 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
486 scvf.unitOuterNormal());
487
488 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
489 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
490 }
491 else if (iBcTypes.hasOnlyDirichlet())
492 {
493 tij[Cache::insideTijIdx] = wIn;
494 tij[Cache::facetTijIdx] = -wIn;
495 }
496 else
497 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
498
499 return tij;
500 }
501};
502
503} // end namespace Dumux
504
505#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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:849
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 (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:127
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:162
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:212
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:157
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:310
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:391
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:440
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:396
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.