Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
More...
template<class Scalar, class GridGeometry, class FluxVariables>
class Dumux::ForchheimerVelocity< Scalar, GridGeometry, FluxVariables >
This file contains the calculation of the Forchheimer velocity for a given Darcy velocity
|
template<class ElementVolumeVariables > |
static DimWorldVector | velocity (const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity) |
|
template<class Problem , class ElementVolumeVariables , bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value, std::enable_if_t< scalarPerm, int > = 0> |
static DimWorldMatrix | calculateHarmonicMeanSqrtPermeability (const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf) |
| Returns the harmonic mean of \(\sqrt{K_0}\) and \(\sqrt{K_1}\). More...
|
|
template<class Problem , class ElementVolumeVariables , bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value, std::enable_if_t<!scalarPerm, int > = 0> |
static DimWorldMatrix | calculateHarmonicMeanSqrtPermeability (const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf) |
| Returns the harmonic mean of \(\sqrt{\mathbf{K_0}}\) and \(\sqrt{\mathbf{K_1}}\). More...
|
|
template<class Scalar , class GridGeometry , class FluxVariables >
template<class Problem , class ElementVolumeVariables , bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value, std::enable_if_t< scalarPerm, int > = 0>
static DimWorldMatrix Dumux::ForchheimerVelocity< Scalar, GridGeometry, FluxVariables >::calculateHarmonicMeanSqrtPermeability |
( |
const Problem & |
problem, |
|
|
const ElementVolumeVariables & |
elemVolVars, |
|
|
const SubControlVolumeFace & |
scvf |
|
) |
| |
|
inlinestatic |
This is a specialization for scalar-valued permeabilities which returns a tensor with identical diagonal entries.
template<class Scalar , class GridGeometry , class FluxVariables >
template<class Problem , class ElementVolumeVariables , bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value, std::enable_if_t<!scalarPerm, int > = 0>
static DimWorldMatrix Dumux::ForchheimerVelocity< Scalar, GridGeometry, FluxVariables >::calculateHarmonicMeanSqrtPermeability |
( |
const Problem & |
problem, |
|
|
const ElementVolumeVariables & |
elemVolVars, |
|
|
const SubControlVolumeFace & |
scvf |
|
) |
| |
|
inlinestatic |
This is a specialization for tensor-valued permeabilities.