Hagen–Poiseuille-type flux law to describe the advective flux for pore-network models.
More...
#include <dumux/flux/porenetwork/advection.hh>
template<class ScalarT, class... TransmissibilityLawTypes>
class Dumux::PoreNetwork::CreepingFlow< ScalarT, TransmissibilityLawTypes >
Hagen–Poiseuille-type flux law to describe the advective flux for pore-network models.
|
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class SubControlVolumeFace , class ElemFluxVarsCache > |
static Scalar | flux (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElemFluxVarsCache &elemFluxVarsCache) |
| Returns the advective flux of a fluid phase across the given sub-control volume face (corresponding to a pore throat). More...
|
|
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class FluxVariablesCache > |
static Scalar | calculateTransmissibility (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx) |
| Returns the throat conductivity. More...
|
|
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class FluxVariablesCache > |
static std::array< Scalar, 2 > | calculateTransmissibilities (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache) |
|
◆ Scalar
template<class ScalarT , class... TransmissibilityLawTypes>
◆ Transmissibility
template<class ScalarT , class... TransmissibilityLawTypes>
Export the transmissibility law.
◆ calculateTransmissibilities()
template<class ScalarT , class... TransmissibilityLawTypes>
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class FluxVariablesCache >
static std::array< Scalar, 2 > Dumux::PoreNetwork::CreepingFlow< ScalarT, TransmissibilityLawTypes >::calculateTransmissibilities |
( |
const Problem & |
problem, |
|
|
const Element & |
element, |
|
|
const FVElementGeometry & |
fvGeometry, |
|
|
const ElementVolumeVariables & |
elemVolVars, |
|
|
const typename FVElementGeometry::SubControlVolumeFace & |
scvf, |
|
|
const FluxVariablesCache & |
fluxVarsCache |
|
) |
| |
|
inlinestatic |
◆ calculateTransmissibility()
template<class ScalarT , class... TransmissibilityLawTypes>
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class FluxVariablesCache >
static Scalar Dumux::PoreNetwork::CreepingFlow< ScalarT, TransmissibilityLawTypes >::calculateTransmissibility |
( |
const Problem & |
problem, |
|
|
const Element & |
element, |
|
|
const FVElementGeometry & |
fvGeometry, |
|
|
const typename FVElementGeometry::SubControlVolumeFace & |
scvf, |
|
|
const ElementVolumeVariables & |
elemVolVars, |
|
|
const FluxVariablesCache & |
fluxVarsCache, |
|
|
const int |
phaseIdx |
|
) |
| |
|
inlinestatic |
Returns the throat conductivity.
◆ flux()
template<class ScalarT , class... TransmissibilityLawTypes>
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class SubControlVolumeFace , class ElemFluxVarsCache >
static Scalar Dumux::PoreNetwork::CreepingFlow< ScalarT, TransmissibilityLawTypes >::flux |
( |
const Problem & |
problem, |
|
|
const Element & |
element, |
|
|
const FVElementGeometry & |
fvGeometry, |
|
|
const ElementVolumeVariables & |
elemVolVars, |
|
|
const SubControlVolumeFace & |
scvf, |
|
|
const int |
phaseIdx, |
|
|
const ElemFluxVarsCache & |
elemFluxVarsCache |
|
) |
| |
|
inlinestatic |
Returns the advective flux of a fluid phase across the given sub-control volume face (corresponding to a pore throat).
- Note
- The flux is given in N*m, and can be converted into a volume flux (m^3/s) or mass flux (kg/s) by applying an upwind scheme for the mobility (1/viscosity) or the product of density and mobility, respectively.
The documentation for this class was generated from the following file: