This is a medium model that is similar to Buildings.Media.PerfectGases.MoistAir but in this model, the air must not be saturated. If the air is saturated, an assert will be triggered and the simulation stops with an error message. If this happens, use the medium model Buildings.Media.PerfectGases.MoistAir instead of this one.
This medium model has been added to allow an explicit computation of the function T_phX so that it is once differentiable in h with a continuous derivative. This allows obtaining an analytic expression for the Jacobian, and therefore simplifies the computation of initial conditions that can be numerically challenging for thermal fluid systems.
Name | Description |
---|---|
Water=1 | Index of water (in substanceNames, massFractions X, etc.) |
Air=2 | Index of air (in substanceNames, massFractions X, etc.) |
k_mair=steam.MM/dryair.MM | ratio of molar weights |
dryair=Buildings.Media.PerfectGases.Common.SingleGasData.Air | |
steam=Buildings.Media.PerfectGases.Common.SingleGasData.H2O | |
BaseProperties | |
Xsaturation | Steam water mass fraction of saturation boundary in kg_water/kg_moistair |
setState_pTX | Thermodynamic state as function of p, T and composition X |
setState_phX | Thermodynamic state as function of p, h and composition X |
setState_dTX | Thermodynamic state as function of d, T and composition X |
gasConstant | Gas constant |
saturationPressureLiquid | Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned |
sublimationPressureIce | Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned |
saturationPressure | Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy) |
pressure | Gas pressure |
temperature | Gas temperature |
density | Gas density |
specificEntropy | Specific entropy (liquid part neglected, mixing entropy included) |
enthalpyOfVaporization | Enthalpy of vaporization of water |
HeatCapacityOfWater | Specific heat capacity of water (liquid only) which is constant |
enthalpyOfLiquid | Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature |
der_enthalpyOfLiquid | Temperature derivative of enthalpy of liquid per unit mass of liquid |
enthalpyOfCondensingGas | Enthalpy of steam per unit mass of steam |
der_enthalpyOfCondensingGas | Derivative of enthalpy of steam per unit mass of steam |
enthalpyOfGas | Enthalpy of gas mixture per unit mass of gas mixture |
enthalpyOfDryAir | Enthalpy of dry air per unit mass of dry air |
der_enthalpyOfDryAir | Derivative of enthalpy of dry air per unit mass of dry air |
specificHeatCapacityCp | Specific heat capacity of gas mixture at constant pressure |
specificHeatCapacityCv | Specific heat capacity of gas mixture at constant volume |
dynamicViscosity | dynamic viscosity of dry air |
thermalConductivity | Thermal conductivity of dry air as a polynomial in the temperature |
h_pTX | Compute specific enthalpy from pressure, temperature and mass fraction |
specificEnthalpy | Specific enthalpy |
specificInternalEnergy | Specific internal energy |
specificGibbsEnergy | Specific Gibbs energy |
specificHelmholtzEnergy | Specific Helmholtz energy |
T_phX | Compute temperature from specific enthalpy and mass fraction |
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;
Type | Name | Default | Description |
---|---|---|---|
Advanced | |||
Boolean | preferredMediumStates | false | = true if StateSelect.prefer shall be used for the independent property variables of the medium |
redeclare replaceable model extends BaseProperties( T(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default), p(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default), Xi(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default)) /* p, T, X = X[Water] are used as preferred states, since only then all other quantities can be computed in a recursive sequence. If other variables are selected as states, static state selection is no longer possible and non-linear algebraic equations occur. */ MassFraction x_water "mass of total water/mass of dry air"; Real phi "relative humidity"; protected constant SI.MolarMass[2] MMX = {steam.MM,dryair.MM} "Molar masses of components"; MassFraction X_liquid "Mass fraction of liquid water. Need to be zero."; MassFraction X_steam "Mass fraction of steam water"; MassFraction X_air "Mass fraction of air"; MassFraction X_sat "Steam water mass fraction of saturation boundary in kg_water/kg_moistair"; MassFraction x_sat "Steam water mass content of saturation boundary in kg_water/kg_dryair"; AbsolutePressure p_steam_sat "Partial saturation pressure of steam"; equation assert(T >= 200.0 and T <= 423.15, " Temperature T is not in the allowed range 200.0 K <= (T =" + String(T) + " K) <= 423.15 K required from medium model \"" + mediumName + "\"."); MM = 1/(Xi[Water]/MMX[Water]+(1.0-Xi[Water])/MMX[Air]); p_steam_sat = min(saturationPressure(T),0.999*p); X_sat = min(p_steam_sat * k_mair/max(100*Modelica.Constants.eps, p - p_steam_sat)*(1 - Xi[Water]), 1.0) "Water content at saturation with respect to actual water content"; X_liquid = max(Xi[Water] - X_sat, 0.0); assert(Xi[Water] < X_sat, "The medium model '" + mediumName + "' must not be saturated.\n" + "To model a saturated medium, use 'Buildings.Media.PerfectGases.MoistAir' instead of this medium."); X_steam = Xi[Water]; X_air = 1-Xi[Water]; h = specificEnthalpy_pTX(p,T,Xi); R = dryair.R*(1 - X_steam) + steam.R*X_steam; // u = h - R*T; d = p/(R*T); /* Note, u and d are computed under the assumption that the volume of the liquid water is neglible with respect to the volume of air and of steam */ state.p = p; state.T = T; state.X = X; // this x_steam is water load / dry air!!!!!!!!!!! x_sat = k_mair*p_steam_sat/max(100*Modelica.Constants.eps,p - p_steam_sat); x_water = Xi[Water]/max(X_air,100*Modelica.Constants.eps); phi = p/p_steam_sat*Xi[Water]/(Xi[Water] + k_mair*X_air); end BaseProperties;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | shermodynamic state |
Type | Name | Description |
---|---|---|
MassFraction | X_sat | steam mass fraction of sat. boundary [kg/kg] |
function Xsaturation = Buildings.Media.PerfectGases.MoistAir.Xsaturation "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | Mass fractions [kg/kg] |
Type | Name | Description |
---|---|---|
ThermodynamicState | state |
redeclare function setState_pTX "Thermodynamic state as function of p, T and composition X" extends Buildings.Media.PerfectGases.MoistAir.setState_pTX; end setState_pTX;
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
MassFraction | X[:] | Mass fractions [kg/kg] |
Type | Name | Description |
---|---|---|
ThermodynamicState | state |
redeclare function setState_phX "Thermodynamic state as function of p, h and composition X" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input SpecificEnthalpy h "Specific enthalpy"; input MassFraction X[:] "Mass fractions"; output ThermodynamicState state; algorithm state := if size(X,1) == nX then ThermodynamicState(p=p,T=T_phX(p,h,X),X=X) else ThermodynamicState(p=p,T=T_phX(p,h,X), X=cat(1,X,{1-sum(X)})); end setState_phX;
Type | Name | Default | Description |
---|---|---|---|
Density | d | density [kg/m3] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | Mass fractions [kg/kg] |
Type | Name | Description |
---|---|---|
ThermodynamicState | state |
redeclare function setState_dTX "Thermodynamic state as function of d, T and composition X" extends Buildings.Media.PerfectGases.MoistAir.setState_dTX; end setState_dTX;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | thermodynamic state |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | R | mixture gas constant [J/(kg.K)] |
redeclare function gasConstant "Gas constant" extends Buildings.Media.PerfectGases.MoistAir.gasConstant; end gasConstant;
Type | Name | Default | Description |
---|---|---|---|
Temperature | Tsat | saturation temperature [K] |
Type | Name | Description |
---|---|---|
AbsolutePressure | psat | saturation pressure [Pa] |
function saturationPressureLiquid = Buildings.Media.PerfectGases.MoistAir.saturationPressureLiquid "Saturation curve valid for 273.16 <= T <= 373.16. Outside of these limits a (less accurate) result is returned";
Type | Name | Default | Description |
---|---|---|---|
Temperature | Tsat | sublimation temperature [K] |
Type | Name | Description |
---|---|---|
AbsolutePressure | psat | sublimation pressure [Pa] |
function sublimationPressureIce = Buildings.Media.PerfectGases.MoistAir.sublimationPressureIce "Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned";
Type | Name | Default | Description |
---|---|---|---|
Temperature | Tsat | saturation temperature [K] |
Type | Name | Description |
---|---|---|
AbsolutePressure | psat | saturation pressure [Pa] |
redeclare function extends saturationPressure "Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)" algorithm psat := Buildings.Utilities.Math.spliceFunction(saturationPressureLiquid(Tsat),sublimationPressureIce(Tsat),Tsat-273.16,1.0); end saturationPressure;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
AbsolutePressure | p | Pressure [Pa] |
redeclare function pressure "Gas pressure" extends Buildings.Media.PerfectGases.MoistAir.pressure; end pressure;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
Temperature | T | Temperature [K] |
redeclare function temperature "Gas temperature" extends Buildings.Media.PerfectGases.MoistAir.temperature; end temperature;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
Density | d | Density [kg/m3] |
redeclare function density "Gas density" extends Buildings.Media.PerfectGases.MoistAir.density; end density;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificEntropy | s | Specific entropy [J/(kg.K)] |
redeclare function specificEntropy "Specific entropy (liquid part neglected, mixing entropy included)" extends Buildings.Media.PerfectGases.MoistAir.specificEntropy; end specificEntropy;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | r0 | vaporization enthalpy [J/kg] |
redeclare function extends enthalpyOfVaporization "Enthalpy of vaporization of water" algorithm r0 := 2501014.5; end enthalpyOfVaporization;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | [K] |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | cp_fl | [J/(kg.K)] |
function HeatCapacityOfWater "Specific heat capacity of water (liquid only) which is constant" extends Modelica.Icons.Function; input Temperature T; output SpecificHeatCapacity cp_fl; algorithm cp_fl := 4186; end HeatCapacityOfWater;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | liquid enthalpy [J/kg] |
redeclare replaceable function extends enthalpyOfLiquid "Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature" annotation(derivative=der_enthalpyOfLiquid); algorithm h := (T - 273.15)*4186; end enthalpyOfLiquid;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] | |
Temperature | der_T | temperature derivative [K] |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | der_h | derivative of liquid enthalpy [J/(kg.K)] |
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;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | steam enthalpy [J/kg] |
redeclare function enthalpyOfCondensingGas "Enthalpy of steam per unit mass of steam" annotation(derivative=der_enthalpyOfCondensingGas); extends Modelica.Icons.Function; input Temperature T "temperature"; output SpecificEnthalpy h "steam enthalpy"; algorithm h := (T-273.15) * steam.cp + enthalpyOfVaporization(T); end enthalpyOfCondensingGas;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] | |
Temperature | der_T | temperature derivative [K] |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | der_h | derivative of steam enthalpy [J/(kg.K)] |
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;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] | |
MassFraction | X[:] | vector of mass fractions [kg/kg] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | liquid enthalpy [J/kg] |
redeclare replaceable function extends enthalpyOfGas "Enthalpy of gas mixture per unit mass of gas mixture" algorithm h := enthalpyOfCondensingGas(T)*X[Water] + enthalpyOfDryAir(T)*(1.0-X[Water]); end enthalpyOfGas;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | dry air enthalpy [J/kg] |
replaceable function enthalpyOfDryAir "Enthalpy of dry air per unit mass of dry air" annotation(derivative=der_enthalpyOfDryAir); extends Modelica.Icons.Function; input Temperature T "temperature"; output SpecificEnthalpy h "dry air enthalpy"; algorithm h := (T - 273.15)*dryair.cp; end enthalpyOfDryAir;
Type | Name | Default | Description |
---|---|---|---|
Temperature | T | temperature [K] | |
Temperature | der_T | temperature derivative [K] |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | der_h | derivative of dry air enthalpy [J/(kg.K)] |
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;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | cp | Specific heat capacity at constant pressure [J/(kg.K)] |
redeclare replaceable function extends specificHeatCapacityCp "Specific heat capacity of gas mixture at constant pressure" algorithm cp := dryair.cp*(1-state.X[Water]) +steam.cp*state.X[Water]; end specificHeatCapacityCp;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificHeatCapacity | cv | Specific heat capacity at constant volume [J/(kg.K)] |
redeclare replaceable function extends specificHeatCapacityCv "Specific heat capacity of gas mixture at constant volume" algorithm cv:= dryair.cv*(1-state.X[Water]) +steam.cv*state.X[Water]; end specificHeatCapacityCv;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
DynamicViscosity | eta | Dynamic viscosity [Pa.s] |
redeclare function extends dynamicViscosity "dynamic viscosity of dry air" algorithm eta := 1.85E-5; end dynamicViscosity;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
ThermalConductivity | lambda | Thermal conductivity [W/(m.K)] |
redeclare function extends thermalConductivity "Thermal conductivity of dry air as a polynomial in the temperature" algorithm lambda := Polynomials_Temp.evaluate({(-4.8737307422969E-008), 7.67803133753502E-005, 0.0241814385504202}, Cv.to_degC(state.T)); end thermalConductivity;
Type | Name | Default | Description |
---|---|---|---|
Pressure | p | Pressure [Pa] | |
Temperature | T | Temperature [K] | |
MassFraction | X[nX] | Mass fractions of moist air [1] |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Specific enthalpy at p, T, X [J/kg] |
function h_pTX "Compute specific enthalpy from pressure, temperature and mass fraction" extends Modelica.Icons.Function; input SI.Pressure p "Pressure"; input SI.Temperature T "Temperature"; input SI.MassFraction X[nX] "Mass fractions of moist air"; output SI.SpecificEnthalpy h "Specific enthalpy at p, T, X"; protected SI.AbsolutePressure p_steam_sat "Partial saturation pressure of steam"; SI.MassFraction x_sat "steam water mass fraction of saturation boundary"; SI.SpecificEnthalpy hDryAir "Enthalpy of dry air"; algorithm p_steam_sat :=saturationPressure(T); x_sat :=k_mair*p_steam_sat/(p - p_steam_sat); assert(X[Water] < x_sat/(1 + x_sat), "The medium model '" + mediumName + "' must not be saturated.\n" + "To model a saturated medium, use 'Buildings.Media.PerfectGases.MoistAir' instead of this medium."); hDryAir := (T - 273.15)*dryair.cp; h := hDryAir * (1 - X[Water]) + ((T-273.15) * steam.cp + 2501014.5) * X[Water]; end h_pTX;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificEnthalpy | h | Specific enthalpy [J/kg] |
redeclare function extends specificEnthalpy "Specific enthalpy" algorithm h := h_pTX(state.p, state.T, state.X); end specificEnthalpy;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificEnergy | u | Specific internal energy [J/kg] |
redeclare function extends specificInternalEnergy "Specific internal energy" extends Modelica.Icons.Function; algorithm u := h_pTX(state.p,state.T,state.X) - gasConstant(state)*state.T; end specificInternalEnergy;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificEnergy | g | Specific Gibbs energy [J/kg] |
redeclare function extends specificGibbsEnergy "Specific Gibbs energy" extends Modelica.Icons.Function; algorithm g := h_pTX(state.p,state.T,state.X) - state.T*specificEntropy(state); end specificGibbsEnergy;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state |
Type | Name | Description |
---|---|---|
SpecificEnergy | f | Specific Helmholtz energy [J/kg] |
redeclare function extends specificHelmholtzEnergy "Specific Helmholtz energy" extends Modelica.Icons.Function; algorithm f := h_pTX(state.p,state.T,state.X) - gasConstant(state)*state.T - state.T*specificEntropy(state); end specificHelmholtzEnergy;
Type | Name | Default | Description |
---|---|---|---|
AbsolutePressure | p | Pressure [Pa] | |
SpecificEnthalpy | h | specific enthalpy [J/kg] | |
MassFraction | X[:] | mass fractions of composition [kg/kg] |
Type | Name | Description |
---|---|---|
Temperature | T | temperature [K] |
function T_phX "Compute temperature from specific enthalpy and mass fraction" input AbsolutePressure p "Pressure"; input SpecificEnthalpy h "specific enthalpy"; input MassFraction[:] X "mass fractions of composition"; output Temperature T "temperature"; protected SI.AbsolutePressure p_steam_sat "Partial saturation pressure of steam"; SI.MassFraction x_sat "steam water mass fraction of saturation boundary"; algorithm T := (h - 2501014.5 * X[Water])/((1 - X[Water])*dryair.cp + X[Water] * steam.cp); // check for saturation p_steam_sat :=saturationPressure(T); x_sat :=k_mair*p_steam_sat/(p - p_steam_sat); assert(X[Water] < x_sat/(1 + x_sat), "The medium model '" + mediumName + "' must not be saturated.\n" + "To model a saturated medium, use 'Buildings.Media.PerfectGases.MoistAir' instead of this medium."); end T_phX;