Buildings.Media.PerfectGases.MoistAir

Information


This is a medium model that is similar to Modelica.Media.Air.MoistAir but it has a constant specific heat capacity.

In particular, the medium is thermally perfect, i.e.,

In addition, the gas is calorically perfect, i.e., the specific heat capacities at constant pressure and constant volume are both constant (Bower 1998).

References

Bower, William B. A primer in fluid mechanics: Dynamics of flows in one space dimension. CRC Press. 1998.

Package Content

NameDescription
Water=1Index of water (in substanceNames, massFractions X, etc.)
Air=2Index of air (in substanceNames, massFractions X, etc.)
k_mair=steam.MM/dryair.MMratio of molar weights
dryair=Buildings.Media.PerfectGases.Common.SingleGasData.Air 
steam=Buildings.Media.PerfectGases.Common.SingleGasData.H2O 
Buildings.Media.PerfectGases.MoistAir.BaseProperties BaseProperties  
Buildings.Media.PerfectGases.MoistAir.Xsaturation Xsaturation Steam water mass fraction of saturation boundary in kg_water/kg_moistair
Buildings.Media.PerfectGases.MoistAir.setState_pTX setState_pTX Thermodynamic state as function of p, T and composition X
Buildings.Media.PerfectGases.MoistAir.setState_phX setState_phX Thermodynamic state as function of p, h and composition X
Buildings.Media.PerfectGases.MoistAir.setState_dTX setState_dTX Thermodynamic state as function of d, T and composition X
Buildings.Media.PerfectGases.MoistAir.gasConstant gasConstant Gas constant (computation neglects liquid fraction)
Buildings.Media.PerfectGases.MoistAir.saturationPressureLiquid saturationPressureLiquid Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned
Buildings.Media.PerfectGases.MoistAir.sublimationPressureIce sublimationPressureIce Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned
Buildings.Media.PerfectGases.MoistAir.saturationPressure saturationPressure Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)
Buildings.Media.PerfectGases.MoistAir.pressure pressure Gas pressure
Buildings.Media.PerfectGases.MoistAir.temperature temperature Gas temperature
Buildings.Media.PerfectGases.MoistAir.density density Gas density
Buildings.Media.PerfectGases.MoistAir.specificEntropy specificEntropy Specific entropy (liquid part neglected, mixing entropy included)
Buildings.Media.PerfectGases.MoistAir.enthalpyOfVaporization enthalpyOfVaporization Enthalpy of vaporization of water
Buildings.Media.PerfectGases.MoistAir.HeatCapacityOfWater HeatCapacityOfWater Specific heat capacity of water (liquid only) which is constant
Buildings.Media.PerfectGases.MoistAir.enthalpyOfLiquid enthalpyOfLiquid Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature
Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfLiquid der_enthalpyOfLiquid Temperature derivative of enthalpy of liquid per unit mass of liquid
Buildings.Media.PerfectGases.MoistAir.enthalpyOfCondensingGas enthalpyOfCondensingGas Enthalpy of steam per unit mass of steam
Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfCondensingGas der_enthalpyOfCondensingGas Derivative of enthalpy of steam per unit mass of steam
Buildings.Media.PerfectGases.MoistAir.enthalpyOfGas enthalpyOfGas Enthalpy of gas mixture per unit mass of gas mixture
Buildings.Media.PerfectGases.MoistAir.enthalpyOfDryAir enthalpyOfDryAir Enthalpy of dry air per unit mass of dry air
Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfDryAir der_enthalpyOfDryAir Derivative of enthalpy of dry air per unit mass of dry air
Buildings.Media.PerfectGases.MoistAir.specificHeatCapacityCp specificHeatCapacityCp Specific heat capacity of gas mixture at constant pressure
Buildings.Media.PerfectGases.MoistAir.specificHeatCapacityCv specificHeatCapacityCv Specific heat capacity of gas mixture at constant volume
Buildings.Media.PerfectGases.MoistAir.dynamicViscosity dynamicViscosity dynamic viscosity of dry air
Buildings.Media.PerfectGases.MoistAir.thermalConductivity thermalConductivity Thermal conductivity of dry air as a polynomial in the temperature
Buildings.Media.PerfectGases.MoistAir.h_pTX h_pTX Compute specific enthalpy from pressure, temperature and mass fraction
Buildings.Media.PerfectGases.MoistAir.specificEnthalpy specificEnthalpy Specific enthalpy
Buildings.Media.PerfectGases.MoistAir.specificInternalEnergy specificInternalEnergy Specific internal energy
Buildings.Media.PerfectGases.MoistAir.specificGibbsEnergy specificGibbsEnergy Specific Gibbs energy
Buildings.Media.PerfectGases.MoistAir.specificHelmholtzEnergy specificHelmholtzEnergy Specific Helmholtz energy
Buildings.Media.PerfectGases.MoistAir.T_phX T_phX Compute temperature from specific enthalpy and mass fraction

Types and constants

  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;


Buildings.Media.PerfectGases.MoistAir.BaseProperties Buildings.Media.PerfectGases.MoistAir.BaseProperties

Parameters

TypeNameDefaultDescription
Advanced
BooleanpreferredMediumStatesfalse= true if StateSelect.prefer shall be used for the independent property variables of the medium

Modelica definition

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);
  //                
  u = h - R*T;
  d = p/(R*T);
  /* 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;

Buildings.Media.PerfectGases.MoistAir.Xsaturation

Steam water mass fraction of saturation boundary in kg_water/kg_moistair

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate shermodynamic state

Outputs

TypeNameDescription
MassFractionX_satsteam mass fraction of sat. boundary [kg/kg]

Modelica definition

function Xsaturation = Modelica.Media.Air.MoistAir.Xsaturation 
  "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";

Buildings.Media.PerfectGases.MoistAir.setState_pTX Buildings.Media.PerfectGases.MoistAir.setState_pTX

Thermodynamic state as function of p, T and composition X

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
TemperatureT Temperature [K]
MassFractionX[:] Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_pTX 
  "Thermodynamic state as function of p, T and composition X" 
    extends Modelica.Media.Air.MoistAir.setState_pTX;
end setState_pTX;

Buildings.Media.PerfectGases.MoistAir.setState_phX Buildings.Media.PerfectGases.MoistAir.setState_phX

Thermodynamic state as function of p, h and composition X

Information


Function to set the state for given pressure, enthalpy and species concentration.
This function needed to be reimplemented in order for the medium model to use
the implementation of T_phX provided by this package as opposed to the 
implementation provided by its parent package.

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEnthalpyh Specific enthalpy [J/kg]
MassFractionX[:] Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

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)}));
end setState_phX;

Buildings.Media.PerfectGases.MoistAir.setState_dTX Buildings.Media.PerfectGases.MoistAir.setState_dTX

Thermodynamic state as function of d, T and composition X

Inputs

TypeNameDefaultDescription
Densityd density [kg/m3]
TemperatureT Temperature [K]
MassFractionX[:] Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_dTX 
  "Thermodynamic state as function of d, T and composition X" 
   extends Modelica.Media.Air.MoistAir.setState_dTX;
end setState_dTX;

Buildings.Media.PerfectGases.MoistAir.gasConstant Buildings.Media.PerfectGases.MoistAir.gasConstant

Gas constant (computation neglects liquid fraction)

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state

Outputs

TypeNameDescription
SpecificHeatCapacityRmixture gas constant [J/(kg.K)]

Modelica definition

redeclare function gasConstant 
  "Gas constant (computation neglects liquid fraction)" 
   extends Modelica.Media.Air.MoistAir.gasConstant;
end gasConstant;

Buildings.Media.PerfectGases.MoistAir.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

Inputs

TypeNameDefaultDescription
TemperatureTsat saturation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatsaturation pressure [Pa]

Modelica definition

function saturationPressureLiquid = 
    Modelica.Media.Air.MoistAir.saturationPressureLiquid 
  "Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned";

Buildings.Media.PerfectGases.MoistAir.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

Inputs

TypeNameDefaultDescription
TemperatureTsat sublimation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatsublimation pressure [Pa]

Modelica definition

function sublimationPressureIce = 
    Modelica.Media.Air.MoistAir.sublimationPressureIce 
  "Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned";

Buildings.Media.PerfectGases.MoistAir.saturationPressure Buildings.Media.PerfectGases.MoistAir.saturationPressure

Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)

Inputs

TypeNameDefaultDescription
TemperatureTsat saturation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatsaturation pressure [Pa]

Modelica definition

redeclare function extends saturationPressure 
  "Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)" 
  
algorithm 
  psat := Buildings.Utilities.Math.spliceFunction(saturationPressureLiquid(Tsat),sublimationPressureIce(Tsat),Tsat-273.16,1.0);
end saturationPressure;

Buildings.Media.PerfectGases.MoistAir.pressure Buildings.Media.PerfectGases.MoistAir.pressure

Gas pressure

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function pressure "Gas pressure" 
   extends Modelica.Media.Air.MoistAir.pressure;
end pressure;

Buildings.Media.PerfectGases.MoistAir.temperature Buildings.Media.PerfectGases.MoistAir.temperature

Gas temperature

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function temperature "Gas temperature" 
   extends Modelica.Media.Air.MoistAir.temperature;
end temperature;

Buildings.Media.PerfectGases.MoistAir.density Buildings.Media.PerfectGases.MoistAir.density

Gas density

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function density "Gas density" 
   extends Modelica.Media.Air.MoistAir.density;
end density;

Buildings.Media.PerfectGases.MoistAir.specificEntropy Buildings.Media.PerfectGases.MoistAir.specificEntropy

Specific entropy (liquid part neglected, mixing entropy included)

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEntropysSpecific entropy [J/(kg.K)]

Modelica definition

redeclare function specificEntropy 
  "Specific entropy (liquid part neglected, mixing entropy included)" 
   extends Modelica.Media.Air.MoistAir.specificEntropy;
end specificEntropy;

Buildings.Media.PerfectGases.MoistAir.enthalpyOfVaporization Buildings.Media.PerfectGases.MoistAir.enthalpyOfVaporization

Enthalpy of vaporization of water

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyr0vaporization enthalpy [J/kg]

Modelica definition

redeclare function extends enthalpyOfVaporization 
  "Enthalpy of vaporization of water" 
algorithm 
 r0 := 2501014.5;
end enthalpyOfVaporization;

Buildings.Media.PerfectGases.MoistAir.HeatCapacityOfWater Buildings.Media.PerfectGases.MoistAir.HeatCapacityOfWater

Specific heat capacity of water (liquid only) which is constant

Inputs

TypeNameDefaultDescription
TemperatureT [K]

Outputs

TypeNameDescription
SpecificHeatCapacitycp_fl[J/(kg.K)]

Modelica definition

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.PerfectGases.MoistAir.enthalpyOfLiquid Buildings.Media.PerfectGases.MoistAir.enthalpyOfLiquid

Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhliquid enthalpy [J/kg]

Modelica definition

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.PerfectGases.MoistAir.der_enthalpyOfLiquid Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfLiquid

Temperature derivative of enthalpy of liquid per unit mass of liquid

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]
Temperatureder_T temperature derivative [K]

Outputs

TypeNameDescription
SpecificHeatCapacityder_hderivative of liquid enthalpy [J/(kg.K)]

Modelica definition

replaceable function der_enthalpyOfLiquid 
  "Temperature derivative of enthalpy of liquid per unit mass of liquid" 
  extends Modelica.Icons.Function;
  input Temperature T "temperature";
  input Temperature der_T "temperature derivative";
  output SpecificHeatCapacity der_h "derivative of liquid enthalpy";
algorithm 
  der_h := 4186;
end der_enthalpyOfLiquid;

Buildings.Media.PerfectGases.MoistAir.enthalpyOfCondensingGas Buildings.Media.PerfectGases.MoistAir.enthalpyOfCondensingGas

Enthalpy of steam per unit mass of steam

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhsteam enthalpy [J/kg]

Modelica definition

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.PerfectGases.MoistAir.der_enthalpyOfCondensingGas Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfCondensingGas

Derivative of enthalpy of steam per unit mass of steam

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]
Temperatureder_T temperature derivative [K]

Outputs

TypeNameDescription
SpecificHeatCapacityder_hderivative of steam enthalpy [J/(kg.K)]

Modelica definition

replaceable function der_enthalpyOfCondensingGas 
  "Derivative of enthalpy of steam per unit mass of steam" 
  extends Modelica.Icons.Function;
  input Temperature T "temperature";
  input Temperature der_T "temperature derivative";
  output SpecificHeatCapacity der_h "derivative of steam enthalpy";
algorithm 
  der_h := steam.cp;
end der_enthalpyOfCondensingGas;

Buildings.Media.PerfectGases.MoistAir.enthalpyOfGas Buildings.Media.PerfectGases.MoistAir.enthalpyOfGas

Enthalpy of gas mixture per unit mass of gas mixture

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]
MassFractionX[:] vector of mass fractions [kg/kg]

Outputs

TypeNameDescription
SpecificEnthalpyhliquid enthalpy [J/kg]

Modelica definition

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.PerfectGases.MoistAir.enthalpyOfDryAir Buildings.Media.PerfectGases.MoistAir.enthalpyOfDryAir

Enthalpy of dry air per unit mass of dry air

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhdry air enthalpy [J/kg]

Modelica definition

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.PerfectGases.MoistAir.der_enthalpyOfDryAir Buildings.Media.PerfectGases.MoistAir.der_enthalpyOfDryAir

Derivative of enthalpy of dry air per unit mass of dry air

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]
Temperatureder_T temperature derivative [K]

Outputs

TypeNameDescription
SpecificHeatCapacityder_hderivative of dry air enthalpy [J/(kg.K)]

Modelica definition

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 Temperature der_T "temperature derivative";
  output SpecificHeatCapacity der_h "derivative of dry air enthalpy";
algorithm 
  der_h := dryair.cp;
end der_enthalpyOfDryAir;

Buildings.Media.PerfectGases.MoistAir.specificHeatCapacityCp Buildings.Media.PerfectGases.MoistAir.specificHeatCapacityCp

Specific heat capacity of gas mixture at constant pressure

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificHeatCapacitycpSpecific heat capacity at constant pressure [J/(kg.K)]

Modelica definition

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.PerfectGases.MoistAir.specificHeatCapacityCv Buildings.Media.PerfectGases.MoistAir.specificHeatCapacityCv

Specific heat capacity of gas mixture at constant volume

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificHeatCapacitycvSpecific heat capacity at constant volume [J/(kg.K)]

Modelica definition

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.PerfectGases.MoistAir.dynamicViscosity Buildings.Media.PerfectGases.MoistAir.dynamicViscosity

dynamic viscosity of dry air

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

redeclare function extends dynamicViscosity 
  "dynamic viscosity of dry air" 
algorithm 
  eta := 1.85E-5;
end dynamicViscosity;

Buildings.Media.PerfectGases.MoistAir.thermalConductivity Buildings.Media.PerfectGases.MoistAir.thermalConductivity

Thermal conductivity of dry air as a polynomial in the temperature

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
ThermalConductivitylambdaThermal conductivity [W/(m.K)]

Modelica definition

redeclare function extends thermalConductivity 
  "Thermal conductivity of dry air as a polynomial in the temperature" 
algorithm 
  lambda := Polynomials_Temp.evaluate({(-4.8737307422969E-008), 7.67803133753502E-005, 0.0241814385504202},
   Cv.to_degC(state.T));
end thermalConductivity;

Buildings.Media.PerfectGases.MoistAir.h_pTX Buildings.Media.PerfectGases.MoistAir.h_pTX

Compute specific enthalpy from pressure, temperature and mass fraction

Inputs

TypeNameDefaultDescription
Pressurep Pressure [Pa]
TemperatureT Temperature [K]
MassFractionX[nX] Mass fractions of moist air [1]

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy at p, T, X [J/kg]

Modelica definition

function h_pTX 
  "Compute specific enthalpy from pressure, temperature and mass fraction" 
  extends Modelica.Icons.Function;
  input SI.Pressure p "Pressure";
  input SI.Temperature T "Temperature";
  input SI.MassFraction X[nX] "Mass fractions of moist air";
  output SI.SpecificEnthalpy h "Specific enthalpy at p, T, X";
  
protected 
  SI.AbsolutePressure p_steam_sat "Partial saturation pressure of steam";
  SI.MassFraction x_sat "steam water mass fraction of saturation boundary";
  SI.MassFraction X_liquid "mass fraction of liquid water";
  SI.MassFraction X_steam "mass fraction of steam water";
  SI.MassFraction X_air "mass fraction of air";
  SI.SpecificEnthalpy hDryAir "Enthalpy of dry air";
algorithm 
  p_steam_sat :=saturationPressure(T);
  x_sat    :=k_mair*p_steam_sat/(p - p_steam_sat);
  X_liquid :=max(X[Water] - x_sat/(1 + x_sat), 0.0);
  X_steam  :=X[Water] - X_liquid;
  X_air    :=1 - X[Water];
  
/* THIS DOES NOT WORK --------------------------    
  h := enthalpyOfDryAir(T) * X_air + 
       Modelica.Media.Air.MoistAir.enthalpyOfCondensingGas(T) * X_steam + enthalpyOfLiquid(T)*X_liquid;
--------------------------------- */
  
/* THIS WORKS!!!! +++++++++++++++++++++
  h := (T - 273.15)*dryair.cp * X_air + 
       Modelica.Media.Air.MoistAir.enthalpyOfCondensingGas(T) * X_steam + enthalpyOfLiquid(T)*X_liquid;
 +++++++++++++++++++++*/
  
  hDryAir := (T - 273.15)*dryair.cp;
  h := hDryAir * X_air +
       ((T-273.15) * steam.cp + 2501014.5) * X_steam +
       (T - 273.15)*4186*X_liquid;
end h_pTX;

Buildings.Media.PerfectGases.MoistAir.specificEnthalpy Buildings.Media.PerfectGases.MoistAir.specificEnthalpy

Specific enthalpy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

redeclare function extends specificEnthalpy "Specific enthalpy" 
algorithm 
  h := h_pTX(state.p, state.T, state.X);
end specificEnthalpy;

Buildings.Media.PerfectGases.MoistAir.specificInternalEnergy Buildings.Media.PerfectGases.MoistAir.specificInternalEnergy

Specific internal energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergyuSpecific internal energy [J/kg]

Modelica definition

redeclare function extends specificInternalEnergy 
  "Specific internal energy" 
  extends Modelica.Icons.Function;
algorithm 
  u := h_pTX(state.p,state.T,state.X) - gasConstant(state)*state.T;
end specificInternalEnergy;

Buildings.Media.PerfectGases.MoistAir.specificGibbsEnergy Buildings.Media.PerfectGases.MoistAir.specificGibbsEnergy

Specific Gibbs energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

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.PerfectGases.MoistAir.specificHelmholtzEnergy Buildings.Media.PerfectGases.MoistAir.specificHelmholtzEnergy

Specific Helmholtz energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

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;

Buildings.Media.PerfectGases.MoistAir.T_phX

Compute temperature from specific enthalpy and mass fraction

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEnthalpyh specific enthalpy [J/kg]
MassFractionX[:] mass fractions of composition [kg/kg]

Outputs

TypeNameDescription
TemperatureTtemperature [K]

Modelica definition

function T_phX 
  "Compute temperature from specific enthalpy and mass fraction" 
  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "specific enthalpy";
  input MassFraction[:] X "mass fractions of composition";
  output Temperature T "temperature";
protected 
package Internal 
    "Solve h(data,T) for T with given h (use only indirectly via temperature_phX)" 
  extends Modelica.Media.Common.OneNonLinearEquation;
  redeclare record extends f_nonlinear_Data 
      "Data to be passed to non-linear function" 
    extends Buildings.Media.PerfectGases.Common.DataRecord;
  end f_nonlinear_Data;
    
  redeclare function extends f_nonlinear 
  algorithm 
      y := h_pTX(p,x,X);
  end f_nonlinear;
    
  // Dummy definition has to be added for current Dymola
  redeclare function extends solve 
  end solve;
end Internal;
  
algorithm 
 /* The function call below has been changed from 
      Internal.solve(h, 200, 6000, p, X[1:nXi], steam);
    to  
      Internal.solve(h, 200, 6000, p, X, steam);
    The reason is that when running the problem
       Buildings.Media.PerfectGases.Examples.MoistAirComparison
    then an assertion is triggered because the vector X had the wrong
    dimension. The above example verifies that T(h(T)) = T.
 */
  T := Internal.solve(h, 200, 6000, p, X, steam);
end T_phX;

HTML-documentation generated by Dymola Tue Sep 30 14:24:47 2008.