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 |
ThermodynamicState | ThermodynamicState record for moist air |
BaseProperties | Base properties (p, d, T, h, u, R, MM and X and Xi) of a medium |
density | Gas density |
dynamicViscosity | Return the dynamic viscosity of dry air |
enthalpyOfCondensingGas | Enthalpy of steam per unit mass of steam |
enthalpyOfGas | Enthalpy of gas mixture per unit mass of gas mixture |
enthalpyOfLiquid | Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature |
enthalpyOfNonCondensingGas | Enthalpy of non-condensing gas per unit mass of steam |
enthalpyOfVaporization | Enthalpy of vaporization of water |
gasConstant | Return ideal gas constant as a function from thermodynamic state, only valid for phi<1 |
pressure | Returns pressure of ideal gas as a function of the thermodynamic state record |
isobaricExpansionCoefficient | Isobaric expansion coefficient beta |
isothermalCompressibility | Isothermal compressibility factor |
saturationPressure | Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy) |
specificEntropy | Return the specific entropy, only valid for phi<1 |
density_derp_T | Return the partial derivative of density with respect to pressure at constant temperature |
density_derT_p | Return the partial derivative of density with respect to temperature at constant pressure |
density_derX | Return the partial derivative of density with respect to mass fractions at constant pressure and temperature |
specificHeatCapacityCp | Specific heat capacity of gas mixture at constant pressure |
specificHeatCapacityCv | Specific heat capacity of gas mixture at constant volume |
setState_dTX | Return thermodynamic state as function of density d, temperature T and composition X |
setState_phX | Return thermodynamic state as function of pressure p, specific enthalpy h and composition X |
setState_pTX | Return thermodynamic state as function of p, T and composition X or Xi |
setState_psX | Return the thermodynamic state as function of p, s and composition X or Xi |
specificEnthalpy | Compute specific enthalpy from pressure, temperature and mass fraction |
specificEnthalpy_pTX | Specific enthalpy |
specificGibbsEnergy | Specific Gibbs energy |
specificHelmholtzEnergy | Specific Helmholtz energy |
isentropicEnthalpy | Return the isentropic enthalpy |
specificInternalEnergy | Specific internal energy |
temperature | Return temperature of ideal gas as a function of the thermodynamic state record |
molarMass | Return the molar mass |
temperature_phX | Compute temperature from specific enthalpy and mass fraction |
thermalConductivity | Thermal conductivity of dry air as a polynomial in the temperature |
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 |
der_enthalpyOfLiquid | Temperature derivative of enthalpy of liquid per unit mass of liquid |
der_enthalpyOfCondensingGas | Derivative of enthalpy of steam per unit mass of steam |
enthalpyOfDryAir | Enthalpy of dry air per unit mass of dry air |
der_enthalpyOfDryAir | Derivative of enthalpy of dry air per unit mass of dry air |
der_enthalpyOfNonCondensingGas | Derivative of enthalpy of non-condensing gas per unit mass of steam |
der_specificHeatCapacityCp | Derivative of specific heat capacity of gas mixture at constant pressure |
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 |
moleToMassFractions | Return mass fractions X from mole fractions |
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 |
FluidConstants | Critical, triple, molecular and other standard data of fluid |
setSmoothState | Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b |
prandtlNumber | Return the Prandtl number |
heatCapacity_cp | Alias for deprecated name |
heatCapacity_cv | Alias for deprecated name |
isentropicExponent | Return isentropic exponent |
velocityOfSound | Return velocity of sound |
beta | Alias for isobaricExpansionCoefficient for user convenience |
kappa | Alias of isothermalCompressibility for user convenience |
density_derp_h | Return density derivative w.r.t. pressure at const specific enthalpy |
density_derh_p | Return density derivative w.r.t. specific enthalpy at constant pressure |
specificEntropy_pTX | Return specific enthalpy from p, T, and X or Xi |
density_pTX | Return density from p, T, and X or Xi |
density_phX | Return density from p, h, and X or Xi |
temperature_psX | Return temperature from p,s, and X or Xi |
density_psX | Return density from p, s, and X or Xi |
specificEnthalpy_psX | Return specific enthalpy from p, s, and X or Xi |
MassFlowRate | Type for mass flow rate with medium specific attributes |
AbsolutePressure | Type for absolute pressure with medium specific attributes |
Density | Type for density with medium specific attributes |
DynamicViscosity | Type for dynamic viscosity with medium specific attributes |
EnthalpyFlowRate | Type for enthalpy flow rate with medium specific attributes |
MassFraction | Type for mass fraction with medium specific attributes |
MoleFraction | Type for mole fraction with medium specific attributes |
MolarMass | Type for molar mass with medium specific attributes |
MolarVolume | Type for molar volume with medium specific attributes |
IsentropicExponent | Type for isentropic exponent with medium specific attributes |
SpecificEnergy | Type for specific energy with medium specific attributes |
SpecificInternalEnergy | Type for specific internal energy with medium specific attributes |
SpecificEnthalpy | Type for specific enthalpy with medium specific attributes |
SpecificEntropy | Type for specific entropy with medium specific attributes |
SpecificHeatCapacity | Type for specific heat capacity with medium specific attributes |
SurfaceTension | Type for surface tension with medium specific attributes |
Temperature | Type for temperature with medium specific attributes |
ThermalConductivity | Type for thermal conductivity with medium specific attributes |
PrandtlNumber | Type for Prandtl number with medium specific attributes |
VelocityOfSound | Type for velocity of sound with medium specific attributes |
ExtraProperty | Type for unspecified, mass-specific property transported by flow |
CumulativeExtraProperty | Type for conserved integral of unspecified, mass specific property |
ExtraPropertyFlowRate | Type for flow rate of unspecified, mass-specific property |
IsobaricExpansionCoefficient | Type for isobaric expansion coefficient with medium specific attributes |
DipoleMoment | Type for dipole moment with medium specific attributes |
DerDensityByPressure | Type for partial derivative of density with respect to pressure with medium specific attributes |
DerDensityByEnthalpy | Type for partial derivative of density with respect to enthalpy with medium specific attributes |
DerEnthalpyByPressure | Type for partial derivative of enthalpy with respect to pressure with medium specific attributes |
DerDensityByTemperature | Type for partial derivative of density with respect to temperature with medium specific attributes |
DerTemperatureByPressure | Type for partial derivative of temperature with respect to pressure with medium specific attributes |
SaturationProperties | Saturation properties of two phase medium |
FluidLimits | Validity limits for fluid model |
FixedPhase | Phase of the fluid: 1 for 1-phase, 2 for two-phase, 0 for not known, e.g., interactive use |
Basic | The most basic version of a record used in several degrees of detail |
IdealGas | The ideal gas version of a record used in several degrees of detail |
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
ThermodynamicState record for moist air
Information
Extends from (Thermodynamic state variables).
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Advanced | |||
Boolean | preferredMediumStates | false | = true if StateSelect.prefer shall be used for the independent property variables of the medium |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Outputs
Type | Name | Description |
---|---|---|
Density | d | Density [kg/m3] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
DynamicViscosity | eta | Dynamic viscosity [Pa.s] |
Modelica definition
Buildings.Media.Air.enthalpyOfCondensingGas
Enthalpy of steam per unit mass of steam
Information
Extends from Modelica.Icons.Function (Icon for functions).
Inputs
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | steam enthalpy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] | |
MassFraction | X[:] | Vector of mass fractions [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Specific enthalpy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Liquid enthalpy [J/kg] |
Modelica definition
Buildings.Media.Air.enthalpyOfNonCondensingGas
Enthalpy of non-condensing gas per unit mass of steam
Information
Extends from Modelica.Icons.Function (Icon for functions).
Inputs
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | enthalpy [J/kg] |
Modelica definition
Buildings.Media.Air.enthalpyOfVaporization
Enthalpy of vaporization of water
Information
Extends from (Return vaporization enthalpy of condensing fluid).
Inputs
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | r0 | Vaporization enthalpy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state |
Outputs
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | R | Mixture gas constant [J/(kg.K)] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
AbsolutePressure | p | Pressure [Pa] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
IsobaricExpansionCoefficient | beta | Isobaric expansion coefficient [1/K] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
IsothermalCompressibility | kappa | Isothermal compressibility [1/Pa] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | Tsat | Saturation temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
AbsolutePressure | psat | Saturation pressure [Pa] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificEntropy | s | Specific entropy [J/(kg.K)] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
DerDensityByPressure | ddpT | Density derivative w.r.t. pressure [s2/m2] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
DerDensityByTemperature | ddTp | Density derivative w.r.t. temperature [kg/(m3.K)] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
Density | dddX[nX] | Derivative of density w.r.t. mass fraction [kg/m3] |
Modelica definition
Buildings.Media.Air.specificHeatCapacityCp
Specific heat capacity of gas mixture at constant pressure
Information
Extends from (Return specific heat capacity at constant pressure).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | cp | Specific heat capacity at constant pressure [J/(kg.K)] |
Modelica definition
Buildings.Media.Air.specificHeatCapacityCv
Specific heat capacity of gas mixture at constant volume
Information
Extends from (Return specific heat capacity at constant volume).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | cv | Specific heat capacity at constant volume [J/(kg.K)] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Density | d | Density [kg/m3] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | reference_X | Mass fractions [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
ThermodynamicState | state | Thermodynamic state |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
MassFraction | X[:] | reference_X | Mass fractions [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | reference_X | Mass fractions [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
SpecificEntropy | s | Specific entropy [J/(kg.K)] | |
MassFraction | X[:] | reference_X | Mass fractions [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Modelica definition
Buildings.Media.Air.specificEnthalpy
Compute specific enthalpy from pressure, temperature and mass fraction
Information
Extends from (Return specific enthalpy).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Specific enthalpy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Pressure | p | Pressure [Pa] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | Mass fractions of moist air [1] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Specific enthalpy at p, T, X [J/kg] |
Modelica definition
Buildings.Media.Air.specificGibbsEnergy
Specific Gibbs energy
Information
Extends from (Return specific Gibbs energy).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnergy | g | Specific Gibbs energy [J/kg] |
Modelica definition
Buildings.Media.Air.specificHelmholtzEnergy
Specific Helmholtz energy
Information
Extends from (Return specific Helmholtz energy).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnergy | f | Specific Helmholtz energy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p_downstream | Downstream pressure [Pa] | |
ThermodynamicState | refState | Reference state for entropy |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h_is | Isentropic enthalpy [J/kg] |
Modelica definition
Buildings.Media.Air.specificInternalEnergy
Specific internal energy
Information
Extends from Modelica.Icons.Function (Icon for functions), (Return specific internal energy).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnergy | u | Specific internal energy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
Temperature | T | Temperature [K] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
MolarMass | MM | Mixture molar mass [kg/mol] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
SpecificEnthalpy | h | specific enthalpy [J/kg] | |
MassFraction | X[:] | mass fractions of composition [kg/kg] |
Outputs
Type | Name | Description |
---|---|---|
Temperature | T | temperature [K] |
Modelica definition
Buildings.Media.Air.thermalConductivity
Thermal conductivity of dry air as a polynomial in the temperature
Information
Extends from (Return thermal conductivity).
Inputs
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state record |
Outputs
Type | Name | Description |
---|---|---|
ThermalConductivity | lambda | Thermal conductivity [W/(m.K)] |
Modelica definition
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
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] | |
Real | der_T | Temperature derivative |
Outputs
Type | Name | Description |
---|---|---|
Real | der_h | Derivative of liquid enthalpy |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] | |
Real | der_T | Temperature derivative |
Outputs
Type | Name | Description |
---|---|---|
Real | der_h | Derivative of steam enthalpy |
Modelica definition
Buildings.Media.Air.enthalpyOfDryAir
Enthalpy of dry air per unit mass of dry air
Information
Extends from Modelica.Icons.Function (Icon for functions).
Inputs
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] |
Outputs
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Dry air enthalpy [J/kg] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] | |
Real | der_T | Temperature derivative |
Outputs
Type | Name | Description |
---|---|---|
Real | der_h | Derivative of dry air enthalpy |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | Temperature [K] | |
Real | der_T | Temperature derivative |
Outputs
Type | Name | Description |
---|---|---|
Real | der_h | Derivative of steam enthalpy |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state | |
ThermodynamicState | der_state | Derivative of thermodynamic state |
Outputs
Type | Name | Description |
---|---|---|
Real | der_cp | Derivative of specific heat capacity [J/(kg.K.s)] |
Modelica definition
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
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | Thermodynamic state | |
ThermodynamicState | der_state | Derivative of thermodynamic state |
Outputs
Type | Name | Description |
---|---|---|
Real | der_cv | Derivative of specific heat capacity [J/(kg.K.s)] |
Modelica definition
Buildings.Media.Air.BaseProperties.InputAbsolutePressure
Pressure as input signal connector
Modelica definition
Buildings.Media.Air.BaseProperties.InputMassFraction
Mass fraction as input signal connector
Modelica definition
Buildings.Media.Air.BaseProperties.InputSpecificEnthalpy
Specific enthalpy as input signal connector