LBL logo

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated

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

Information

This is a medium model that is identical to Buildings.Media.GasesPTDecoupled.MoistAir, but in this model, the air must not be saturated. If the air is saturated, use the medium model Buildings.Media.GasesPTDecoupled.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 thermo-fluid systems.

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

Extends from Modelica.Media.Interfaces.PartialCondensingGases (Base class for mixtures of condensing and non-condensing gases).

Package Content

NameDescription
Water=1Index of water (in substanceNames, massFractions X, etc.)
Air=2Index of air (in substanceNames, massFractions X, etc.)
k_mair=steam.MM/dryair.MMratio of molar weights
dryair=Buildings.Media.PerfectGases.Common.SingleGasData.Air 
steam=Buildings.Media.PerfectGases.Common.SingleGasData.H2O 
pStp=101325Pressure for which dStp is defined
dStp=1.2Fluid density at pressure pStp
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.ThermodynamicState ThermodynamicState ThermodynamicState record for moist air
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.BaseProperties BaseProperties  
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.Xsaturation Xsaturation Steam water mass fraction of saturation boundary in kg_water/kg_moistair
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_pTX setState_pTX Thermodynamic state as function of p, T and composition X
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_phX setState_phX Thermodynamic state as function of p, h and composition X
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_dTX setState_dTX Thermodynamic state as function of d, T and composition X
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.gasConstant gasConstant Gas constant (computation neglects liquid fraction)
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid saturationPressureLiquid Return saturation pressure of water as a function of temperature T in the range of 273.16 to 373.16 K
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid_der saturationPressureLiquid_der Time derivative of saturationPressureLiquid
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.sublimationPressureIce sublimationPressureIce Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressure saturationPressure Saturation curve valid for 223.16 <= T <= 373.16 (and slightly outside with less accuracy)
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.pressure pressure Gas pressure
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.temperature temperature Gas temperature
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.density density Gas density
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEntropy specificEntropy Specific entropy (liquid part neglected, mixing entropy included)
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfVaporization enthalpyOfVaporization Enthalpy of vaporization of water
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.HeatCapacityOfWater HeatCapacityOfWater Specific heat capacity of water (liquid only) which is constant
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfLiquid enthalpyOfLiquid Enthalpy of liquid (per unit mass of liquid) which is linear in the temperature
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfLiquid der_enthalpyOfLiquid Temperature derivative of enthalpy of liquid per unit mass of liquid
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfCondensingGas enthalpyOfCondensingGas Enthalpy of steam per unit mass of steam
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfCondensingGas der_enthalpyOfCondensingGas Derivative of enthalpy of steam per unit mass of steam
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfGas enthalpyOfGas Enthalpy of gas mixture per unit mass of gas mixture
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfDryAir enthalpyOfDryAir Enthalpy of dry air per unit mass of dry air
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfDryAir der_enthalpyOfDryAir Derivative of enthalpy of dry air per unit mass of dry air
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHeatCapacityCp specificHeatCapacityCp Specific heat capacity of gas mixture at constant pressure
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_specificHeatCapacityCp der_specificHeatCapacityCp Derivative of specific heat capacity of gas mixture at constant pressure
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHeatCapacityCv specificHeatCapacityCv Specific heat capacity of gas mixture at constant volume
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_specificHeatCapacityCv der_specificHeatCapacityCv Derivative of specific heat capacity of gas mixture at constant volume
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.dynamicViscosity dynamicViscosity dynamic viscosity of dry air
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.thermalConductivity thermalConductivity Thermal conductivity of dry air as a polynomial in the temperature
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEnthalpy specificEnthalpy Specific enthalpy
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificInternalEnergy specificInternalEnergy Specific internal energy
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificGibbsEnergy specificGibbsEnergy Specific Gibbs energy
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHelmholtzEnergy specificHelmholtzEnergy Specific Helmholtz energy
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX h_pTX Compute specific enthalpy from pressure, temperature and mass fraction
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.T_phX T_phX Compute temperature from specific enthalpy and mass fraction
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfNonCondensingGas enthalpyOfNonCondensingGas Enthalpy of non-condensing gas per unit mass
Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfNonCondensingGas der_enthalpyOfNonCondensingGas Derivative of enthalpy of non-condensing gas per unit mass
Inherited
Modelica.Media.Interfaces.PartialMixtureMedium.FluidConstants FluidConstants extended fluid constants
fluidConstantsconstant 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
ThermoStatesEnumeration type for independent variables
mediumName="unusablePartialMedium"Name of the medium
substanceNames={mediumName}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= 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=101325Reference pressure of Medium: default 1 atmosphere
reference_T=298.15Reference temperature of Medium: default 25 deg Celsius
reference_X=fill(1/nX, nX)Default mass fractions of medium
p_default=101325Default 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_XDefault value for mass fractions of medium (for initialization)
nS=size(substanceNames, 1)Number of substances
nX=nSNumber of mass fractions
nXi=if fixedX then 0 else if reducedX then nS - 1 else nSNumber 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.setState_psX setState_psX Return thermodynamic state as function of p, s and composition X or Xi
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.isentropicEnthalpy isentropicEnthalpy Return isentropic enthalpy
Modelica.Media.Interfaces.PartialMedium.velocityOfSound velocityOfSound Return velocity of sound
Modelica.Media.Interfaces.PartialMedium.isobaricExpansionCoefficient isobaricExpansionCoefficient Return overall the isobaric expansion coefficient beta
Modelica.Media.Interfaces.PartialMedium.beta beta alias for isobaricExpansionCoefficient for user convenience
Modelica.Media.Interfaces.PartialMedium.isothermalCompressibility isothermalCompressibility Return overall the isothermal compressibility factor
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.density_derp_T density_derp_T Return density derivative w.r.t. pressure at const temperature
Modelica.Media.Interfaces.PartialMedium.density_derT_p density_derT_p Return density derivative w.r.t. temperature at constant pressure
Modelica.Media.Interfaces.PartialMedium.density_derX density_derX Return density derivative w.r.t. mass fraction
Modelica.Media.Interfaces.PartialMedium.molarMass molarMass Return the molar mass of the medium
Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_pTX specificEnthalpy_pTX Return specific enthalpy from p, T, and X or Xi
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.temperature_phX temperature_phX Return temperature from p, h, 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
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
MassFlowRate Type for mass 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 resect to pressure with medium specific attributes
DerDensityByEnthalpy Type for partial derivative of density with resect to enthalpy with medium specific attributes
DerEnthalpyByPressure Type for partial derivative of enthalpy with resect to pressure with medium specific attributes
DerDensityByTemperature Type for partial derivative of density with resect to temperature with medium specific attributes
Modelica.Media.Interfaces.PartialMedium.Choices Choices Types, constants to define menu choices

Types and constants

  constant Integer Water=1 
  "Index of water (in substanceNames, massFractions X, etc.)";

  constant Integer Air=2 
  "Index of air (in substanceNames, massFractions X, etc.)";

  constant Real k_mair =  steam.MM/dryair.MM "ratio of molar weights";

  constant Buildings.Media.PerfectGases.Common.DataRecord dryair=
        Buildings.Media.PerfectGases.Common.SingleGasData.Air;

  constant Buildings.Media.PerfectGases.Common.DataRecord steam=
        Buildings.Media.PerfectGases.Common.SingleGasData.H2O;

  constant AbsolutePressure pStp = 101325 "Pressure for which dStp is defined";

  constant Density dStp = 1.2 "Fluid density at pressure pStp";


Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.ThermodynamicState Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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.GasesPTDecoupled.MoistAirUnsaturated.BaseProperties Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.BaseProperties

Information

Extends from (Base properties (p, d, T, h, u, R, MM and, if applicable, X and Xi) of a medium).

Parameters

TypeNameDefaultDescription
BooleanstandardOrderComponentstrueif true, and reducedX = true, the last element of X will be computed from the other ones
Advanced
BooleanpreferredMediumStatesfalse= true if StateSelect.prefer shall be used for the independent property variables of the medium

Modelica definition

redeclare replaceable model extends BaseProperties(
  T(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default),
  p(stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default),
  Xi(each stateSelect=if preferredMediumStates then StateSelect.prefer else StateSelect.default),
  final standardOrderComponents=true)

  /* p, T, X = X[Water] are used as preferred states, since only then all
     other quantities can be computed in a recursive sequence. 
     If other variables are selected as states, static state selection
     is no longer possible and non-linear algebraic equations occur.
      */
  MassFraction x_water "Mass of total water/mass of dry air";
  Real phi "Relative humidity";

protected 
  constant SI.MolarMass[2] MMX = {steam.MM,dryair.MM} 
    "Molar masses of components";

  //    MassFraction X_liquid "Mass fraction of liquid water";
  MassFraction X_steam "Mass fraction of steam water";
  MassFraction X_air "Mass fraction of air";
  MassFraction X_sat 
    "Steam water mass fraction of saturation boundary in kg_water/kg_moistair";
  MassFraction x_sat 
    "Steam water mass content of saturation boundary in kg_water/kg_dryair";
  AbsolutePressure p_steam_sat "Partial saturation pressure of steam";

equation 
  assert(T >= 200.0 and T <= 423.15, "
Temperature T is not in the allowed range
200.0 K <= (T ="
             + String(T) + " K) <= 423.15 K
required from medium model \""   + mediumName + "\".");

/*
    assert(Xi[Water] <= X_sat, "The medium model '" + mediumName + "' must not be saturated.\n"
     + "To model a saturated medium, use 'Buildings.Media.GasesPTDecoupled.MoistAir' instead of this medium.\n"
     + " T         = " + String(T) + "\n"
     + " X_sat     = " + String(X_sat) + "\n"
     + " Xi[Water] = " + String(Xi[Water]) + "\n"
     + " phi       = " + String(phi) + "\n"
     + " p         = " + String(p));
  */

  MM = 1/(Xi[Water]/MMX[Water]+(1.0-Xi[Water])/MMX[Air]);

  p_steam_sat = min(saturationPressure(T),0.999*p);
  X_sat = min(p_steam_sat * k_mair/max(100*Modelica.Constants.eps, p - p_steam_sat)*(1 - Xi[Water]), 1.0) 
    "Water content at saturation with respect to actual water content";
  //    X_liquid = max(Xi[Water] - X_sat, 0.0);
  //    X_steam  = Xi[Water]-X_liquid;

  X_steam  = Xi[Water]; // There is no liquid in this medium model
  X_air    = 1-Xi[Water];

  h = specificEnthalpy_pTX(p,T,Xi);
  R = dryair.R*(1 - Xi[Water]) + steam.R*Xi[Water];

  // 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;

  //    d = p/(R*T);
  d/dStp = p/pStp;

  /* Note, u and d are computed under the assumption that the volume of the liquid
         water is neglible with respect to the volume of air and of steam
      */
  state.p = p;
  state.T = T;
  state.X = X;

  // this x_steam is water load / dry air!!!!!!!!!!!
  x_sat    = k_mair*p_steam_sat/max(100*Modelica.Constants.eps,p - p_steam_sat);
  x_water = Xi[Water]/max(X_air,100*Modelica.Constants.eps);
  phi = p/p_steam_sat*Xi[Water]/(Xi[Water] + k_mair*X_air);
end BaseProperties;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.Xsaturation

Steam water mass fraction of saturation boundary in kg_water/kg_moistair

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_pTX Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_pTX

Thermodynamic state as function of p, T and composition X

Information

Extends from Buildings.Media.PerfectGases.MoistAir.setState_pTX (Thermodynamic state as function of p, T and composition X).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_phX Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_phX

Thermodynamic state as function of p, h and composition X

Information

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

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

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_phX 
  "Thermodynamic state as function of p, h and composition X"
extends Modelica.Icons.Function;
input AbsolutePressure p "Pressure";
input SpecificEnthalpy h "Specific enthalpy";
input MassFraction X[:] "Mass fractions";
output ThermodynamicState state;
algorithm 
state := if size(X,1) == nX then 
       ThermodynamicState(p=p,T=T_phX(p,h,X),X=X) else 
      ThermodynamicState(p=p,T=T_phX(p,h,cat(1,X,{1-sum(X)})), X=cat(1,X,{1-sum(X)}));
  //    ThermodynamicState(p=p,T=T_phX(p,h,X), X=cat(1,X,{1-sum(X)}));
end setState_phX;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_dTX Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.setState_dTX

Thermodynamic state as function of d, T and composition X

Information

Extends from Buildings.Media.PerfectGases.MoistAir.setState_dTX (Thermodynamic state as function of d, T and composition X).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.gasConstant Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.gasConstant

Gas constant (computation neglects liquid fraction)

Information

Extends from Buildings.Media.PerfectGases.MoistAir.gasConstant (Gas constant (computation neglects liquid fraction)).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state

Outputs

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

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid

Return saturation pressure of water as a function of temperature T in the range of 273.16 to 373.16 K

Information

Saturation pressure of water above the triple point temperature is computed from temperature. It's range of validity is between 273.16 and 373.16 K. Outside these limits a less accurate result is returned.

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

Inputs

TypeNameDefaultDescription
TemperatureTsat saturation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatsaturation pressure [Pa]

Modelica definition

function saturationPressureLiquid 
  "Return saturation pressure of water as a function of temperature T in the range of 273.16 to 373.16 K"
  annotation(derivative=Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid_der);

  extends Modelica.Icons.Function;
  input SI.Temperature Tsat "saturation temperature";
  output SI.AbsolutePressure psat "saturation pressure";
algorithm 
  psat := 611.657*Modelica.Math.exp(17.2799 - 4102.99/(Tsat - 35.719));
end saturationPressureLiquid;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid_der Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressureLiquid_der

Time derivative of saturationPressureLiquid

Information

Derivative function of saturationPressureLiquid

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

Inputs

TypeNameDefaultDescription
TemperatureTsat Saturation temperature [K]
RealdTsat Saturation temperature derivative [K/s]

Outputs

TypeNameDescription
Realpsat_derSaturation pressure [Pa/s]

Modelica definition

function saturationPressureLiquid_der 
  "Time derivative of saturationPressureLiquid"

  extends Modelica.Icons.Function;
  input SI.Temperature Tsat "Saturation temperature";
  input Real dTsat(unit="K/s") "Saturation temperature derivative";
  output Real psat_der(unit="Pa/s") "Saturation pressure";
algorithm 
  psat_der:=611.657*Modelica.Math.exp(17.2799 - 4102.99/(Tsat - 35.719))*4102.99*dTsat/(Tsat - 35.719)/(Tsat - 35.719);

end saturationPressureLiquid_der;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.sublimationPressureIce Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.sublimationPressureIce

Saturation curve valid for 223.16 <= T <= 273.16. Outside of these limits a (less accurate) result is returned

Inputs

TypeNameDefaultDescription
TemperatureTsat sublimation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatsublimation pressure [Pa]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.saturationPressure Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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.Math.Functions.spliceFunction(
                                                  saturationPressureLiquid(Tsat),sublimationPressureIce(Tsat),Tsat-273.16,1.0);
end saturationPressure;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.pressure Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.pressure

Gas pressure

Information

Extends from Buildings.Media.PerfectGases.MoistAir.pressure (Gas pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.temperature Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.temperature

Gas temperature

Information

Extends from Buildings.Media.PerfectGases.MoistAir.temperature (Gas temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.density Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.density

Gas density

Information

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.GasesPTDecoupled.MoistAirUnsaturated.specificEntropy Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEntropy

Specific entropy (liquid part neglected, mixing entropy included)

Information

Extends from Buildings.Media.PerfectGases.MoistAir.specificEntropy (Specific entropy (liquid part neglected, mixing entropy included)).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfVaporization Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 := 2501014.5;
end enthalpyOfVaporization;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.HeatCapacityOfWater Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.HeatCapacityOfWater

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

Information

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

Inputs

TypeNameDefaultDescription
TemperatureT [K]

Outputs

TypeNameDescription
SpecificHeatCapacitycp_fl[J/(kg.K)]

Modelica definition

function HeatCapacityOfWater 
  "Specific heat capacity of water (liquid only) which is constant"
  extends Modelica.Icons.Function;
  input Temperature T;
  output SpecificHeatCapacity cp_fl;
algorithm 
  cp_fl := 4186;
end HeatCapacityOfWater;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfLiquid Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 - 273.15)*4186;
end enthalpyOfLiquid;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfLiquid Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 := 4186*der_T;
end der_enthalpyOfLiquid;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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-273.15) * steam.cp + Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfVaporization(T);
end enthalpyOfCondensingGas;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 [K/s]

Outputs

TypeNameDescription
Realder_hderivative of steam enthalpy [J/(kg.s)]

Modelica definition

replaceable function der_enthalpyOfCondensingGas 
  "Derivative of enthalpy of steam per unit mass of steam"
  extends Modelica.Icons.Function;
  input Temperature T "temperature";
  input Real der_T(unit="K/s") "temperature derivative";
  output Real der_h(unit="J/(kg.s)") "derivative of steam enthalpy";
algorithm 
  der_h := steam.cp*der_T;
end der_enthalpyOfCondensingGas;

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfDryAir Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 - 273.15)*dryair.cp;
end enthalpyOfDryAir;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfDryAir Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 [K/s]

Outputs

TypeNameDescription
Realder_hderivative of dry air enthalpy [J/(kg.s)]

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(unit="K/s") "temperature derivative";
  output Real der_h(unit="J/(kg.s)") "derivative of dry air enthalpy";
algorithm 
  der_h := dryair.cp*der_T;
end der_enthalpyOfDryAir;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHeatCapacityCp Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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.GasesPTDecoupled.MoistAirUnsaturated.der_specificHeatCapacityCp

Derivative of specific heat capacity of gas mixture at constant pressure

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  
ThermodynamicStateder_state  

Outputs

TypeNameDescription
Realder_cp[J/(kg.K.s)]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHeatCapacityCv Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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.GasesPTDecoupled.MoistAirUnsaturated.der_specificHeatCapacityCv

Derivative of specific heat capacity of gas mixture at constant volume

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  
ThermodynamicStateder_state  

Outputs

TypeNameDescription
Realder_cv[J/(kg.K.s)]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.dynamicViscosity Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.dynamicViscosity

dynamic viscosity of dry air

Information

Extends from (Return dynamic viscosity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

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

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEnthalpy Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEnthalpy

Specific enthalpy

Information

Extends from (Return specific enthalpy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

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

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificInternalEnergy Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.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 := Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX(state.p,state.T,state.X) - pStp/dStp;
end specificInternalEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificGibbsEnergy Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificGibbsEnergy

Specific Gibbs energy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

redeclare function extends specificGibbsEnergy 
  "Specific Gibbs energy"
  extends Modelica.Icons.Function;
algorithm 
  g := Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX(state.p,state.T,state.X) - state.T*specificEntropy(state);
end specificGibbsEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHelmholtzEnergy Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificHelmholtzEnergy

Specific Helmholtz energy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

redeclare function extends specificHelmholtzEnergy 
  "Specific Helmholtz energy"
  extends Modelica.Icons.Function;
algorithm 
  f := Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX(state.p,state.T,state.X)
         - gasConstant(state)*state.T
         - state.T*Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.specificEntropy(state);
end specificHelmholtzEnergy;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.h_pTX

Compute specific enthalpy from pressure, temperature and mass fraction

Information

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

Inputs

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

Outputs

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

Modelica definition

function h_pTX 
  "Compute specific enthalpy from pressure, temperature and mass fraction"
  extends Modelica.Icons.Function;

  input SI.Pressure p "Pressure";
  input SI.Temperature T "Temperature";
  input SI.MassFraction X[nX] "Mass fractions of moist air";
  output SI.SpecificEnthalpy h "Specific enthalpy at p, T, X";
protected 
  SI.AbsolutePressure p_steam_sat "Partial saturation pressure of steam";
  SI.MassFraction x_sat "steam water mass fraction of saturation boundary";
  SI.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]-0.001 < x_sat/(1 + x_sat), "The medium model '" + mediumName + "' must not be saturated.\n"
     + "To model a saturated medium, use 'Buildings.Media.GasesPTDecoupled.MoistAir' instead of this medium.\n"
     + " T         = " + String(T) + "\n"
     + " x_sat     = " + String(x_sat) + "\n"
     + " X[Water] = "  + String(X[Water]) + "\n"
     + " phi       = " + String(X[Water]/((x_sat)/(1+x_sat))) + "\n"
     + " p         = " + String(p));
     */
  h := (T - 273.15)*dryair.cp * (1 - X[Water]) + ((T-273.15) * steam.cp + 2501014.5) * X[Water];
end h_pTX;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.T_phX Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.T_phX

Compute temperature from specific enthalpy and mass fraction

Information

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

function T_phX 
  "Compute temperature from specific enthalpy and mass fraction"
  extends Modelica.Icons.Function;

  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "specific enthalpy";
  input MassFraction[:] X "mass fractions of composition";
  output Temperature T "temperature";
protected 
  SI.AbsolutePressure p_steam_sat "Partial saturation pressure of steam";
  SI.MassFraction x_sat "steam water mass fraction of saturation boundary";
algorithm 
  T := 273.15 + (h-2501014.5 * X[Water])/(dryair.cp * (1 - X[Water])+steam.cp*X[Water]);
  // Check for saturation
  p_steam_sat :=saturationPressure(T);
  x_sat    :=k_mair*p_steam_sat/(p - p_steam_sat);
/*  
  assert(X[Water]-0.001 < x_sat/(1 + x_sat), "The medium model '" + mediumName + "' must not be saturated.\n"
     + "To model a saturated medium, use 'Buildings.Media.GasesPTDecoupled.MoistAir' instead of this medium.\n"
     + " T         = " + String(T) + "\n"
     + " x_sat     = " + String(x_sat) + "\n"
     + " X[Water] = " + String(X[Water]) + "\n"
     + " phi       = " + String(X[Water]/((x_sat)/(1+x_sat))) + "\n"
     + " p         = " + String(p));
*/
end T_phX;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfNonCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.enthalpyOfNonCondensingGas

Enthalpy of non-condensing gas per unit mass

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"
  annotation(derivative=der_enthalpyOfNonCondensingGas);
  extends Modelica.Icons.Function;

  input Temperature T "temperature";
  output SpecificEnthalpy h "enthalpy";
algorithm 
  h := enthalpyOfDryAir(T);
end enthalpyOfNonCondensingGas;

Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfNonCondensingGas Buildings.Media.GasesPTDecoupled.MoistAirUnsaturated.der_enthalpyOfNonCondensingGas

Derivative of enthalpy of non-condensing gas per unit mass

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"
  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;

Automatically generated Thu Oct 24 15:10:44 2013.