Linear capillary pressure and relative permeability <-> saturation relations.
More...
#include <dumux/material/fluidmatrixinteractions/2p/linearmaterial.hh>
template<class ScalarT, class ParamsT = LinearMaterialParams<ScalarT>>
class Dumux::LinearMaterial< ScalarT, ParamsT >
Linear capillary pressure and relative permeability <-> saturation relations.
The entry pressure is reached at \(\mathrm{ \overline{S}_w = 1}\), the maximum capillary pressure is observed at \(\mathrm{ \overline{S}_w = 0}\).
For general info: EffToAbsLaw
- See also
- LinearMaterialParams
|
using | Params = ParamsT |
|
using | Scalar = typename Params::Scalar |
|
◆ Params
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
◆ Scalar
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
◆ dpc_dswe()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
Returns the partial derivative of the capillary pressure w.r.t. the effective saturation.
This is equivalent to \(\mathrm{ \frac{\partial p_C}{\partial \overline{S}_w} = - (p_{C,max} - p_{C,min}) }\)
- Parameters
-
swe | Effective saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\) conversion from absolute saturation happened in EffToAbsLaw. |
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
- Returns
- Partial derivative of \(\mathrm{[p_c]}\) w.r.t. effective saturation according to linear material relation.
◆ dsw_dpc()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
Returns the partial derivative of the effective saturation to the capillary pressure.
- Parameters
-
pc | Capillary pressure \(\mathrm{[p_C]}\) in \(\mathrm{[Pa]}\). |
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
- Returns
- Partial derivative of effective saturation w.r.t. \(\mathrm{[p_c]}\) according to linear relation.
◆ endPointPc()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
- Parameters
-
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
◆ krn()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
The relative permeability for the nonwetting phase.
- Parameters
-
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
swe | Effective saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\) conversion from absolute saturation happened in EffToAbsLaw. |
- Returns
- Relative permeability of the nonwetting phase calculated as linear relation.
◆ krw()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
The relative permeability for the wetting phase.
- Parameters
-
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
swe | Effective saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\) conversion from absolute saturation happened in EffToAbsLaw. |
- Returns
- Relative permeability of the wetting phase calculated as linear relation.
◆ pc()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
The linear capillary pressure-saturation curve.
This material law is linear: \(\mathrm{ p_C = (1 - \overline{S}_w ) (p_{C,max} - p_{C,entry}) + p_{C,entry} }\)
- Parameters
-
swe | Effective saturation of the wetting phase \(\overline{S}_w\) conversion from absolute saturation happened in EffToAbsLaw. |
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
- Returns
- Capillary pressure calculated by linear constitutive relation.
◆ sw()
template<class ScalarT , class ParamsT = LinearMaterialParams<ScalarT>>
The saturation-capillary pressure curve.
This is the inverse of the capillary pressure-saturation curve: \(\mathrm{ S_w = 1 - \frac{p_C - p_{C,entry}}{p_{C,max} - p_{C,entry}} }\)
- Parameters
-
pc | Capillary pressure \(\mathrm{[p_C]}\) in \(\mathrm{[Pa]}\). |
params | A container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too. |
- Returns
- Effective wetting phase saturation calculated as inverse of the linear constitutive relation.
The documentation for this class was generated from the following file: