Buildings.Media.Air

Package with moist air model that decouples pressure and temperature

Information

This medium package models moist air using a gas law in which pressure and temperature are independent, which often leads to significantly faster and more robust computations. The specific heat capacities at constant pressure and at constant volume are constant. The air is assumed to be not saturated.

This medium uses the gas law

ρ/ρstp = p/pstp,

where pstd and ρstp are constant reference temperature and density, rathern than the ideal gas law

ρ = p ⁄(R T),

where R is the gas constant and T is the temperature.

This formulation often leads to smaller systems of nonlinear equations because equations for pressure and temperature are decoupled. Therefore, if air inside a control volume such as room air is heated, it does not increase its specific volume. Consequently, merely heating or cooling a control volume does not affect the air flow calculations in a duct network that may be connected to that volume. Note that multizone air exchange simulation in which buoyancy drives the air flow is still possible as the models in Buildings.Airflow.Multizone compute the mass density using the function Buildings.Utilities.Psychrometrics.Functions.density_pTX in which density is a function of temperature.

Note that models in this package implement the equation for the internal energy as

u = h - pstp ⁄ ρstp,

where u is the internal energy per unit mass, h is the enthalpy per unit mass, pstp is the static pressure and ρstp is the mass density at standard pressure and temperature. The reason for this implementation is that in general,

h = u + p v,

from which follows that

u = h - p v = h - p ⁄ ρ = h - pstp ⁄ ρstd,

because p ⁄ ρ = pstp ⁄ ρstp in this medium model.

The enthalpy is computed using the convention that h=0 if T=0 °C and no water vapor is present.

Extends from Modelica.Media.Interfaces.PartialCondensingGases (Base class for mixtures of condensing and non-condensing gases), Modelica.Icons.Package (Icon for standard packages).

Package Content

Name Description
Water=1 Index of water (in substanceNames, massFractions X, etc.)
Air=2 Index of air (in substanceNames, massFractions X, etc.)
pStp=reference_p Pressure for which fluid density is defined
dStp=1.2 Fluid density at pressure pStp
Buildings.Media.Air.ThermodynamicState ThermodynamicState ThermodynamicState record for moist air
Buildings.Media.Air.BaseProperties BaseProperties Base properties (p, d, T, h, u, R, MM and X and Xi) of a medium
Buildings.Media.Air.density density Gas density
Buildings.Media.Air.dynamicViscosity dynamicViscosity Return the dynamic viscosity of dry air
Buildings.Media.Air.enthalpyOfCondensingGas enthalpyOfCondensingGas Enthalpy of steam per unit mass of steam
Buildings.Media.Air.enthalpyOfGas enthalpyOfGas Enthalpy of gas mixture per unit mass of gas mixture
Buildings.Media.Air.enthalpyOfLiquid enthalpyOfLiquid Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature
Buildings.Media.Air.enthalpyOfNonCondensingGas enthalpyOfNonCondensingGas Enthalpy of non-condensing gas per unit mass of steam
Buildings.Media.Air.enthalpyOfVaporization enthalpyOfVaporization Enthalpy of vaporization of water
Buildings.Media.Air.gasConstant gasConstant Return ideal gas constant as a function from thermodynamic state, only valid for phi<1
Buildings.Media.Air.pressure pressure Returns pressure of ideal gas as a function of the thermodynamic state record
Buildings.Media.Air.isobaricExpansionCoefficient isobaricExpansionCoefficient Isobaric expansion coefficient beta
Buildings.Media.Air.isothermalCompressibility isothermalCompressibility Isothermal compressibility factor
Buildings.Media.Air.saturationPressure saturationPressure Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)
Buildings.Media.Air.specificEntropy specificEntropy Return the specific entropy, only valid for phi<1
Buildings.Media.Air.density_derp_T density_derp_T Return the partial derivative of density with respect to pressure at constant temperature
Buildings.Media.Air.density_derT_p density_derT_p Return the partial derivative of density with respect to temperature at constant pressure
Buildings.Media.Air.density_derX density_derX Return the partial derivative of density with respect to mass fractions at constant pressure and temperature
Buildings.Media.Air.specificHeatCapacityCp specificHeatCapacityCp Specific heat capacity of gas mixture at constant pressure
Buildings.Media.Air.specificHeatCapacityCv specificHeatCapacityCv Specific heat capacity of gas mixture at constant volume
Buildings.Media.Air.setState_dTX setState_dTX Return thermodynamic state as function of density d, temperature T and composition X
Buildings.Media.Air.setState_phX setState_phX Return thermodynamic state as function of pressure p, specific enthalpy h and composition X
Buildings.Media.Air.setState_pTX setState_pTX Return thermodynamic state as function of p, T and composition X or Xi
Buildings.Media.Air.setState_psX setState_psX Return the thermodynamic state as function of p, s and composition X or Xi
Buildings.Media.Air.specificEnthalpy specificEnthalpy Compute specific enthalpy from pressure, temperature and mass fraction
Buildings.Media.Air.specificEnthalpy_pTX specificEnthalpy_pTX Specific enthalpy
Buildings.Media.Air.specificGibbsEnergy specificGibbsEnergy Specific Gibbs energy
Buildings.Media.Air.specificHelmholtzEnergy specificHelmholtzEnergy Specific Helmholtz energy
Buildings.Media.Air.isentropicEnthalpy isentropicEnthalpy Return the isentropic enthalpy
Buildings.Media.Air.specificInternalEnergy specificInternalEnergy Specific internal energy
Buildings.Media.Air.temperature temperature Return temperature of ideal gas as a function of the thermodynamic state record
Buildings.Media.Air.molarMass molarMass Return the molar mass
Buildings.Media.Air.temperature_phX temperature_phX Compute temperature from specific enthalpy and mass fraction
Buildings.Media.Air.thermalConductivity thermalConductivity Thermal conductivity of dry air as a polynomial in the temperature
Buildings.Media.Air.GasProperties GasProperties Coefficient data record for properties of perfect gases
dryair Dry air properties
steam Steam properties
k_mair=steam.MM/dryair.MM Ratio of molar weights
MMX={steam.MM,dryair.MM} Molar masses of components
h_fg=Buildings.Utilities.Psychrometrics.Constants.h_fg Latent heat of evaporation of water
cpWatLiq=Buildings.Utilities.Psychrometrics.Constants.cpWatLiq Specific heat capacity of liquid water
Buildings.Media.Air.der_enthalpyOfLiquid der_enthalpyOfLiquid Temperature derivative of enthalpy of liquid per unit mass of liquid
Buildings.Media.Air.der_enthalpyOfCondensingGas der_enthalpyOfCondensingGas Derivative of enthalpy of steam per unit mass of steam
Buildings.Media.Air.enthalpyOfDryAir enthalpyOfDryAir Enthalpy of dry air per unit mass of dry air
Buildings.Media.Air.der_enthalpyOfDryAir der_enthalpyOfDryAir Derivative of enthalpy of dry air per unit mass of dry air
Buildings.Media.Air.der_enthalpyOfNonCondensingGas der_enthalpyOfNonCondensingGas Derivative of enthalpy of non-condensing gas per unit mass of steam
Buildings.Media.Air.der_specificHeatCapacityCp der_specificHeatCapacityCp Derivative of specific heat capacity of gas mixture at constant pressure
Buildings.Media.Air.der_specificHeatCapacityCv der_specificHeatCapacityCv Derivative of specific heat capacity of gas mixture at constant volume
Inherited
fluidConstants={Modelica.Media.IdealGases.Common.FluidData.H2O,Modelica.Media.IdealGases.Common.FluidData.N2} Constant data for the fluid
Modelica.Media.Interfaces.PartialMixtureMedium.moleToMassFractions moleToMassFractions Return mass fractions X from mole fractions
Modelica.Media.Interfaces.PartialMixtureMedium.massToMoleFractions massToMoleFractions Return mole fractions from mass fractions X
ThermoStates=Modelica.Media.Interfaces.Choices.IndependentVariables.pTX Enumeration type for independent variables
mediumName="Air" Name of the medium
substanceNames={"water","air"} 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=false = 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=273.15 Reference temperature of Medium: default 25 deg Celsius
reference_X={0.01,0.99} 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)
C_default=fill(0, nC) Default value for trace substances 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
Modelica.Media.Interfaces.PartialMedium.FluidConstants FluidConstants Critical, triple, molecular and other standard data of fluid
Modelica.Media.Interfaces.PartialMedium.setSmoothState setSmoothState Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b
Modelica.Media.Interfaces.PartialMedium.prandtlNumber prandtlNumber Return the Prandtl number
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp heatCapacity_cp Alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cv heatCapacity_cv Alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.isentropicExponent isentropicExponent Return isentropic exponent
Modelica.Media.Interfaces.PartialMedium.velocityOfSound velocityOfSound Return velocity of sound
Modelica.Media.Interfaces.PartialMedium.beta beta Alias for isobaricExpansionCoefficient for user convenience
Modelica.Media.Interfaces.PartialMedium.kappa kappa Alias of isothermalCompressibility for user convenience
Modelica.Media.Interfaces.PartialMedium.density_derp_h density_derp_h Return density derivative w.r.t. pressure at const specific enthalpy
Modelica.Media.Interfaces.PartialMedium.density_derh_p density_derh_p Return density derivative w.r.t. specific enthalpy at constant pressure
Modelica.Media.Interfaces.PartialMedium.specificEntropy_pTX specificEntropy_pTX Return specific enthalpy from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_pTX density_pTX Return density from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_phX density_phX Return density from p, h, and X or Xi
Modelica.Media.Interfaces.PartialMedium.temperature_psX temperature_psX Return temperature from p,s, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_psX density_psX Return density from p, s, and X or Xi
Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_psX specificEnthalpy_psX Return specific enthalpy from p, s, and X or Xi
Modelica.Media.Interfaces.PartialMedium.MassFlowRate MassFlowRate Type for mass flow rate with medium specific attributes
Modelica.Media.Interfaces.Types.AbsolutePressure AbsolutePressure Type for absolute pressure with medium specific attributes
Modelica.Media.Interfaces.Types.Density Density Type for density with medium specific attributes
Modelica.Media.Interfaces.Types.DynamicViscosity DynamicViscosity Type for dynamic viscosity with medium specific attributes
Modelica.Media.Interfaces.Types.EnthalpyFlowRate EnthalpyFlowRate Type for enthalpy flow rate with medium specific attributes
Modelica.Media.Interfaces.Types.MassFraction MassFraction Type for mass fraction with medium specific attributes
Modelica.Media.Interfaces.Types.MoleFraction MoleFraction Type for mole fraction with medium specific attributes
Modelica.Media.Interfaces.Types.MolarMass MolarMass Type for molar mass with medium specific attributes
Modelica.Media.Interfaces.Types.MolarVolume MolarVolume Type for molar volume with medium specific attributes
Modelica.Media.Interfaces.Types.IsentropicExponent IsentropicExponent Type for isentropic exponent with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEnergy SpecificEnergy Type for specific energy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificInternalEnergy SpecificInternalEnergy Type for specific internal energy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEnthalpy SpecificEnthalpy Type for specific enthalpy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEntropy SpecificEntropy Type for specific entropy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificHeatCapacity SpecificHeatCapacity Type for specific heat capacity with medium specific attributes
Modelica.Media.Interfaces.Types.SurfaceTension SurfaceTension Type for surface tension with medium specific attributes
Modelica.Media.Interfaces.Types.Temperature Temperature Type for temperature with medium specific attributes
Modelica.Media.Interfaces.Types.ThermalConductivity ThermalConductivity Type for thermal conductivity with medium specific attributes
Modelica.Media.Interfaces.Types.PrandtlNumber PrandtlNumber Type for Prandtl number with medium specific attributes
Modelica.Media.Interfaces.Types.VelocityOfSound VelocityOfSound Type for velocity of sound with medium specific attributes
Modelica.Media.Interfaces.Types.ExtraProperty ExtraProperty Type for unspecified, mass-specific property transported by flow
Modelica.Media.Interfaces.Types.CumulativeExtraProperty CumulativeExtraProperty Type for conserved integral of unspecified, mass specific property
Modelica.Media.Interfaces.Types.ExtraPropertyFlowRate ExtraPropertyFlowRate Type for flow rate of unspecified, mass-specific property
Modelica.Media.Interfaces.Types.IsobaricExpansionCoefficient IsobaricExpansionCoefficient Type for isobaric expansion coefficient with medium specific attributes
Modelica.Media.Interfaces.Types.DipoleMoment DipoleMoment Type for dipole moment with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByPressure DerDensityByPressure Type for partial derivative of density with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByEnthalpy DerDensityByEnthalpy Type for partial derivative of density with respect to enthalpy with medium specific attributes
Modelica.Media.Interfaces.Types.DerEnthalpyByPressure DerEnthalpyByPressure Type for partial derivative of enthalpy with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByTemperature DerDensityByTemperature Type for partial derivative of density with respect to temperature with medium specific attributes
Modelica.Media.Interfaces.Types.DerTemperatureByPressure DerTemperatureByPressure Type for partial derivative of temperature with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.SaturationProperties SaturationProperties Saturation properties of two phase medium
Modelica.Media.Interfaces.Types.FluidLimits FluidLimits Validity limits for fluid model
Modelica.Media.Interfaces.Types.FixedPhase FixedPhase Phase of the fluid: 1 for 1-phase, 2 for two-phase, 0 for not known, e.g., interactive use
Modelica.Media.Interfaces.Types.Basic Basic The most basic version of a record used in several degrees of detail
Modelica.Media.Interfaces.Types.IdealGas IdealGas The ideal gas version of a record used in several degrees of detail
Modelica.Media.Interfaces.Types.TwoPhase TwoPhase The two phase fluid version of a record used in several degrees of detail

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 AbsolutePressure pStp = reference_p
    "Pressure for which fluid density is defined";
  constant Density dStp = 1.2 "Fluid density at pressure pStp";
  constant GasProperties dryair(
    R =    Modelica.Media.IdealGases.Common.SingleGasesData.Air.R,
    MM =   Modelica.Media.IdealGases.Common.SingleGasesData.Air.MM,
    cp =   Buildings.Utilities.Psychrometrics.Constants.cpAir,
    cv =   Buildings.Utilities.Psychrometrics.Constants.cpAir
             -Modelica.Media.IdealGases.Common.SingleGasesData.Air.R)
    "Dry air properties";
  constant GasProperties steam(
    R =    Modelica.Media.IdealGases.Common.SingleGasesData.H2O.R,
    MM =   Modelica.Media.IdealGases.Common.SingleGasesData.H2O.MM,
    cp =   Buildings.Utilities.Psychrometrics.Constants.cpSte,
    cv =   Buildings.Utilities.Psychrometrics.Constants.cpSte
             -Modelica.Media.IdealGases.Common.SingleGasesData.H2O.R)
    "Steam properties";
  constant Real k_mair =  steam.MM/dryair.MM "Ratio of molar weights";
  constant Modelica.SIunits.MolarMass[2] MMX={steam.MM,dryair.MM}
    "Molar masses of components";
  constant Modelica.SIunits.SpecificEnergy h_fg=
    Buildings.Utilities.Psychrometrics.Constants.h_fg
    "Latent heat of evaporation of water";
  constant Modelica.SIunits.SpecificHeatCapacity cpWatLiq=
    Buildings.Utilities.Psychrometrics.Constants.cpWatLiq
    "Specific heat capacity of liquid water";

Buildings.Media.Air.ThermodynamicState Buildings.Media.Air.ThermodynamicState

ThermodynamicState record for moist air

Information

Extends from (Thermodynamic state variables).

Modelica definition

redeclare record extends ThermodynamicState "ThermodynamicState record for moist air" end ThermodynamicState;

Buildings.Media.Air.BaseProperties Buildings.Media.Air.BaseProperties

Base properties (p, d, T, h, u, R, MM and X and Xi) of a medium

Information

Model with basic thermodynamic properties.

This model provides equation for the following thermodynamic properties:

Variable Unit Description
T K temperature
p Pa absolute pressure
d kg/m3 density
h J/kg specific enthalpy
u J/kg specific internal energy
Xi[nXi] kg/kg independent mass fractions m_i/m
R J/kg.K gas constant
M kg/mol molar mass

Parameters

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

Modelica definition

redeclare model BaseProperties "Base properties (p, d, T, h, u, R, MM and X and Xi) of a medium" parameter Boolean preferredMediumStates=false "= true if StateSelect.prefer shall be used for the independent property variables of the medium"; final parameter Boolean standardOrderComponents=true "If true, and reducedX = true, the last element of X will be computed from the other ones"; InputAbsolutePressure p "Absolute pressure of medium"; InputMassFraction[1] Xi( start=reference_X[1:1], each stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default) "Structurally independent mass fractions"; InputSpecificEnthalpy h "Specific enthalpy of medium"; Modelica.Media.Interfaces.Types.Density d "Density of medium"; Modelica.Media.Interfaces.Types.Temperature T( stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default) "Temperature of medium"; Modelica.Media.Interfaces.Types.MassFraction[2] X(start=reference_X) "Mass fractions (= (component mass)/total mass m_i/m)"; Modelica.Media.Interfaces.Types.SpecificInternalEnergy u "Specific internal energy of medium"; Modelica.Media.Interfaces.Types.SpecificHeatCapacity R "Gas constant (of mixture if applicable)"; Modelica.Media.Interfaces.Types.MolarMass MM "Molar mass (of mixture or single fluid)"; ThermodynamicState state "Thermodynamic state record for optional functions"; Modelica.SIunits.Conversions.NonSIunits.Temperature_degC T_degC= Modelica.SIunits.Conversions.to_degC(T) "Temperature of medium in [degC]"; Modelica.SIunits.Conversions.NonSIunits.Pressure_bar p_bar= Modelica.SIunits.Conversions.to_bar(p) "Absolute pressure of medium in [bar]"; // Local connector definition, used for equation balancing check connector InputAbsolutePressure = input Modelica.SIunits.AbsolutePressure "Pressure as input signal connector"; connector InputSpecificEnthalpy = input Modelica.SIunits.SpecificEnthalpy "Specific enthalpy as input signal connector"; connector InputMassFraction = input Modelica.SIunits.MassFraction "Mass fraction as input signal connector"; // Declarations for Air only protected Modelica.SIunits.TemperatureDifference dT(start=T_default-reference_T) "Temperature difference used to compute enthalpy"; equation MM = 1/(X[1]/steam.MM+(X[2])/dryair.MM); dT = T - reference_T; h = dT*dryair.cp * X[2] + (dT * steam.cp + h_fg) * X[1]; R = dryair.R*X[2] + steam.R*X[1]; // 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; // In this medium model, the density depends only // on temperature, but not on pressure. // d = p/(R*T); d/dStp = p/pStp; state.p = p; state.T = T; state.X = X; X[1] = Xi[1]; X[2] = 1 - X[1]; // Assertions to test for bounds assert(noEvent(X[1] >= -1.e-5) and noEvent(X[1] <= 1 + 1.e-5), "Mass fraction X[1] = " + String(X[1]) + " of substance water" + "\nof medium \"Buildings.Media.Air\" is not in the range 0..1"); assert(noEvent(T >= 200.0), "In " + getInstanceName() + ": Temperature T exceeded its minimum allowed value of -73.15 degC (200 Kelvin) as required from medium model \"Buildings.Media.Air\"."); assert(noEvent(T <= 423.15), "In " + getInstanceName() + ": Temperature T exceeded its maximum allowed value of 150 degC (423.15 Kelvin) as required from medium model \"Buildings.Media.Air\"."); assert(noEvent(p >= 0.0), "Pressure (= " + String(p) + " Pa) of medium \"Buildings.Media.Air\" is negative\n(Temperature = " + String(T) + " K)"); end BaseProperties;

Buildings.Media.Air.density Buildings.Media.Air.density

Gas density

Information

Density is computed from pressure, temperature and composition in the thermodynamic state record applying the ideal gas law.

Extends from Modelica.Icons.Function (Icon for functions).

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*dStp/pStp; end density;

Buildings.Media.Air.dynamicViscosity Buildings.Media.Air.dynamicViscosity

Return the dynamic viscosity of dry air

Information

This function returns the dynamic viscosity.

Implementation

The function is based on the 5th order polynomial of Modelica.Media.Air.MoistAir.dynamicViscosity. However, for the typical range of temperatures encountered in building applications, a linear function sufficies. This implementation is therefore the above 5th order polynomial, linearized around 20°C. The relative error of this linearization is 0.4% at -20°C, and less then 0.2% between -5°C and +50°C.

Extends from (Return dynamic viscosity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

redeclare function extends dynamicViscosity "Return the dynamic viscosity of dry air" algorithm eta := 4.89493640395e-08 * state.T + 3.88335940547e-06; end dynamicViscosity;

Buildings.Media.Air.enthalpyOfCondensingGas Buildings.Media.Air.enthalpyOfCondensingGas

Enthalpy of steam per unit mass of steam

Information

Extends from Modelica.Icons.Function (Icon for functions).

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-reference_T) * steam.cp + h_fg; end enthalpyOfCondensingGas;

Buildings.Media.Air.enthalpyOfGas Buildings.Media.Air.enthalpyOfGas

Enthalpy of gas mixture per unit mass of gas mixture

Information

Extends from (Return enthalpy of non-condensing gas mixture).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
MassFractionX[:] Vector of mass fractions [kg/kg]

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific 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.Air.enthalpyOfLiquid Buildings.Media.Air.enthalpyOfLiquid

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

Information

Extends from (Return liquid enthalpy of condensing fluid).

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 - reference_T)*cpWatLiq; end enthalpyOfLiquid;

Buildings.Media.Air.enthalpyOfNonCondensingGas Buildings.Media.Air.enthalpyOfNonCondensingGas

Enthalpy of non-condensing gas per unit mass of steam

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhenthalpy [J/kg]

Modelica definition

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.Air.enthalpyOfVaporization Buildings.Media.Air.enthalpyOfVaporization

Enthalpy of vaporization of water

Information

Extends from (Return vaporization enthalpy of condensing fluid).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyr0Vaporization enthalpy [J/kg]

Modelica definition

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

Buildings.Media.Air.gasConstant Buildings.Media.Air.gasConstant

Return ideal gas constant as a function from thermodynamic state, only valid for phi<1

Information

The ideal gas constant for moist air is computed from thermodynamic state assuming that all water is in the gas phase.

Extends from (Return the gas constant of the mixture (also for liquids)).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state

Outputs

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

Modelica definition

redeclare function extends gasConstant "Return ideal gas constant as a function from thermodynamic state, only valid for phi<1" algorithm R := dryair.R*(1 - state.X[Water]) + steam.R*state.X[Water]; end gasConstant;

Buildings.Media.Air.pressure Buildings.Media.Air.pressure

Returns pressure of ideal gas as a function of the thermodynamic state record

Information

Pressure is returned from the thermodynamic state record input as a simple assignment.

Extends from (Return pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function extends pressure "Returns pressure of ideal gas as a function of the thermodynamic state record" algorithm p := state.p; end pressure;

Buildings.Media.Air.isobaricExpansionCoefficient Buildings.Media.Air.isobaricExpansionCoefficient

Isobaric expansion coefficient beta

Information

This function returns the isobaric expansion coefficient at constant pressure, which is zero for this medium. The isobaric expansion coefficient at constant pressure is

βp = - 1 ⁄ v   (∂ v ⁄ ∂ T)p = 0,

where v is the specific volume, T is the temperature and p is the pressure.

Extends from (Return overall the isobaric expansion coefficient beta).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
IsobaricExpansionCoefficientbetaIsobaric expansion coefficient [1/K]

Modelica definition

redeclare function extends isobaricExpansionCoefficient "Isobaric expansion coefficient beta" algorithm beta := 0; end isobaricExpansionCoefficient;

Buildings.Media.Air.isothermalCompressibility Buildings.Media.Air.isothermalCompressibility

Isothermal compressibility factor

Information

This function returns the isothermal compressibility coefficient. The isothermal compressibility is

κT = -1 ⁄ v   (∂ v ⁄ ∂ p)T = -1 ⁄ p,

where v is the specific volume, T is the temperature and p is the pressure.

Extends from (Return overall the isothermal compressibility factor).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
IsothermalCompressibilitykappaIsothermal compressibility [1/Pa]

Modelica definition

redeclare function extends isothermalCompressibility "Isothermal compressibility factor" algorithm kappa := -1/state.p; end isothermalCompressibility;

Buildings.Media.Air.saturationPressure Buildings.Media.Air.saturationPressure

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

Information

Extends from (Return saturation pressure of condensing fluid).

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.Psychrometrics.Functions.saturationPressure(Tsat); end saturationPressure;

Buildings.Media.Air.specificEntropy Buildings.Media.Air.specificEntropy

Return the specific entropy, only valid for phi<1

Information

This function computes the specific entropy.

The specific entropy of the mixture is obtained from

s = ss + sm,

where ss is the entropy change due to the state change (relative to the reference temperature) and sm is the entropy change due to mixing of the dry air and water vapor.

The entropy change due to change in state is obtained from

ss = cv ln(T/T0) + R ln(v/v0)
= cv ln(T/T0) + R ln(ρ0/ρ)

If we assume ρ = p0/(R T), and because cp = cv + R, we can write

ss = cv ln(T/T0) + R ln(T/T0)
=cp ln(T/T0).

Next, the entropy of mixing is obtained from a reversible isothermal expansion process. Hence,

sm = -R ∑i( Xi ⁄ Mi ln(Yi p/p0)),

where R is the gas constant, X is the mass fraction, M is the molar mass, and Y is the mole fraction.

To obtain the state for a given pressure, entropy and mass fraction, use Buildings.Media.Air.setState_psX.

Limitations

This function is only valid for a relative humidity below 100%.

Extends from (Return specific entropy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends specificEntropy "Return the specific entropy, only valid for phi<1" protected Modelica.SIunits.MoleFraction[2] Y "Molar fraction"; algorithm Y := massToMoleFractions( state.X, {steam.MM,dryair.MM}); s := specificHeatCapacityCp(state) * Modelica.Math.log(state.T/reference_T) - Modelica.Constants.R * sum(state.X[i]/MMX[i]* Modelica.Math.log(max(Y[i], Modelica.Constants.eps)*state.p/reference_p) for i in 1:2); end specificEntropy;

Buildings.Media.Air.density_derp_T Buildings.Media.Air.density_derp_T

Return the partial derivative of density with respect to pressure at constant temperature

Information

This function returns the partial derivative of density with respect to pressure at constant temperature.

Extends from (Return density derivative w.r.t. pressure at const temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DerDensityByPressureddpTDensity derivative w.r.t. pressure [s2/m2]

Modelica definition

redeclare function extends density_derp_T "Return the partial derivative of density with respect to pressure at constant temperature" algorithm ddpT := dStp/pStp; end density_derp_T;

Buildings.Media.Air.density_derT_p Buildings.Media.Air.density_derT_p

Return the partial derivative of density with respect to temperature at constant pressure

Information

This function computes the derivative of density with respect to temperature at constant pressure.

Extends from (Return density derivative w.r.t. temperature at constant pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DerDensityByTemperatureddTpDensity derivative w.r.t. temperature [kg/(m3.K)]

Modelica definition

redeclare function extends density_derT_p "Return the partial derivative of density with respect to temperature at constant pressure" algorithm ddTp := 0; end density_derT_p;

Buildings.Media.Air.density_derX Buildings.Media.Air.density_derX

Return the partial derivative of density with respect to mass fractions at constant pressure and temperature

Information

This function returns the partial derivative of density with respect to mass fraction. This value is zero because in this medium, density is proportional to pressure, but independent of the species concentration.

Extends from (Return density derivative w.r.t. mass fraction).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DensitydddX[nX]Derivative of density w.r.t. mass fraction [kg/m3]

Modelica definition

redeclare function extends density_derX "Return the partial derivative of density with respect to mass fractions at constant pressure and temperature" algorithm dddX := fill(0, nX); end density_derX;

Buildings.Media.Air.specificHeatCapacityCp Buildings.Media.Air.specificHeatCapacityCp

Specific heat capacity of gas mixture at constant pressure

Information

Extends from (Return specific heat capacity at constant pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

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" annotation(derivative=der_specificHeatCapacityCp); algorithm cp := dryair.cp*(1-state.X[Water]) +steam.cp*state.X[Water]; end specificHeatCapacityCp;

Buildings.Media.Air.specificHeatCapacityCv Buildings.Media.Air.specificHeatCapacityCv

Specific heat capacity of gas mixture at constant volume

Information

Extends from (Return specific heat capacity at constant volume).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

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" annotation(derivative=der_specificHeatCapacityCv); algorithm cv:= dryair.cv*(1-state.X[Water]) +steam.cv*state.X[Water]; end specificHeatCapacityCv;

Buildings.Media.Air.setState_dTX Buildings.Media.Air.setState_dTX

Return thermodynamic state as function of density d, temperature T and composition X

Information

The thermodynamic state record is computed from density d, temperature T and composition X.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Densityd Density [kg/m3]
TemperatureT Temperature [K]
MassFractionX[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state

Modelica definition

redeclare function setState_dTX "Return thermodynamic state as function of density d, temperature T and composition X" extends Modelica.Icons.Function; input Density d "Density"; input Temperature T "Temperature"; input MassFraction X[:]=reference_X "Mass fractions"; output ThermodynamicState state "Thermodynamic state"; algorithm // Note that d/dStp = p/pStp, hence p = d*pStp/dStp state := if size(X, 1) == nX then ThermodynamicState(p=d*pStp/dStp, T=T, X=X) else ThermodynamicState(p=d*pStp/dStp, T=T, X=cat(1, X, {1 - sum(X)})); end setState_dTX;

Buildings.Media.Air.setState_phX Buildings.Media.Air.setState_phX

Return thermodynamic state as function of pressure p, specific enthalpy h and composition X

Information

The thermodynamic state record is computed from pressure p, specific enthalpy h and composition X.

Extends from (Return thermodynamic state as function of p, h and composition X or Xi).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_phX "Return thermodynamic state as function of pressure p, specific enthalpy h and composition X" algorithm state := if size(X, 1) == nX then ThermodynamicState(p=p, T=temperature_phX(p, h, X), X=X) else ThermodynamicState(p=p, T=temperature_phX(p, h, X), X=cat(1, X, {1 - sum(X)})); end setState_phX;

Buildings.Media.Air.setState_pTX Buildings.Media.Air.setState_pTX

Return thermodynamic state as function of p, T and composition X or Xi

Information

The thermodynamic state record is computed from pressure p, temperature T and composition X.

Extends from (Return thermodynamic state as function of p, T and composition X or Xi).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_pTX "Return thermodynamic state as function of p, T and composition X or Xi" algorithm state := if size(X, 1) == nX then ThermodynamicState(p=p, T=T, X=X) else ThermodynamicState(p=p, T=T, X=cat(1, X, {1 - sum(X)})); end setState_pTX;

Buildings.Media.Air.setState_psX Buildings.Media.Air.setState_psX

Return the thermodynamic state as function of p, s and composition X or Xi

Information

This function returns the thermodynamic state based on pressure, specific entropy and mass fraction.

The state is computed by symbolically solving Buildings.Media.Air.specificEntropy for temperature.

Extends from (Return thermodynamic state as function of p, s and composition X or Xi).

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEntropys Specific entropy [J/(kg.K)]
MassFractionX[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_psX "Return the thermodynamic state as function of p, s and composition X or Xi" protected Modelica.SIunits.MassFraction[2] X_int "Mass fraction"; Modelica.SIunits.MoleFraction[2] Y "Molar fraction"; Modelica.SIunits.Temperature T "Temperature"; algorithm if size(X, 1) == nX then X_int:=X; else X_int :=cat( 1, X, {1 - sum(X)}); end if; Y := massToMoleFractions( X_int, {steam.MM,dryair.MM}); // The next line is obtained from symbolic solving the // specificEntropy function for T. // In this formulation, we can set T to any value when calling // specificHeatCapacityCp as cp does not depend on T. T := 273.15 * Modelica.Math.exp((s + Modelica.Constants.R * sum(X_int[i]/MMX[i]* Modelica.Math.log(max(Y[i], Modelica.Constants.eps)) for i in 1:2)) / specificHeatCapacityCp(setState_pTX(p=p, T=273.15, X=X_int))); state := ThermodynamicState(p=p, T=T, X=X_int); end setState_psX;

Buildings.Media.Air.specificEnthalpy Buildings.Media.Air.specificEnthalpy

Compute specific enthalpy from pressure, temperature and mass fraction

Information

Extends from (Return specific enthalpy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

redeclare replaceable function extends specificEnthalpy "Compute specific enthalpy from pressure, temperature and mass fraction" algorithm h := (state.T - reference_T)*dryair.cp * (1 - state.X[Water]) + ((state.T-reference_T) * steam.cp + h_fg) * state.X[Water]; end specificEnthalpy;

Buildings.Media.Air.specificEnthalpy_pTX Buildings.Media.Air.specificEnthalpy_pTX

Specific enthalpy

Information

Specific enthalpy as a function of temperature and species concentration. The pressure is input for compatibility with the medium models, but the specific enthalpy is independent of the pressure.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

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

Outputs

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

Modelica definition

redeclare replaceable function specificEnthalpy_pTX "Specific enthalpy" extends Modelica.Icons.Function; input Modelica.SIunits.Pressure p "Pressure"; input Modelica.SIunits.Temperature T "Temperature"; input Modelica.SIunits.MassFraction X[:] "Mass fractions of moist air"; output Modelica.SIunits.SpecificEnthalpy h "Specific enthalpy at p, T, X"; algorithm h := specificEnthalpy(setState_pTX(p, T, X)); end specificEnthalpy_pTX;

Buildings.Media.Air.specificGibbsEnergy Buildings.Media.Air.specificGibbsEnergy

Specific Gibbs energy

Information

Extends from (Return specific Gibbs energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

redeclare replaceable function extends specificGibbsEnergy "Specific Gibbs energy" algorithm g := specificEnthalpy(state) - state.T*specificEntropy(state); end specificGibbsEnergy;

Buildings.Media.Air.specificHelmholtzEnergy Buildings.Media.Air.specificHelmholtzEnergy

Specific Helmholtz energy

Information

Extends from (Return specific Helmholtz energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

redeclare replaceable function extends specificHelmholtzEnergy "Specific Helmholtz energy" algorithm f := specificEnthalpy(state) - gasConstant(state)*state.T - state.T*specificEntropy(state); end specificHelmholtzEnergy;

Buildings.Media.Air.isentropicEnthalpy Buildings.Media.Air.isentropicEnthalpy

Return the isentropic enthalpy

Information

This function computes the specific enthalpy for an isentropic state change from the temperature that corresponds to the state refState to reference_T.

Extends from (Return isentropic enthalpy).

Inputs

TypeNameDefaultDescription
AbsolutePressurep_downstream Downstream pressure [Pa]
ThermodynamicStaterefState Reference state for entropy

Outputs

TypeNameDescription
SpecificEnthalpyh_isIsentropic enthalpy [J/kg]

Modelica definition

redeclare function extends isentropicEnthalpy "Return the isentropic enthalpy" algorithm h_is := specificEnthalpy(setState_psX( p=p_downstream, s=specificEntropy(refState), X=refState.X)); end isentropicEnthalpy;

Buildings.Media.Air.specificInternalEnergy Buildings.Media.Air.specificInternalEnergy

Specific internal energy

Information

Extends from Modelica.Icons.Function (Icon for functions), (Return specific internal energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyuSpecific internal energy [J/kg]

Modelica definition

redeclare function extends specificInternalEnergy "Specific internal energy" extends Modelica.Icons.Function; algorithm u := specificEnthalpy(state) - pStp/dStp; end specificInternalEnergy;

Buildings.Media.Air.temperature Buildings.Media.Air.temperature

Return temperature of ideal gas as a function of the thermodynamic state record

Information

Temperature is returned from the thermodynamic state record input as a simple assignment.

Extends from (Return temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function extends temperature "Return temperature of ideal gas as a function of the thermodynamic state record" algorithm T := state.T; end temperature;

Buildings.Media.Air.molarMass Buildings.Media.Air.molarMass

Return the molar mass

Information

This function returns the molar mass.

Extends from (Return the molar mass of the medium).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
MolarMassMMMixture molar mass [kg/mol]

Modelica definition

redeclare function extends molarMass "Return the molar mass" algorithm MM := 1/(state.X[Water]/MMX[Water]+(1.0-state.X[Water])/MMX[Air]); end molarMass;

Buildings.Media.Air.temperature_phX Buildings.Media.Air.temperature_phX

Compute temperature from specific enthalpy and mass fraction

Information

Temperature as a function of specific enthalpy and species concentration. The pressure is input for compatibility with the medium models, but the temperature is independent of the pressure.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

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

Outputs

TypeNameDescription
TemperatureTtemperature [K]

Modelica definition

redeclare replaceable function temperature_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 := reference_T + (h - h_fg * X[Water]) /((1 - X[Water])*dryair.cp + X[Water] * steam.cp); end temperature_phX;

Buildings.Media.Air.thermalConductivity Buildings.Media.Air.thermalConductivity

Thermal conductivity of dry air as a polynomial in the temperature

Information

Extends from (Return thermal conductivity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

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 := Modelica.Media.Incompressible.TableBased.Polynomials_Temp.evaluate( {(-4.8737307422969E-008), 7.67803133753502E-005, 0.0241814385504202}, Modelica.SIunits.Conversions.to_degC(state.T)); end thermalConductivity;

Buildings.Media.Air.GasProperties Buildings.Media.Air.GasProperties

Coefficient data record for properties of perfect gases

Information

This data record contains the coefficients for perfect gases.

Extends from Modelica.Icons.Record (Icon for records).

Modelica definition

record GasProperties "Coefficient data record for properties of perfect gases" extends Modelica.Icons.Record; Modelica.SIunits.MolarMass MM "Molar mass"; Modelica.SIunits.SpecificHeatCapacity R "Gas constant"; Modelica.SIunits.SpecificHeatCapacity cp "Specific heat capacity at constant pressure"; Modelica.SIunits.SpecificHeatCapacity cv = cp-R "Specific heat capacity at constant volume"; end GasProperties;

Buildings.Media.Air.der_enthalpyOfLiquid Buildings.Media.Air.der_enthalpyOfLiquid

Temperature derivative of enthalpy of liquid per unit mass of liquid

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
Realder_T Temperature derivative

Outputs

TypeNameDescription
Realder_hDerivative of liquid enthalpy

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 Real der_T "Temperature derivative"; output Real der_h "Derivative of liquid enthalpy"; algorithm der_h := cpWatLiq*der_T; end der_enthalpyOfLiquid;

Buildings.Media.Air.der_enthalpyOfCondensingGas Buildings.Media.Air.der_enthalpyOfCondensingGas

Derivative of enthalpy of steam per unit mass of steam

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
Realder_T Temperature derivative

Outputs

TypeNameDescription
Realder_hDerivative of steam enthalpy

Modelica definition

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

Enthalpy of dry air per unit mass of dry air

Information

Extends from Modelica.Icons.Function (Icon for functions).

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 - reference_T)*dryair.cp; end enthalpyOfDryAir;

Buildings.Media.Air.der_enthalpyOfDryAir Buildings.Media.Air.der_enthalpyOfDryAir

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

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
Realder_T Temperature derivative

Outputs

TypeNameDescription
Realder_hDerivative of dry air enthalpy

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 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.Air.der_enthalpyOfNonCondensingGas Buildings.Media.Air.der_enthalpyOfNonCondensingGas

Derivative of enthalpy of non-condensing gas per unit mass of steam

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
Realder_T Temperature derivative

Outputs

TypeNameDescription
Realder_hDerivative of steam enthalpy

Modelica definition

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.Air.der_specificHeatCapacityCp Buildings.Media.Air.der_specificHeatCapacityCp

Derivative of specific heat capacity of gas mixture at constant pressure

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state
ThermodynamicStateder_state Derivative of thermodynamic state

Outputs

TypeNameDescription
Realder_cpDerivative of specific heat capacity [J/(kg.K.s)]

Modelica definition

replaceable function der_specificHeatCapacityCp "Derivative of specific heat capacity of gas mixture at constant pressure" extends Modelica.Icons.Function; input ThermodynamicState state "Thermodynamic state"; input ThermodynamicState der_state "Derivative of thermodynamic state"; output Real der_cp(unit="J/(kg.K.s)") "Derivative of specific heat capacity"; algorithm der_cp := (steam.cp-dryair.cp)*der_state.X[Water]; end der_specificHeatCapacityCp;

Buildings.Media.Air.der_specificHeatCapacityCv Buildings.Media.Air.der_specificHeatCapacityCv

Derivative of specific heat capacity of gas mixture at constant volume

Information

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state
ThermodynamicStateder_state Derivative of thermodynamic state

Outputs

TypeNameDescription
Realder_cvDerivative of specific heat capacity [J/(kg.K.s)]

Modelica definition

replaceable function der_specificHeatCapacityCv "Derivative of specific heat capacity of gas mixture at constant volume" extends Modelica.Icons.Function; input ThermodynamicState state "Thermodynamic state"; input ThermodynamicState der_state "Derivative of thermodynamic state"; output Real der_cv(unit="J/(kg.K.s)") "Derivative of specific heat capacity"; algorithm der_cv := (steam.cv-dryair.cv)*der_state.X[Water]; end der_specificHeatCapacityCv;

Buildings.Media.Air.BaseProperties.InputAbsolutePressure

Pressure as input signal connector

Modelica definition

connector InputAbsolutePressure = input Modelica.SIunits.AbsolutePressure "Pressure as input signal connector";

Buildings.Media.Air.BaseProperties.InputMassFraction

Mass fraction as input signal connector

Modelica definition

connector InputMassFraction = input Modelica.SIunits.MassFraction "Mass fraction as input signal connector";

Buildings.Media.Air.BaseProperties.InputSpecificEnthalpy

Specific enthalpy as input signal connector

Modelica definition

connector InputSpecificEnthalpy = input Modelica.SIunits.SpecificEnthalpy "Specific enthalpy as input signal connector";