This is a medium model that is identical to 
Buildings.Media.PerfectGases.MoistAir, except the 
equation d = p/(R*T) has been replaced with 
d/dStp = p/pStp where 
pStd and dStp are constants for a reference
temperature and density.
This new formulation often leads to smaller systems of nonlinear equations because pressure and temperature are decoupled, at the expense of accuracy.
Extends from Modelica.Media.Interfaces.PartialCondensingGases (Base class for mixtures of condensing and non-condensing gases).
| Name | Description | 
|---|---|
| Water=1 | Index of water (in substanceNames, massFractions X, etc.) | 
| Air=2 | Index of air (in substanceNames, massFractions X, etc.) | 
| k_mair=steam.MM/dryair.MM | ratio of molar weights | 
| dryair=Buildings.Media.PerfectGases.Common.SingleGasData.Air | |
| steam=Buildings.Media.PerfectGases.Common.SingleGasData.H2O | |
| pStp=101325 | Pressure for which dStp is defined | 
| dStp=1.2 | Fluid density at pressure pStp | 
|  BaseProperties | |
|  Xsaturation | Steam water mass fraction of saturation boundary in kg_water/kg_moistair | 
|  setState_pTX | Thermodynamic state as function of p, T and composition X | 
|  setState_phX | Thermodynamic state as function of p, h and composition X | 
|  setState_dTX | Thermodynamic state as function of d, T and composition X | 
|  gasConstant | Gas constant (computation neglects liquid fraction) | 
|  saturationPressureLiquid | Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned | 
|  sublimationPressureIce | Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned | 
|  saturationPressure | Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy) | 
|  pressure | Gas pressure | 
|  temperature | Gas temperature | 
|  density | Gas density | 
|  specificEntropy | Specific entropy (liquid part neglected, mixing entropy included) | 
|  enthalpyOfVaporization | Enthalpy of vaporization of water | 
|  HeatCapacityOfWater | Specific heat capacity of water (liquid only) which is constant | 
|  enthalpyOfLiquid | Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature | 
|  der_enthalpyOfLiquid | Temperature derivative of enthalpy of liquid per unit mass of liquid | 
|  enthalpyOfCondensingGas | Enthalpy of steam per unit mass of steam | 
|  der_enthalpyOfCondensingGas | Derivative of enthalpy of steam per unit mass of steam | 
|  enthalpyOfNonCondensingGas | Enthalpy of non-condensing gas per unit mass of steam | 
|  der_enthalpyOfNonCondensingGas | Derivative of enthalpy of non-condensing gas per unit mass of steam | 
|  enthalpyOfGas | Enthalpy of gas mixture per unit mass of gas mixture | 
|  enthalpyOfDryAir | Enthalpy of dry air per unit mass of dry air | 
|  der_enthalpyOfDryAir | Derivative of enthalpy of dry air per unit mass of dry air | 
|  specificHeatCapacityCp | Specific heat capacity of gas mixture at constant pressure | 
|  specificHeatCapacityCv | Specific heat capacity of gas mixture at constant volume | 
|  dynamicViscosity | dynamic viscosity of dry air | 
|  thermalConductivity | Thermal conductivity of dry air as a polynomial in the temperature | 
|  h_pTX | Compute specific enthalpy from pressure, temperature and mass fraction | 
|  specificEnthalpy | Specific enthalpy | 
|  specificInternalEnergy | Specific internal energy | 
|  specificGibbsEnergy | Specific Gibbs energy | 
|  specificHelmholtzEnergy | Specific Helmholtz energy | 
|  T_phX | Compute temperature from specific enthalpy and mass fraction | 
| Inherited | |
|  ThermodynamicState | thermodynamic state variables | 
|  FluidConstants | extended fluid constants | 
| fluidConstants | constant data for the fluid | 
|  moleToMassFractions | Return mass fractions X from mole fractions | 
|  massToMoleFractions | Return mole fractions from mass fractions X | 
| ThermoStates | Enumeration type for independent variables | 
| mediumName="unusablePartialMedium" | Name of the medium | 
| substanceNames={mediumName} | Names of the mixture substances. Set substanceNames={mediumName} if only one substance. | 
| extraPropertiesNames=fill("", 0) | Names of the additional (extra) transported properties. Set extraPropertiesNames=fill("",0) if unused | 
| singleState | = true, if u and d are not a function of pressure | 
| reducedX=true | = true if medium contains the equation sum(X) = 1.0; set reducedX=true if only one substance (see docu for details) | 
| fixedX=false | = true if medium contains the equation X = reference_X | 
| reference_p=101325 | Reference pressure of Medium: default 1 atmosphere | 
| reference_T=298.15 | Reference temperature of Medium: default 25 deg Celsius | 
| reference_X=fill(1/nX, nX) | Default mass fractions of medium | 
| p_default=101325 | Default value for pressure of medium (for initialization) | 
| T_default=Modelica.SIunits.Conversions.from_degC(20) | Default value for temperature of medium (for initialization) | 
| h_default=specificEnthalpy_pTX(p_default, T_default, X_default) | Default value for specific enthalpy of medium (for initialization) | 
| X_default=reference_X | Default value for mass fractions of medium (for initialization) | 
| nS=size(substanceNames, 1) | Number of substances | 
| nX=nS | Number of mass fractions | 
| nXi=if fixedX then 0 else if reducedX then nS - 1 else nS | Number of structurally independent mass fractions (see docu for details) | 
| nC=size(extraPropertiesNames, 1) | Number of extra (outside of standard mass-balance) transported properties | 
| C_nominal=1.0e-6*ones(nC) | Default for the nominal values for the extra properties | 
|  setState_psX | Return thermodynamic state as function of p, s and composition X or Xi | 
|  setSmoothState | Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b | 
|  prandtlNumber | Return the Prandtl number | 
|  heatCapacity_cp | alias for deprecated name | 
|  heatCapacity_cv | alias for deprecated name | 
|  isentropicExponent | Return isentropic exponent | 
|  isentropicEnthalpy | Return isentropic enthalpy | 
|  velocityOfSound | Return velocity of sound | 
|  isobaricExpansionCoefficient | Return overall the isobaric expansion coefficient beta | 
|  beta | alias for isobaricExpansionCoefficient for user convenience | 
|  isothermalCompressibility | Return overall the isothermal compressibility factor | 
|  kappa | alias of isothermalCompressibility for user convenience | 
|  density_derp_h | Return density derivative w.r.t. pressure at const specific enthalpy | 
|  density_derh_p | Return density derivative w.r.t. specific enthalpy at constant pressure | 
|  density_derp_T | Return density derivative w.r.t. pressure at const temperature | 
|  density_derT_p | Return density derivative w.r.t. temperature at constant pressure | 
|  density_derX | Return density derivative w.r.t. mass fraction | 
|  molarMass | Return the molar mass of the medium | 
|  specificEnthalpy_pTX | Return specific enthalpy from p, T, and X or Xi | 
|  specificEntropy_pTX | Return specific enthalpy from p, T, and X or Xi | 
|  density_pTX | Return density from p, T, and X or Xi | 
|  temperature_phX | Return temperature from p, h, and X or Xi | 
|  density_phX | Return density from p, h, and X or Xi | 
|  temperature_psX | Return temperature from p,s, and X or Xi | 
|  density_psX | Return density from p, s, and X or Xi | 
|  specificEnthalpy_psX | Return specific enthalpy from p, s, and X or Xi | 
| AbsolutePressure | Type for absolute pressure with medium specific attributes | 
| Density | Type for density with medium specific attributes | 
| DynamicViscosity | Type for dynamic viscosity with medium specific attributes | 
| EnthalpyFlowRate | Type for enthalpy flow rate with medium specific attributes | 
| MassFlowRate | Type for mass flow rate with medium specific attributes | 
| MassFraction | Type for mass fraction with medium specific attributes | 
| MoleFraction | Type for mole fraction with medium specific attributes | 
| MolarMass | Type for molar mass with medium specific attributes | 
| MolarVolume | Type for molar volume with medium specific attributes | 
| IsentropicExponent | Type for isentropic exponent with medium specific attributes | 
| SpecificEnergy | Type for specific energy with medium specific attributes | 
| SpecificInternalEnergy | Type for specific internal energy with medium specific attributes | 
| SpecificEnthalpy | Type for specific enthalpy with medium specific attributes | 
| SpecificEntropy | Type for specific entropy with medium specific attributes | 
| SpecificHeatCapacity | Type for specific heat capacity with medium specific attributes | 
| SurfaceTension | Type for surface tension with medium specific attributes | 
| Temperature | Type for temperature with medium specific attributes | 
| ThermalConductivity | Type for thermal conductivity with medium specific attributes | 
| PrandtlNumber | Type for Prandtl number with medium specific attributes | 
| VelocityOfSound | Type for velocity of sound with medium specific attributes | 
| ExtraProperty | Type for unspecified, mass-specific property transported by flow | 
| CumulativeExtraProperty | Type for conserved integral of unspecified, mass specific property | 
| ExtraPropertyFlowRate | Type for flow rate of unspecified, mass-specific property | 
| IsobaricExpansionCoefficient | Type for isobaric expansion coefficient with medium specific attributes | 
| DipoleMoment | Type for dipole moment with medium specific attributes | 
| DerDensityByPressure | Type for partial derivative of density with resect to pressure with medium specific attributes | 
| DerDensityByEnthalpy | Type for partial derivative of density with resect to enthalpy with medium specific attributes | 
| DerEnthalpyByPressure | Type for partial derivative of enthalpy with resect to pressure with medium specific attributes | 
| DerDensityByTemperature | Type for partial derivative of density with resect to temperature with medium specific attributes | 
|  Choices | Types, constants to define menu choices | 
constant Integer Water=1 "Index of water (in substanceNames, massFractions X, etc.)";
constant Integer Air=2 "Index of air (in substanceNames, massFractions X, etc.)";
constant Real k_mair = steam.MM/dryair.MM "ratio of molar weights";
  constant Buildings.Media.PerfectGases.Common.DataRecord dryair=
        Buildings.Media.PerfectGases.Common.SingleGasData.Air;
  constant Buildings.Media.PerfectGases.Common.DataRecord steam=
        Buildings.Media.PerfectGases.Common.SingleGasData.H2O;
constant AbsolutePressure pStp = 101325 "Pressure for which dStp is defined";
constant Density dStp = 1.2 "Fluid density at pressure pStp";
 Buildings.Media.GasesPTDecoupled.MoistAir.BaseProperties
Buildings.Media.GasesPTDecoupled.MoistAir.BaseProperties
| Type | Name | Default | Description | 
|---|---|---|---|
| Boolean | standardOrderComponents | true | if true, and reducedX = true, the last element of X will be computed from the other ones | 
| Advanced | |||
| Boolean | preferredMediumStates | false | = true if StateSelect.prefer shall be used for the independent property variables of the medium | 
redeclare replaceable model extends BaseProperties(
  T(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default),
  p(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default),
  Xi(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default))
  /* p, T, X = X[Water] are used as preferred states, since only then all
     other quantities can be computed in a recursive sequence. 
     If other variables are selected as states, static state selection
     is no longer possible and non-linear algebraic equations occur.
      */
  MassFraction x_water "Mass of total water/mass of dry air";
  Real phi "Relative humidity";
protected 
  constant SI.MolarMass[2] MMX = {steam.MM,dryair.MM} 
    "Molar masses of components";
  MassFraction X_liquid "Mass fraction of liquid water";
  MassFraction X_steam "Mass fraction of steam water";
  MassFraction X_air "Mass fraction of air";
  MassFraction X_sat 
    "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";
  MassFraction x_sat 
    "Steam water mass content of saturation boundary in kg_water/kg_dryair";
  AbsolutePressure p_steam_sat "Partial saturation pressure of steam";
equation 
  assert(T >= 200.0 and T <= 423.15, "
Temperature T is not in the allowed range
200.0 K <= (T ="
             + String(T) + " K) <= 423.15 K
required from medium model \""   + mediumName + "\".");
  MM = 1/(Xi[Water]/MMX[Water]+(1.0-Xi[Water])/MMX[Air]);
  p_steam_sat = min(saturationPressure(T),0.999*p);
  X_sat = min(p_steam_sat * k_mair/max(100*Modelica.Constants.eps, p - p_steam_sat)*(1 - Xi[Water]), 1.0) 
    "Water content at saturation with respect to actual water content";
  X_liquid = max(Xi[Water] - X_sat, 0.0);
  X_steam  = Xi[Water]-X_liquid;
  X_air    = 1-Xi[Water];
  h = specificEnthalpy_pTX(p,T,Xi);
  R = dryair.R*(1 - X_steam/(1 - X_liquid)) + steam.R*X_steam/(1 - X_liquid);
  // Equation for ideal gas, from h=u+p*v and R*T=p*v, from which follows that  u = h-R*T.
  // u = h-R*T;
  // However, in this medium, the gas law is d/dStp=p/pStp, from which follows using h=u+pv that
  // u= h-p*v = h-p/d = h-pStp/dStp
  u = h-pStp/dStp;
  //    d = p/(R*T);
  d/dStp = p/pStp;
  /* Note, u and d are computed under the assumption that the volume of the liquid
         water is neglible with respect to the volume of air and of steam
      */
  state.p = p;
  state.T = T;
  state.X = X;
  // this x_steam is water load / dry air!!!!!!!!!!!
  x_sat    = k_mair*p_steam_sat/max(100*Modelica.Constants.eps,p - p_steam_sat);
  x_water = Xi[Water]/max(X_air,100*Modelica.Constants.eps);
  phi = p/p_steam_sat*Xi[Water]/(Xi[Water] + k_mair*X_air);
end BaseProperties;
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | Thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| MassFraction | X_sat | Steam mass fraction of sat. boundary [kg/kg] | 
function Xsaturation = Buildings.Media.PerfectGases.MoistAir.Xsaturation "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";
 Buildings.Media.GasesPTDecoupled.MoistAir.setState_pTX
Buildings.Media.GasesPTDecoupled.MoistAir.setState_pTX
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| Temperature | T | Temperature [K] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | Thermodynamic state | 
redeclare function setState_pTX 
  "Thermodynamic state as function of p, T and composition X"
    extends Buildings.Media.PerfectGases.MoistAir.setState_pTX;
end setState_pTX;
 Buildings.Media.GasesPTDecoupled.MoistAir.setState_phX
Buildings.Media.GasesPTDecoupled.MoistAir.setState_phX
T_phX provided by this package as opposed to the 
implementation provided by its parent package.
Extends from Modelica.Icons.Function (Icon for functions).
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
| MassFraction | X[:] | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | 
redeclare function setState_phX 
  "Thermodynamic state as function of p, h and composition X"
extends Modelica.Icons.Function;
input AbsolutePressure p "Pressure";
input SpecificEnthalpy h "Specific enthalpy";
input MassFraction X[:] "Mass fractions";
output ThermodynamicState state;
algorithm 
state := if size(X,1) == nX then 
       ThermodynamicState(p=p,T=T_phX(p,h,X),X=X) else 
      ThermodynamicState(p=p,T=T_phX(p,h,X), X=cat(1,X,{1-sum(X)}));
   //   ThermodynamicState(p=p,T=T_phX(p,h,cat(1,X,{1-sum(X)})), X=cat(1,X,{1-sum(X)}));
end setState_phX;
 
 Buildings.Media.GasesPTDecoupled.MoistAir.setState_dTX
Buildings.Media.GasesPTDecoupled.MoistAir.setState_dTX
| Type | Name | Default | Description | 
|---|---|---|---|
| Density | d | density [kg/m3] | |
| Temperature | T | Temperature [K] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | Thermodynamic state | 
redeclare function setState_dTX "Thermodynamic state as function of d, T and composition X" extends Buildings.Media.PerfectGases.MoistAir.setState_dTX; end setState_dTX;
 Buildings.Media.GasesPTDecoupled.MoistAir.gasConstant
Buildings.Media.GasesPTDecoupled.MoistAir.gasConstant
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | R | mixture gas constant [J/(kg.K)] | 
redeclare function gasConstant "Gas constant (computation neglects liquid fraction)" extends Buildings.Media.PerfectGases.MoistAir.gasConstant; end gasConstant;
 Buildings.Media.GasesPTDecoupled.MoistAir.saturationPressureLiquid
Buildings.Media.GasesPTDecoupled.MoistAir.saturationPressureLiquid
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | Tsat | saturation temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | psat | saturation pressure [Pa] | 
function saturationPressureLiquid =
    Buildings.Media.PerfectGases.MoistAir.saturationPressureLiquid 
  "Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned";
 Buildings.Media.GasesPTDecoupled.MoistAir.sublimationPressureIce
Buildings.Media.GasesPTDecoupled.MoistAir.sublimationPressureIce
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | Tsat | sublimation temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | psat | sublimation pressure [Pa] | 
function sublimationPressureIce =
    Buildings.Media.PerfectGases.MoistAir.sublimationPressureIce 
  "Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned";
 Buildings.Media.GasesPTDecoupled.MoistAir.saturationPressure
Buildings.Media.GasesPTDecoupled.MoistAir.saturationPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | Tsat | saturation temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | psat | saturation pressure [Pa] | 
redeclare function extends saturationPressure 
  "Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)"
algorithm 
  psat := Buildings.Utilities.Math.Functions.spliceFunction(
                                                  saturationPressureLiquid(Tsat),sublimationPressureIce(Tsat),Tsat-273.16,1.0);
end saturationPressure;
 
 Buildings.Media.GasesPTDecoupled.MoistAir.pressure
Buildings.Media.GasesPTDecoupled.MoistAir.pressure
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | 
redeclare function pressure "Gas pressure" extends Buildings.Media.PerfectGases.MoistAir.pressure; end pressure;
 Buildings.Media.GasesPTDecoupled.MoistAir.temperature
Buildings.Media.GasesPTDecoupled.MoistAir.temperature
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | Temperature [K] | 
redeclare function temperature "Gas temperature" extends Buildings.Media.PerfectGases.MoistAir.temperature; end temperature;
 Buildings.Media.GasesPTDecoupled.MoistAir.density
Buildings.Media.GasesPTDecoupled.MoistAir.density
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | 
| Type | Name | Description | 
|---|---|---|
| Density | d | Density [kg/m3] | 
redeclare function density "Gas density" extends Modelica.Icons.Function; input ThermodynamicState state; output Density d "Density"; algorithm d :=state.p*dStp/pStp; end density;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificEntropy
Buildings.Media.GasesPTDecoupled.MoistAir.specificEntropy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | 
redeclare function specificEntropy "Specific entropy (liquid part neglected, mixing entropy included)" extends Buildings.Media.PerfectGases.MoistAir.specificEntropy; end specificEntropy;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfVaporization
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfVaporization
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | r0 | vaporization enthalpy [J/kg] | 
redeclare function extends enthalpyOfVaporization "Enthalpy of vaporization of water" algorithm r0 := 2501014.5; end enthalpyOfVaporization;
 Buildings.Media.GasesPTDecoupled.MoistAir.HeatCapacityOfWater
Buildings.Media.GasesPTDecoupled.MoistAir.HeatCapacityOfWater
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | cp_fl | [J/(kg.K)] | 
function HeatCapacityOfWater "Specific heat capacity of water (liquid only) which is constant" extends Modelica.Icons.Function; input Temperature T; output SpecificHeatCapacity cp_fl; algorithm cp_fl := 4186; end HeatCapacityOfWater;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfLiquid
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfLiquid
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | liquid enthalpy [J/kg] | 
redeclare replaceable function extends enthalpyOfLiquid "Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature" annotation(derivative=der_enthalpyOfLiquid); algorithm h := (T - 273.15)*4186;end enthalpyOfLiquid; 
 Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfLiquid
Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfLiquid
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | |
| Real | der_T | temperature derivative | 
| Type | Name | Description | 
|---|---|---|
| Real | der_h | derivative of liquid enthalpy | 
replaceable function der_enthalpyOfLiquid "Temperature derivative of enthalpy of liquid per unit mass of liquid" extends Modelica.Icons.Function; input Temperature T "temperature"; input Real der_T "temperature derivative"; output Real der_h "derivative of liquid enthalpy"; algorithm der_h := 4186*der_T; end der_enthalpyOfLiquid;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfCondensingGas
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfCondensingGas
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | steam enthalpy [J/kg] | 
redeclare function enthalpyOfCondensingGas "Enthalpy of steam per unit mass of steam" annotation(derivative=der_enthalpyOfCondensingGas); extends Modelica.Icons.Function; input Temperature T "temperature"; output SpecificEnthalpy h "steam enthalpy"; algorithm h := (T-273.15) * steam.cp + enthalpyOfVaporization(T);end enthalpyOfCondensingGas; 
 Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfCondensingGas
Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfCondensingGas
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | |
| Real | der_T | temperature derivative | 
| Type | Name | Description | 
|---|---|---|
| Real | der_h | derivative of steam enthalpy | 
replaceable function der_enthalpyOfCondensingGas "Derivative of enthalpy of steam per unit mass of steam" extends Modelica.Icons.Function; input Temperature T "temperature"; input Real der_T "temperature derivative"; output Real der_h "derivative of steam enthalpy"; algorithm der_h := steam.cp*der_T; end der_enthalpyOfCondensingGas;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfNonCondensingGas
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfNonCondensingGas
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | enthalpy [J/kg] | 
redeclare function enthalpyOfNonCondensingGas "Enthalpy of non-condensing gas per unit mass of steam" annotation(derivative=der_enthalpyOfNonCondensingGas); extends Modelica.Icons.Function; input Temperature T "temperature"; output SpecificEnthalpy h "enthalpy"; algorithm h := enthalpyOfDryAir(T);end enthalpyOfNonCondensingGas; 
 Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfNonCondensingGas
Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfNonCondensingGas
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | |
| Real | der_T | temperature derivative | 
| Type | Name | Description | 
|---|---|---|
| Real | der_h | derivative of steam enthalpy | 
replaceable function der_enthalpyOfNonCondensingGas "Derivative of enthalpy of non-condensing gas per unit mass of steam" extends Modelica.Icons.Function; input Temperature T "temperature"; input Real der_T "temperature derivative"; output Real der_h "derivative of steam enthalpy"; algorithm der_h := der_enthalpyOfDryAir(T, der_T); end der_enthalpyOfNonCondensingGas;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfGas
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfGas
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | |
| MassFraction | X[:] | vector of mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | specific enthalpy [J/kg] | 
redeclare replaceable function extends enthalpyOfGas 
  "Enthalpy of gas mixture per unit mass of gas mixture"
algorithm 
  h := enthalpyOfCondensingGas(T)*X[Water]
       + enthalpyOfDryAir(T)*(1.0-X[Water]);
end enthalpyOfGas;
 Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfDryAir
Buildings.Media.GasesPTDecoupled.MoistAir.enthalpyOfDryAir
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | dry air enthalpy [J/kg] | 
replaceable function enthalpyOfDryAir "Enthalpy of dry air per unit mass of dry air" annotation(derivative=der_enthalpyOfDryAir); extends Modelica.Icons.Function; input Temperature T "temperature"; output SpecificEnthalpy h "dry air enthalpy"; algorithm h := (T - 273.15)*dryair.cp;end enthalpyOfDryAir; 
 Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfDryAir
Buildings.Media.GasesPTDecoupled.MoistAir.der_enthalpyOfDryAir
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | |
| Real | der_T | temperature derivative | 
| Type | Name | Description | 
|---|---|---|
| Real | der_h | derivative of dry air enthalpy | 
replaceable function der_enthalpyOfDryAir "Derivative of enthalpy of dry air per unit mass of dry air" extends Modelica.Icons.Function; input Temperature T "temperature"; input Real der_T "temperature derivative"; output Real der_h "derivative of dry air enthalpy"; algorithm der_h := dryair.cp*der_T; end der_enthalpyOfDryAir;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificHeatCapacityCp
Buildings.Media.GasesPTDecoupled.MoistAir.specificHeatCapacityCp
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | cp | Specific heat capacity at constant pressure [J/(kg.K)] | 
redeclare replaceable function extends specificHeatCapacityCp "Specific heat capacity of gas mixture at constant pressure" algorithm cp := dryair.cp*(1-state.X[Water]) +steam.cp*state.X[Water]; end specificHeatCapacityCp;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificHeatCapacityCv
Buildings.Media.GasesPTDecoupled.MoistAir.specificHeatCapacityCv
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | cv | Specific heat capacity at constant volume [J/(kg.K)] | 
redeclare replaceable function extends specificHeatCapacityCv "Specific heat capacity of gas mixture at constant volume" algorithm cv:= dryair.cv*(1-state.X[Water]) +steam.cv*state.X[Water]; end specificHeatCapacityCv;
 Buildings.Media.GasesPTDecoupled.MoistAir.dynamicViscosity
Buildings.Media.GasesPTDecoupled.MoistAir.dynamicViscosity
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| DynamicViscosity | eta | Dynamic viscosity [Pa.s] | 
redeclare function extends dynamicViscosity "dynamic viscosity of dry air" algorithm eta := 1.85E-5; end dynamicViscosity;
 Buildings.Media.GasesPTDecoupled.MoistAir.thermalConductivity
Buildings.Media.GasesPTDecoupled.MoistAir.thermalConductivity
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| ThermalConductivity | lambda | Thermal conductivity [W/(m.K)] | 
redeclare function extends thermalConductivity 
  "Thermal conductivity of dry air as a polynomial in the temperature"
  import Modelica.Media.Incompressible.TableBased.Polynomials_Temp;
algorithm 
  lambda := Polynomials_Temp.evaluate({(-4.8737307422969E-008), 7.67803133753502E-005, 0.0241814385504202},
  Modelica.SIunits.Conversions.to_degC(state.T));
end thermalConductivity;
 Buildings.Media.GasesPTDecoupled.MoistAir.h_pTX
Buildings.Media.GasesPTDecoupled.MoistAir.h_pTX
| Type | Name | Default | Description | 
|---|---|---|---|
| Pressure | p | Pressure [Pa] | |
| Temperature | T | Temperature [K] | |
| MassFraction | X[:] | Mass fractions of moist air [1] | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | Specific enthalpy at p, T, X [J/kg] | 
function h_pTX "Compute specific enthalpy from pressure, temperature and mass fraction" extends Buildings.Media.PerfectGases.MoistAir.h_pTX; end h_pTX;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificEnthalpy
Buildings.Media.GasesPTDecoupled.MoistAir.specificEnthalpy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | 
redeclare function extends specificEnthalpy "Specific enthalpy" algorithm h := h_pTX(state.p, state.T, state.X); end specificEnthalpy;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificInternalEnergy
Buildings.Media.GasesPTDecoupled.MoistAir.specificInternalEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | u | Specific internal energy [J/kg] | 
redeclare function extends specificInternalEnergy "Specific internal energy" extends Modelica.Icons.Function; algorithm u := h_pTX(state.p,state.T,state.X) - pStp/dStp; end specificInternalEnergy;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificGibbsEnergy
Buildings.Media.GasesPTDecoupled.MoistAir.specificGibbsEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | g | Specific Gibbs energy [J/kg] | 
redeclare function extends specificGibbsEnergy "Specific Gibbs energy" extends Modelica.Icons.Function; algorithm g := h_pTX(state.p,state.T,state.X) - state.T*specificEntropy(state); end specificGibbsEnergy;
 Buildings.Media.GasesPTDecoupled.MoistAir.specificHelmholtzEnergy
Buildings.Media.GasesPTDecoupled.MoistAir.specificHelmholtzEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | f | Specific Helmholtz energy [J/kg] | 
redeclare function extends specificHelmholtzEnergy "Specific Helmholtz energy" extends Modelica.Icons.Function; algorithm f := h_pTX(state.p,state.T,state.X) - gasConstant(state)*state.T - state.T*specificEntropy(state); end specificHelmholtzEnergy;
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
| MassFraction | X[:] | Mass fractions of composition [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | Temperature [K] | 
function T_phX "Compute temperature from specific enthalpy and mass fraction" extends Buildings.Media.PerfectGases.MoistAir.T_phX; end T_phX;