Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid

Package with moist air model that decouples pressure and temperature and that has no liquid water

Information


This is a medium model that is identical to Buildings.Media.GasesPTDecoupled.MoistAir, except that the water vapor is not condensing to the liquid phase if Xi[Water] > X_sat. Instead, all water is assumed in vapor phase.

This new formulation often leads to smaller systems of nonlinear equations because it allows to invert the function T_phX analytically.

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.GasesPTDecoupled.MoistAirNoLiquid.BaseProperties BaseProperties  
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.Xsaturation Xsaturation Steam water mass fraction of saturation boundary in kg_water/kg_moistair
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.setState_pTX setState_pTX Thermodynamic state as function of p, T and composition X
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.setState_phX setState_phX Thermodynamic state as function of p, h and composition X
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.setState_dTX setState_dTX Thermodynamic state as function of d, T and composition X
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.gasConstant gasConstant Gas constant (computation neglects liquid fraction)
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.saturationPressureLiquid saturationPressureLiquid Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.sublimationPressureIce sublimationPressureIce Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.saturationPressure saturationPressure Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.pressure pressure Gas pressure
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.temperature temperature Gas temperature
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.density density Gas density
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificEntropy specificEntropy Specific entropy (liquid part neglected, mixing entropy included)
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfVaporization enthalpyOfVaporization Enthalpy of vaporization of water
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.HeatCapacityOfWater HeatCapacityOfWater Specific heat capacity of water (liquid only) which is constant
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfLiquid enthalpyOfLiquid Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfLiquid der_enthalpyOfLiquid Temperature derivative of enthalpy of liquid per unit mass of liquid
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfCondensingGas enthalpyOfCondensingGas Enthalpy of steam per unit mass of steam
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfCondensingGas der_enthalpyOfCondensingGas Derivative of enthalpy of steam per unit mass of steam
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfGas enthalpyOfGas Enthalpy of gas mixture per unit mass of gas mixture
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfDryAir enthalpyOfDryAir Enthalpy of dry air per unit mass of dry air
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfDryAir der_enthalpyOfDryAir Derivative of enthalpy of dry air per unit mass of dry air
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificHeatCapacityCp specificHeatCapacityCp Specific heat capacity of gas mixture at constant pressure
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificHeatCapacityCv specificHeatCapacityCv Specific heat capacity of gas mixture at constant volume
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.dynamicViscosity dynamicViscosity dynamic viscosity of dry air
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.thermalConductivity thermalConductivity Thermal conductivity of dry air as a polynomial in the temperature
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificEnthalpy specificEnthalpy Specific enthalpy
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificInternalEnergy specificInternalEnergy Specific internal energy
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificGibbsEnergy specificGibbsEnergy Specific Gibbs energy
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificHelmholtzEnergy specificHelmholtzEnergy Specific Helmholtz energy
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.h_pTX h_pTX Compute specific enthalpy from pressure, temperature and mass fraction
Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.BaseProperties Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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";
 constant AbsolutePressure pStp = 101325 "Pressure for which dStp is defined";
 constant Density dStp = 1.2 "Fluid density at pressure pStp";
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 - Xi[Water]) + steam.R*Xi[Water];
  //                
  u = h - R*T;
  //    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;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 = Buildings.Media.PerfectGases.MoistAir.Xsaturation 
  "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.setState_pTX Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 Buildings.Media.PerfectGases.MoistAir.setState_pTX;
end setState_pTX;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.setState_phX Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.setState_dTX Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 Buildings.Media.PerfectGases.MoistAir.setState_dTX;
end setState_dTX;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.gasConstant Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 Buildings.Media.PerfectGases.MoistAir.gasConstant;
end gasConstant;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.saturationPressureLiquid Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 = 
    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.MoistAirNoLiquid.sublimationPressureIce Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 = 
    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.MoistAirNoLiquid.saturationPressure Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.pressure Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.pressure

Gas pressure

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function pressure "Gas pressure" 
   extends Buildings.Media.PerfectGases.MoistAir.pressure;
end pressure;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.temperature Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.temperature

Gas temperature

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function temperature "Gas temperature" 
   extends Buildings.Media.PerfectGases.MoistAir.temperature;
end temperature;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.density Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.density

Gas density

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function density "Gas density" 
  extends Modelica.Icons.Function;
  input ThermodynamicState state;
  output Density d "Density";
algorithm 
 d :=state.p*1.2/101325;
end density;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificEntropy Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 Buildings.Media.PerfectGases.MoistAir.specificEntropy;
end specificEntropy;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfVaporization Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.HeatCapacityOfWater Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfLiquid Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfLiquid Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 + Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfVaporization(T);
end enthalpyOfCondensingGas;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfGas Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 := Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfCondensingGas(T)*X[Water]
       + Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfDryAir(T)*(1.0-X[Water]);
end enthalpyOfGas;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.enthalpyOfDryAir Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.der_enthalpyOfDryAir Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.specificHeatCapacityCp Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.specificHeatCapacityCv Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.dynamicViscosity Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.thermalConductivity Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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.GasesPTDecoupled.MoistAirNoLiquid.specificEnthalpy Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificEnthalpy

Specific enthalpy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificInternalEnergy Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 := Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.h_pTX(state.p,state.T,state.X) - gasConstant(state)*state.T;
end specificInternalEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificGibbsEnergy Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 := Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.h_pTX(state.p,state.T,state.X) - state.T*specificEntropy(state);
end specificGibbsEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificHelmholtzEnergy Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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 := Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.h_pTX(state.p,state.T,state.X)
         - gasConstant(state)*state.T
         - state.T*Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.specificEntropy(state);
end specificHelmholtzEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.h_pTX Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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";
algorithm 
  h := (T - 273.15)*dryair.cp * (1 - X[Water]) + ((T-273.15) * steam.cp + 2501014.5) * X[Water];
end h_pTX;

Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.T_phX Buildings.Media.GasesPTDecoupled.MoistAirNoLiquid.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" 
  extends Modelica.Icons.Function;
  
  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "specific enthalpy";
  input MassFraction[:] X "mass fractions of composition";
  output Temperature T "temperature";
algorithm 
  T := 273.15 + (h-2501014.5 * X[Water])/(dryair.cp * (1 - X[Water])+steam.cp*X[Water]);
end T_phX;

HTML-documentation generated by Dymola Fri Oct 31 16:24:04 2008.