version 3.8
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// 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_FACET_COUPLING_FICKS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH
14
15#include <array>
16#include <cmath>
17
18#include <dune/common/float_cmp.hh>
19#include <dune/common/fvector.hh>
20
21#include <dumux/common/math.hh>
24
28
31
32namespace Dumux {
33
35template<class TypeTag,
36 ReferenceSystemFormulation referenceSystem,
37 bool isNetwork>
39
46template<class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
48 CCTpfaFacetCouplingFicksLawImpl< TypeTag, referenceSystem,
51
56template<class TypeTag, ReferenceSystemFormulation referenceSystem>
57class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/false>
58: public FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
59{
62
64 using FVElementGeometry = typename GridGeometry::LocalView;
65 using SubControlVolume = typename GridGeometry::SubControlVolume;
66 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
67 using Extrusion = Extrusion_t<GridGeometry>;
68
69 using GridView = typename GridGeometry::GridView;
70 using Element = typename GridView::template Codim<0>::Entity;
71
75
76 static const int numPhases = ModelTraits::numFluidPhases();
77 static const int numComponents = ModelTraits::numFluidComponents();
78
80 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
81
85 class FacetCouplingFicksLawCache
86 {
87 public:
89 using Filler = typename ParentType::Cache::Filler;
90
94 static constexpr int insideTijIdx = 0;
95 static constexpr int outsideTijIdx = 1;
96 static constexpr int facetTijIdx = 2;
97
99 using DiffusionTransmissibilityContainer = std::array<Scalar, 3>;
100
102 template< class Problem, class ElementVolumeVariables >
103 void updateDiffusion(const Problem& problem,
104 const Element& element,
105 const FVElementGeometry& fvGeometry,
106 const ElementVolumeVariables& elemVolVars,
107 const SubControlVolumeFace &scvf,
108 unsigned int phaseIdx,
109 unsigned int compIdx)
110 {
111 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
112 }
113
118 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
119 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
120
122 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
123 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
124
126 Scalar diffusionTijOutside(unsigned int phaseIdx, unsigned int compIdx) const
127 { return tij_[phaseIdx][compIdx][outsideTijIdx]; }
128
130 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
131 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
132
133 private:
134 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
135 };
136
137public:
139 using Cache = FacetCouplingFicksLawCache;
140
143 static constexpr DiscretizationMethod discMethod{};
144
147 { return referenceSystem; }
148
150 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
151 static ComponentFluxVector flux(const Problem& problem,
152 const Element& element,
153 const FVElementGeometry& fvGeometry,
154 const ElementVolumeVariables& elemVolVars,
155 const SubControlVolumeFace& scvf,
156 int phaseIdx,
157 const ElementFluxVarsCache& elemFluxVarsCache)
158 {
159 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
160 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
161
162 ComponentFluxVector componentFlux(0.0);
163 for (int compIdx = 0; compIdx < numComponents; compIdx++)
164 {
165 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
166 continue;
167
168 // get inside/outside volume variables
169 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
170 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
171 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
172 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
173
174 // the inside and outside mass/mole fractions fractions
175 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
176 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
177 const Scalar xOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
178
179 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
180 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
181 const Scalar rho = 0.5*(rhoInside + rhoFacet);
182
183 componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
184 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
185
186 if (!scvf.boundary())
187 componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
188 componentFlux[compIdx] *= rho;
189
190 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
191 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
192 }
193
194 return componentFlux;
195 }
196
197 // The flux variables cache has to be bound to an element prior to flux calculations
198 // During the binding, the transmissibility will be computed and stored using the method below.
199 template< class Problem, class ElementVolumeVariables >
200 static typename Cache::DiffusionTransmissibilityContainer
201 calculateTransmissibility(const Problem& problem,
202 const Element& element,
203 const FVElementGeometry& fvGeometry,
204 const ElementVolumeVariables& elemVolVars,
205 const SubControlVolumeFace& scvf,
206 unsigned int phaseIdx, unsigned int compIdx)
207 {
208 typename Cache::DiffusionTransmissibilityContainer tij;
209 if (!problem.couplingManager().isCoupled(element, scvf))
210 {
212 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
213 return tij;
214 }
215
217 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
218
219 const auto insideScvIdx = scvf.insideScvIdx();
220 const auto& insideScv = fvGeometry.scv(insideScvIdx);
221 const auto& insideVolVars = elemVolVars[insideScvIdx];
222 const auto wIn = Extrusion::area(fvGeometry, scvf)
223 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
224 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
225 insideVolVars.extrusionFactor());
226
227 // proceed depending on the interior BC types used
228 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
229
230 // neumann-coupling
231 if (iBcTypes.hasOnlyNeumann())
232 {
233 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
234 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
235 /facetVolVars.extrusionFactor()
236 *vtmv(scvf.unitOuterNormal(),
237 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
238 scvf.unitOuterNormal());
239
240 // The fluxes across this face and the outside face can be expressed in matrix form:
241 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
242 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
243 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
244 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
245 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
246 // that allow the description of the fluxes as functions of the cell and Dirichlet mass/mole fractions only.
247 if (!scvf.boundary())
248 {
249 const auto outsideScvIdx = scvf.outsideScvIdx();
250 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
251 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
252 *computeTpfaTransmissibility(fvGeometry, scvf, fvGeometry.scv(outsideScvIdx),
253 outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
254 outsideVolVars.extrusionFactor());
255
256 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
257 {
258 // optimized implementation: factorization obtained using sympy
259 // see CCTpfaFacetCouplingDarcysLaw for more details
260 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
261 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
262 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
263 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
264 }
265 else
266 {
267 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
268 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
269 tij[Cache::outsideTijIdx] = 0.0;
270 }
271 }
272 else
273 {
274 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
275 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
276 tij[Cache::outsideTijIdx] = 0.0;
277 }
278 }
279 else if (iBcTypes.hasOnlyDirichlet())
280 {
281 tij[Cache::insideTijIdx] = wIn;
282 tij[Cache::outsideTijIdx] = 0.0;
283 tij[Cache::facetTijIdx] = -wIn;
284 }
285 else
286 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
287
288 return tij;
289 }
290};
291
296template<class TypeTag, ReferenceSystemFormulation referenceSystem>
297class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, /*isNetwork*/true>
298: public FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
299{
302
304 using FVElementGeometry = typename GridGeometry::LocalView;
305 using SubControlVolume = typename GridGeometry::SubControlVolume;
306 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
307 using Extrusion = Extrusion_t<GridGeometry>;
308
309 using GridView = typename GridGeometry::GridView;
310 using Element = typename GridView::template Codim<0>::Entity;
311
315
316 static const int numPhases = ModelTraits::numFluidPhases();
317 static const int numComponents = ModelTraits::numFluidComponents();
318
320 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
321
325 class FacetCouplingFicksLawCache
326 {
327 public:
329 using Filler = typename ParentType::Cache::Filler;
330
334 static constexpr int insideTijIdx = 0;
335 static constexpr int facetTijIdx = 1;
336
338 using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
339
341 template< class Problem, class ElementVolumeVariables >
342 void updateDiffusion(const Problem& problem,
343 const Element& element,
344 const FVElementGeometry& fvGeometry,
345 const ElementVolumeVariables& elemVolVars,
346 const SubControlVolumeFace &scvf,
347 unsigned int phaseIdx,
348 unsigned int compIdx)
349 {
350 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
351 }
352
357 Scalar diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
358 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
359
361 Scalar diffusionTijInside(unsigned int phaseIdx, unsigned int compIdx) const
362 { return tij_[phaseIdx][compIdx][insideTijIdx]; }
363
365 Scalar diffusionTijFacet(unsigned int phaseIdx, unsigned int compIdx) const
366 { return tij_[phaseIdx][compIdx][facetTijIdx]; }
367
368 private:
369 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
370 };
371
372public:
374 using Cache = FacetCouplingFicksLawCache;
375
378 static constexpr DiscretizationMethod discMethod{};
379
382 { return referenceSystem; }
383
385 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
386 static ComponentFluxVector flux(const Problem& problem,
387 const Element& element,
388 const FVElementGeometry& fvGeometry,
389 const ElementVolumeVariables& elemVolVars,
390 const SubControlVolumeFace& scvf,
391 int phaseIdx,
392 const ElementFluxVarsCache& elemFluxVarsCache)
393 {
394 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
395 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
396
397 ComponentFluxVector componentFlux(0.0);
398 for (int compIdx = 0; compIdx < numComponents; compIdx++)
399 {
400 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
401 continue;
402
403 // get inside/outside volume variables
404 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
405 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
406 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
407
408 // the inside and outside mass/mole fractions fractions
409 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
410 const Scalar xFacet = massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
411
412 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
413 const Scalar rhoFacet = massOrMolarDensity(facetVolVars, referenceSystem, phaseIdx);
414 const Scalar rho = 0.5*(rhoInside + rhoFacet);
415
416 componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
417 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
418
419 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
420 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
421 }
422
423 return componentFlux;
424 }
425
426 // The flux variables cache has to be bound to an element prior to flux calculations
427 // During the binding, the transmissibility will be computed and stored using the method below.
428 template< class Problem, class ElementVolumeVariables >
429 static typename Cache::DiffusionTransmissibilityContainer
430 calculateTransmissibility(const Problem& problem,
431 const Element& element,
432 const FVElementGeometry& fvGeometry,
433 const ElementVolumeVariables& elemVolVars,
434 const SubControlVolumeFace& scvf,
435 unsigned int phaseIdx, unsigned int compIdx)
436 {
437 typename Cache::DiffusionTransmissibilityContainer tij;
438 if (!problem.couplingManager().isCoupled(element, scvf))
439 {
441 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
442 return tij;
443 }
444
446 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
447
448 // On surface grids only xi = 1.0 can be used, as the coupling condition
449 // for xi != 1.0 does not generalize for surface grids where the normal
450 // vectors of the inside/outside elements have different orientations.
451 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
452 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
453
454 const auto insideScvIdx = scvf.insideScvIdx();
455 const auto& insideScv = fvGeometry.scv(insideScvIdx);
456 const auto& insideVolVars = elemVolVars[insideScvIdx];
457 const auto wIn = Extrusion::area(fvGeometry, scvf)
458 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
459 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
460 insideVolVars.extrusionFactor());
461
462 // proceed depending on the interior BC types used
463 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
464
465 // neumann-coupling
466 if (iBcTypes.hasOnlyNeumann())
467 {
468 // Here we use the square root of the facet extrusion factor
469 // as an approximate average distance from scvf ip to facet center
470 using std::sqrt;
471 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
472 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
473 /sqrt(facetVolVars.extrusionFactor())
474 *vtmv(scvf.unitOuterNormal(),
475 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
476 scvf.unitOuterNormal());
477
478 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
479 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
480 }
481 else if (iBcTypes.hasOnlyDirichlet())
482 {
483 tij[Cache::insideTijIdx] = wIn;
484 tij[Cache::facetTijIdx] = -wIn;
485 }
486 else
487 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
488
489 return tij;
490 }
491};
492
493} // end namespace Dumux
494
495#endif
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:59
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:151
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:201
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:146
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:299
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:381
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:430
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:386
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:38
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:39
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:118
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.
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
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
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:851
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
Define some often used mathematical functions.
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...