LBL logo

Buildings.Media.GasesPTDecoupled.SimpleAir

Package with dry air model that decouples pressure and temperature

Information

This medium model is identical to Modelica.Media.Air.SimpleAir, except the equation d = p/(R*T) has been replaced with d/dStp = p/pStp where pStd and dStp are constants for a reference temperature and density.

This new formulation often leads to smaller systems of nonlinear equations because pressure and temperature are decoupled, at the expense of accuracy.

As in Modelica.Media.Air.SimpleAir, the specific enthalpy h and specific internal energy u are only a function of temperature T and all other provided medium quantities are constant.

Extends from Buildings.Media.Interfaces.PartialSimpleIdealGasMedium (Medium model of Ideal gas with constant cp and cv. All other quantities, e.g., transport properties, are constant.).

Package Content

NameDescription
fluidConstants=FluidConstants(iupacName={"simple air"}, casRegistryNumber={"not a real substance"}, chemicalFormula={"N2, O2"}, structureFormula={"N2, O2"}, molarMass=Modelica.Media.IdealGases.Common.SingleGasesData.N2.MM)constant data for the fluid
pStp=101325Pressure for which dStp is defined
dStp=1.2Fluid density at pressure pStp
Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties BaseProperties Basic medium properties
Buildings.Media.GasesPTDecoupled.SimpleAir.setState_dTX setState_dTX Return thermodynamic state from d, T, and X or Xi
Buildings.Media.GasesPTDecoupled.SimpleAir.density density return density of ideal gas
Buildings.Media.GasesPTDecoupled.SimpleAir.specificInternalEnergy specificInternalEnergy Return specific internal energy
Buildings.Media.GasesPTDecoupled.SimpleAir.specificEntropy specificEntropy Return specific entropy
Buildings.Media.GasesPTDecoupled.SimpleAir.enthalpyOfCondensingGas enthalpyOfCondensingGas Enthalpy of steam per unit mass of steam
Buildings.Media.GasesPTDecoupled.SimpleAir.saturationPressure saturationPressure Return saturation pressure of condensing fluid
Inherited
cp_constConstant specific heat capacity at constant pressure
cv_const=cp_const - R_gasConstant specific heat capacity at constant volume
R_gasmedium specific gas constant
MM_constMolar mass
eta_constConstant dynamic viscosity
lambda_constConstant thermal conductivity
T_minMinimum temperature valid for medium model
T_maxMaximum temperature valid for medium model
T0=reference_TZero enthalpy temperature
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.ThermodynamicState ThermodynamicState Thermodynamic state of ideal gas
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.FluidConstants FluidConstants fluid constants
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_pTX setState_pTX Return thermodynamic state from p, T, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_phX setState_phX Return thermodynamic state from p, h, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_psX setState_psX Return thermodynamic state from p, s, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setSmoothState setSmoothState Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.pressure pressure Return pressure of ideal gas
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature temperature Return temperature of ideal gas
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy specificEnthalpy Return specific enthalpy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificGibbsEnergy specificGibbsEnergy Return specific Gibbs energy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHelmholtzEnergy specificHelmholtzEnergy Return specific Helmholtz energy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.dynamicViscosity dynamicViscosity Return dynamic viscosity
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.thermalConductivity thermalConductivity Return thermal conductivity
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCp specificHeatCapacityCp Return specific heat capacity at constant pressure
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCv specificHeatCapacityCv Return specific heat capacity at constant volume
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isentropicExponent isentropicExponent Return isentropic exponent
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.velocityOfSound velocityOfSound Return velocity of sound
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy_pTX specificEnthalpy_pTX Return specific enthalpy from p, T, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature_phX temperature_phX Return temperature from p, h, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_phX density_phX Return density from p, h, and X or Xi
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isentropicEnthalpy isentropicEnthalpy Return isentropic enthalpy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isobaricExpansionCoefficient isobaricExpansionCoefficient Returns overall the isobaric expansion coefficient beta
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isothermalCompressibility isothermalCompressibility Returns overall the isothermal compressibility factor
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_derp_T density_derp_T Returns the partial derivative of density with respect to pressure at constant temperature
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_derT_p density_derT_p Returns the partial derivative of density with respect to temperature at constant pressure
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_derX density_derX Returns the partial derivative of density with respect to mass fractions at constant pressure and temperature
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.molarMass molarMass Returns the molar mass of the medium
Modelica.Media.Interfaces.PartialPureSubstance.setState_pT setState_pT Return thermodynamic state from p and T
Modelica.Media.Interfaces.PartialPureSubstance.setState_ph setState_ph Return thermodynamic state from p and h
Modelica.Media.Interfaces.PartialPureSubstance.setState_ps setState_ps Return thermodynamic state from p and s
Modelica.Media.Interfaces.PartialPureSubstance.setState_dT setState_dT Return thermodynamic state from d and T
Modelica.Media.Interfaces.PartialPureSubstance.density_ph density_ph Return density from p and h
Modelica.Media.Interfaces.PartialPureSubstance.temperature_ph temperature_ph Return temperature from p and h
Modelica.Media.Interfaces.PartialPureSubstance.pressure_dT pressure_dT Return pressure from d and T
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_dT specificEnthalpy_dT Return specific enthalpy from d and T
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_ps specificEnthalpy_ps Return specific enthalpy from p and s
Modelica.Media.Interfaces.PartialPureSubstance.temperature_ps temperature_ps Return temperature from p and s
Modelica.Media.Interfaces.PartialPureSubstance.density_ps density_ps Return density from p and s
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_pT specificEnthalpy_pT Return specific enthalpy from p and T
Modelica.Media.Interfaces.PartialPureSubstance.density_pT density_pT Return density from p and T
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.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.beta beta alias for isobaricExpansionCoefficient for user convenience
Modelica.Media.Interfaces.PartialMedium.kappa kappa alias of isothermalCompressibility for user convenience
Modelica.Media.Interfaces.PartialMedium.density_derp_h density_derp_h Return density derivative w.r.t. pressure at const specific enthalpy
Modelica.Media.Interfaces.PartialMedium.density_derh_p density_derh_p Return density derivative w.r.t. specific enthalpy at constant pressure
Modelica.Media.Interfaces.PartialMedium.specificEntropy_pTX specificEntropy_pTX Return specific enthalpy from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_pTX density_pTX Return density from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.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 FluidConstants[nS] fluidConstants=
    FluidConstants(iupacName={"simple air"},
                   casRegistryNumber={"not a real substance"},
                   chemicalFormula={"N2, O2"},
                   structureFormula={"N2, O2"},
                   molarMass=Modelica.Media.IdealGases.Common.SingleGasesData.N2.MM) 
  "constant data for the fluid";

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

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


Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties

Basic medium properties

Parameters

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

Modelica definition

redeclare replaceable model BaseProperties "Basic medium properties"
   // declarations from Modelica.Media.Interfaces.PartialMedium
   InputAbsolutePressure p "Absolute pressure of medium";
   InputMassFraction[nXi] Xi(start=reference_X[1:nXi]) 
    "Structurally independent mass fractions";
   InputSpecificEnthalpy h "Specific enthalpy of medium";
   Density d "Density of medium";
   Temperature T "Temperature of medium";
   MassFraction[nX] X(start=reference_X) 
    "Mass fractions (= (component mass)/total mass  m_i/m)";
   SpecificInternalEnergy u "Specific internal energy of medium";
   SpecificHeatCapacity R "Gas constant (of mixture if applicable)";
   MolarMass MM "Molar mass (of mixture or single fluid)";
   ThermodynamicState state "thermodynamic state record for optional functions";
   parameter Boolean preferredMediumStates=false 
    "= true if StateSelect.prefer shall be used for the independent property variables of the medium";
   final parameter Boolean standardOrderComponents = true 
    "if true, and reducedX = true, the last element of X will be computed from the other ones";
   SI.Conversions.NonSIunits.Temperature_degC T_degC=
       Modelica.SIunits.Conversions.to_degC(T) 
    "Temperature of medium in [degC]";
   SI.Conversions.NonSIunits.Pressure_bar p_bar=
    Modelica.SIunits.Conversions.to_bar(p) 
    "Absolute pressure of medium in [bar]";

   // Local connector definition, used for equation balancing check
   connector InputAbsolutePressure = input SI.AbsolutePressure 
    "Pressure as input signal connector";
   connector InputSpecificEnthalpy = input SI.SpecificEnthalpy 
    "Specific enthalpy as input signal connector";
   connector InputMassFraction = input SI.MassFraction 
    "Mass fraction as input signal connector";

   // own declarations

equation 
   if standardOrderComponents then
     Xi = X[1:nXi];

       if fixedX then
         X = reference_X;
       end if;
       if reducedX and not fixedX then
         X[nX] = 1 - sum(Xi);
       end if;
       for i in 1:nX loop
         assert(X[i] >= -1.e-5 and X[i] <= 1 + 1.e-5, "Mass fraction X[" +
                String(i) + "] = " + String(X[i]) + "of substance "
                + substanceNames[i] + "\nof medium " + mediumName + " is not in the range 0..1");
       end for;

   end if;

   assert(p >= 0.0, "Pressure (= " + String(p) + " Pa) of medium \"" +
     mediumName + "\" is negative\n(Temperature = " + String(T) + " K)");

   // new medium equations
   h = specificEnthalpy_pTX(p,T,X);
   // 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;
   R = R_gas;
   //    d = p/(R*T);
   d/dStp = p/pStp;

   MM = MM_const;
   state.T = T;
   state.p = p;
end BaseProperties;

Buildings.Media.GasesPTDecoupled.SimpleAir.setState_dTX Buildings.Media.GasesPTDecoupled.SimpleAir.setState_dTX

Return thermodynamic state from d, T, and X or Xi

Information

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

Inputs

TypeNameDefaultDescription
Densityd Density [kg/m3]
TemperatureT Temperature [K]
MassFractionX[:]fill(0, 0)Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_dTX 
  "Return thermodynamic state from d, T, and X or Xi"
   extends Modelica.Icons.Function;
   input Density d "Density";
   input Temperature T "Temperature";
   input MassFraction X[:] = fill(0,0) "Mass fractions";
   output ThermodynamicState state;
algorithm 
   state := ThermodynamicState(p=d/dStp*pStp,T=T);
end setState_dTX;

Buildings.Media.GasesPTDecoupled.SimpleAir.density Buildings.Media.GasesPTDecoupled.SimpleAir.density

return density of ideal gas

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function density "return density of ideal gas"
   extends Modelica.Icons.Function;
   input ThermodynamicState state "Thermodynamic state record";
   output Density d "Density";
algorithm 
   d := dStp*state.p/pStp;
end density;

Buildings.Media.GasesPTDecoupled.SimpleAir.specificInternalEnergy Buildings.Media.GasesPTDecoupled.SimpleAir.specificInternalEnergy

Return specific internal energy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyuSpecific internal energy [J/kg]

Modelica definition

redeclare function specificInternalEnergy 
  "Return specific internal energy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEnergy u "Specific internal energy";
algorithm 
  u := specificEnthalpy(state) - pStp/dStp;
end specificInternalEnergy;

Buildings.Media.GasesPTDecoupled.SimpleAir.specificEntropy Buildings.Media.GasesPTDecoupled.SimpleAir.specificEntropy

Return specific entropy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare replaceable function specificEntropy 
  "Return specific entropy"
   extends Modelica.Icons.Function;
   input ThermodynamicState state "Thermodynamic state record";
   output SpecificEntropy s "Specific entropy";
algorithm 
   s := cp_const*Modelica.Math.log(state.T/T0);// - R_gas*Modelica.Math.log(state.p/reference_p);
end specificEntropy;

Buildings.Media.GasesPTDecoupled.SimpleAir.enthalpyOfCondensingGas Buildings.Media.GasesPTDecoupled.SimpleAir.enthalpyOfCondensingGas

Enthalpy of steam per unit mass of steam

Information

Dummy function that returns 0.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhSteam enthalpy [J/kg]

Modelica definition

replaceable function enthalpyOfCondensingGas 
  "Enthalpy of steam per unit mass of steam"
  extends Modelica.Icons.Function;
  input Temperature T "Temperature";
  output SpecificEnthalpy h "Steam enthalpy";
algorithm 
  h := 0;
end enthalpyOfCondensingGas;

Buildings.Media.GasesPTDecoupled.SimpleAir.saturationPressure Buildings.Media.GasesPTDecoupled.SimpleAir.saturationPressure

Return saturation pressure of condensing fluid

Information

Dummy function that returns 0.

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

Inputs

TypeNameDefaultDescription
TemperatureTsat Saturation temperature [K]

Outputs

TypeNameDescription
AbsolutePressurepsatSaturation pressure [Pa]

Modelica definition

replaceable function saturationPressure 
  "Return saturation pressure of condensing fluid"
  extends Modelica.Icons.Function;
  input Temperature Tsat "Saturation temperature";
  output AbsolutePressure psat "Saturation pressure";
algorithm 
  psat := 0;
end saturationPressure;

Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties.InputAbsolutePressure

Pressure as input signal connector

Modelica definition

connector InputAbsolutePressure = input SI.AbsolutePressure 
  "Pressure as input signal connector";

Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties.InputMassFraction

Mass fraction as input signal connector

Modelica definition

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

Buildings.Media.GasesPTDecoupled.SimpleAir.BaseProperties.InputSpecificEnthalpy

Specific enthalpy as input signal connector

Modelica definition

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

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