Buildings.Media.Interfaces.PartialSimpleIdealGasMedium

Medium model of Ideal gas with constant cp and cv. All other quantities, e.g. transport properties, are constant.

Information


This is the most simple incompressible medium model, where specific enthalpy h and specific internal energy u are only a function of temperature T and all other provided medium quantities are assumed to be constant.

This partial package is identical to Modelica.Media.Interfaces.PartialSimpleIdealGasMedium but it fixes two bugs in the original implementation. Namely, in the functions Modelica.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy and Modelica.Media.Interfaces.PartialSimpleIdealGasMedium.specificInternalEnergy the computation of the enthalpy or internal energy does not use the reference temperature T0. With this implementation, initial states of fluid volumes can be far away from the steady-state value because of this inconsistent implementation. Because these two functions are not replaceable, we needed to make a copy of the whole package, and reimplement these functions. When the Buildings library is upgraded to to Modelica 3.0.0, it should be safe to remove this package as the bug is fixed since Modelica 3.0.0 (but not in Modelica 2.2.2).

Package Content

NameDescription
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
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.BaseProperties BaseProperties Base properties
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.setState_dTX setState_dTX Return thermodynamic state from d, T, and X or Xi
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.density density return density of ideal gas
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy specificEnthalpy Return specific enthalpy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificInternalEnergy specificInternalEnergy Return specific internal energy
Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEntropy specificEntropy Return specific entropy
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 dya_constnamic 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 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

Types and constants

  constant SpecificHeatCapacity cp_const 
  "Constant specific heat capacity at constant pressure";

  constant SpecificHeatCapacity cv_const= cp_const/R_gas 
  "Constant specific heat capacity at constant volume";

  constant SpecificHeatCapacity R_gas "medium specific gas constant";

  constant MolarMass MM_const "Molar mass";

  constant DynamicViscosity eta_const "Constant dynamic viscosity";

  constant ThermalConductivity lambda_const "Constant thermal conductivity";

  constant Temperature T_min "Minimum temperature valid for medium model";

  constant Temperature T_max "Maximum temperature valid for medium model";

  constant Temperature T0= reference_T "Zero enthalpy temperature";


Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.ThermodynamicState Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.ThermodynamicState

thermodynamic state

Modelica definition

redeclare replaceable record extends ThermodynamicState 
  "thermodynamic state" 
  AbsolutePressure p "Absolute pressure of medium";
  Temperature T "Temperature of medium";
end ThermodynamicState;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.BaseProperties Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.BaseProperties

Base 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 extends BaseProperties(
        T(stateSelect=StateSelect.prefer)) "Base properties" 
equation 
      assert(T >= T_min and T <= T_max, "
Temperature T (= " + String(T) + " K) is not
in the allowed range (" + String(T_min) + " K <= T <= " + String(T_max)
         + " K)
required from medium model \"" + mediumName + "\".
");
  h = specificEnthalpy_pTX(p,T,X);
  u = h-R*T;
  R = R_gas;
  d = p/(R*T);
  MM = MM_const;
  state.T = T;
  state.p = p;
  
end BaseProperties;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_pTX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_pTX

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

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
TemperatureT Temperature [K]
MassFractionX[:]fill(0, 0)Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_pTX 
  "Return thermodynamic state from p, T, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input Temperature T "Temperature";
  input MassFraction X[:] = fill(0,0) "Mass fractions";
  output ThermodynamicState state;
algorithm 
  state := ThermodynamicState(p=p,T=T);
end setState_pTX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_phX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_phX

Return thermodynamic state from p, h, and X or Xi

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEnthalpyh Specific enthalpy [J/kg]
MassFractionX[:]fill(0, 0)Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare function setState_phX 
  "Return thermodynamic state from p, h, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "Specific enthalpy";
  input MassFraction X[:] = fill(0,0) "Mass fractions";
  output ThermodynamicState state;
algorithm 
  state := ThermodynamicState(p=p,T=T0+h/cp_const);
end setState_phX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_psX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_psX

Return thermodynamic state from p, s, and X or Xi

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEntropys Specific entropy [J/(kg.K)]
MassFractionX[:]fill(0, 0)Mass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestate 

Modelica definition

redeclare replaceable function setState_psX 
  "Return thermodynamic state from p, s, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEntropy s "Specific entropy";
  input MassFraction X[:] = fill(0,0) "Mass fractions";
  output ThermodynamicState state;
algorithm 
  state := ThermodynamicState(p=p,T=Modelica.Math.exp(s/cp_const + Modelica.Math.log(reference_T))
                              + R_gas*Modelica.Math.log(p/reference_p));
end setState_psX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_dTX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.setState_dTX

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

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*R_gas*T,T=T);
end setState_dTX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.pressure Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.pressure

return pressure of ideal gas

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function extends pressure "return pressure of ideal gas" 
algorithm 
  p := state.p;
end pressure;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature

return temperature of ideal gas

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function extends temperature 
  "return temperature of ideal gas" 
algorithm 
  T := state.T;
end temperature;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density

return density of ideal gas

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function extends density "return density of ideal gas" 
algorithm 
  d := state.p/(R_gas*state.T);
end density;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy

Return specific enthalpy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

redeclare function extends specificEnthalpy 
  "Return specific enthalpy" 
    extends Modelica.Icons.Function;
algorithm 
  h := cp_const*(state.T-T0); // this is one of the bug fixes
end specificEnthalpy;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificInternalEnergy Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificInternalEnergy

Return specific internal energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergyuSpecific internal energy [J/kg]

Modelica definition

redeclare function extends specificInternalEnergy 
  "Return specific internal energy" 
  extends Modelica.Icons.Function;
algorithm 
  u := cp_const*(state.T-T0) - R_gas*state.T; // this is one of the bug fixes
end specificInternalEnergy;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEntropy Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEntropy

Return specific entropy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

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

Modelica definition

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

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificGibbsEnergy Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificGibbsEnergy

Return specific Gibbs energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

redeclare function extends specificGibbsEnergy 
  "Return specific Gibbs energy" 
  extends Modelica.Icons.Function;
algorithm 
  g := cp_const*state.T - state.T*specificEntropy(state);
end specificGibbsEnergy;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHelmholtzEnergy Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHelmholtzEnergy

Return specific Helmholtz energy

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

redeclare function extends specificHelmholtzEnergy 
  "Return specific Helmholtz energy" 
  extends Modelica.Icons.Function;
algorithm 
  f := cp_const*state.T - R_gas*state.T - state.T*specificEntropy(state);
end specificHelmholtzEnergy;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.dynamicViscosity Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.dynamicViscosity

Return dya_constnamic viscosity

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

redeclare function extends dynamicViscosity 
  "Return dya_constnamic viscosity" 
algorithm 
  eta := eta_const;
end dynamicViscosity;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.thermalConductivity Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.thermalConductivity

Return thermal conductivity

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
ThermalConductivitylambdaThermal conductivity [W/(m.K)]

Modelica definition

redeclare function extends thermalConductivity 
  "Return thermal conductivity" 
algorithm 
  lambda := lambda_const;
end thermalConductivity;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCp Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCp

Return specific heat capacity at constant pressure

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificHeatCapacitycpSpecific heat capacity at constant pressure [J/(kg.K)]

Modelica definition

redeclare function extends specificHeatCapacityCp 
  "Return specific heat capacity at constant pressure" 
algorithm 
  cp := cp_const;
end specificHeatCapacityCp;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCv Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificHeatCapacityCv

Return specific heat capacity at constant volume

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
SpecificHeatCapacitycvSpecific heat capacity at constant volume [J/(kg.K)]

Modelica definition

redeclare function extends specificHeatCapacityCv 
  "Return specific heat capacity at constant volume" 
algorithm 
  cv := cv_const;
end specificHeatCapacityCv;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isentropicExponent Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.isentropicExponent

Return isentropic exponent

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
IsentropicExponentgammaIsentropic exponent [1]

Modelica definition

redeclare function extends isentropicExponent 
  "Return isentropic exponent" 
algorithm 
  gamma := cp_const/cv_const;
end isentropicExponent;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.velocityOfSound Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.velocityOfSound

Velocity of sound

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate  

Outputs

TypeNameDescription
VelocityOfSoundaVelocity of sound [m/s]

Modelica definition

redeclare function extends velocityOfSound "Velocity of sound " 
algorithm 
  a := sqrt(cp_const/cv_const*R_gas*state.T);
end velocityOfSound;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy_pTX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.specificEnthalpy_pTX

Return specific enthalpy from p, T, and X or Xi

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
TemperatureT Temperature [K]
MassFractionX[nX] Mass fractions [kg/kg]

Outputs

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

Modelica definition

redeclare function specificEnthalpy_pTX 
  "Return specific enthalpy from p, T, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input Temperature T "Temperature";
  input MassFraction X[nX] "Mass fractions";
  output SpecificEnthalpy h "Specific enthalpy at p, T, X";
algorithm 
  h := cp_const*(T-T0);
end specificEnthalpy_pTX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature_phX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.temperature_phX

Return temperature from p, h, and X or Xi

Inputs

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

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function temperature_phX 
  "Return temperature from p, h, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "Specific enthalpy";
  input MassFraction X[nX] "Mass fractions";
  output Temperature T "Temperature";
algorithm 
  T := h/cp_const + T0;
end temperature_phX;

Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_phX Buildings.Media.Interfaces.PartialSimpleIdealGasMedium.density_phX

Return density from p, h, and X or Xi

Inputs

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

Outputs

TypeNameDescription
Densityddensity [kg/m3]

Modelica definition

redeclare function density_phX 
  "Return density from p, h, and X or Xi" 
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEnthalpy h "Specific enthalpy";
  input MassFraction X[nX] "Mass fractions";
  output Density d "density";
algorithm 
  d := density(setState_phX(p,h,X));
end density_phX;

HTML-documentation generated by Dymola Tue Sep 30 14:24:46 2008.