3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Classes | Public Types | Static Public Member Functions | List of all members
Dumux::PoreNetwork::NonCreepingFlow< ScalarT, TransmissibilityLawTypes > Class Template Reference

#include <dumux/flux/porenetwork/advection.hh>

Classes

class  Transmissibility
 Inherit transmissibility from creeping flow transmissibility law to cache non-creeping flow-related parameters. More...
 

Public Types

using Scalar = ScalarT
 Export the Scalar type. More...
 
using TransmissibilityCreepingFlow = Detail::Transmissibility< TransmissibilityLawTypes... >
 Export the creeping flow transmissibility law. More...
 

Static Public Member Functions

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)
 
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)
 

Member Typedef Documentation

◆ Scalar

template<class ScalarT , class... TransmissibilityLawTypes>
using Dumux::PoreNetwork::NonCreepingFlow< ScalarT, TransmissibilityLawTypes >::Scalar = ScalarT

Export the Scalar type.

◆ TransmissibilityCreepingFlow

template<class ScalarT , class... TransmissibilityLawTypes>
using Dumux::PoreNetwork::NonCreepingFlow< ScalarT, TransmissibilityLawTypes >::TransmissibilityCreepingFlow = Detail::Transmissibility<TransmissibilityLawTypes...>

Export the creeping flow transmissibility law.

Member Function Documentation

◆ calculateTransmissibilities()

template<class ScalarT , class... TransmissibilityLawTypes>
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class FluxVariablesCache >
static std::array< Scalar, 2 > Dumux::PoreNetwork::NonCreepingFlow< 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::NonCreepingFlow< 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

◆ flux()

template<class ScalarT , class... TransmissibilityLawTypes>
template<class Problem , class Element , class FVElementGeometry , class ElementVolumeVariables , class SubControlVolumeFace , class ElemFluxVarsCache >
static Scalar Dumux::PoreNetwork::NonCreepingFlow< 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.

El-Zehairy et al.(2019): Eq.(14): A0Q^2 + B0Q + C0 = 0 where Q is volumetric flowrate, A0 = [Kc + Ke] * rho /[2aij^2], Kc and Ke are contraction and expansion coefficient respectively, aij is cross sectional area of the throat, B0 = 1/K0 (K0 is trnasmissibility used in paper) and C0 = -deltaP In our implementation viscosity of fluid will be included using upwinding later. Thus, creepingFlowTransmissibility = mu * K0 and the volumetric flowrate calculated here q = mu * Q. Substitution of new variables into El-Zehairy et al.(2019), Eq.(14) gives: Aq^2 + Bq + C = 0 where A = A0 / mu^2, B = B0/mu = 1/creepingFlowTransmissibility and C = C0 attention: the q, volumetric flowrate, calculated here is always positive and its sign needs to be determined based on flow direction this approach is taken to prevent the term under the square root (discriminant) becoming negative

give the volume flowrate proper sign based on flow direction.


The documentation for this class was generated from the following file: