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.).
Name | Description |
---|---|
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=101325 | Pressure for which dStp is defined |
dStp=1.2 | Fluid density at pressure pStp |
BaseProperties | Basic medium properties |
setState_dTX | Return thermodynamic state from d, T, and X or Xi |
density | return density of ideal gas |
specificEntropy | Return specific entropy |
Inherited | |
cp_const | Constant specific heat capacity at constant pressure |
cv_const=cp_const - R_gas | Constant specific heat capacity at constant volume |
R_gas | medium specific gas constant |
MM_const | Molar mass |
eta_const | Constant dynamic viscosity |
lambda_const | Constant thermal conductivity |
T_min | Minimum temperature valid for medium model |
T_max | Maximum temperature valid for medium model |
T0=reference_T | Zero enthalpy temperature |
ThermodynamicState | Thermodynamic state of ideal gas |
FluidConstants | fluid constants |
setState_pTX | Return thermodynamic state from p, T, and X or Xi |
setState_phX | Return thermodynamic state from p, h, and X or Xi |
setState_psX | Return thermodynamic state from p, s, and X or Xi |
pressure | Return pressure of ideal gas |
temperature | Return temperature of ideal gas |
specificEnthalpy | Return specific enthalpy |
specificInternalEnergy | Return specific internal energy |
specificGibbsEnergy | Return specific Gibbs energy |
specificHelmholtzEnergy | Return specific Helmholtz energy |
dynamicViscosity | Return dynamic viscosity |
thermalConductivity | Return thermal conductivity |
specificHeatCapacityCp | Return specific heat capacity at constant pressure |
specificHeatCapacityCv | Return specific heat capacity at constant volume |
isentropicExponent | Return isentropic exponent |
velocityOfSound | Return velocity of sound |
specificEnthalpy_pTX | Return specific enthalpy from p, T, and X or Xi |
temperature_phX | Return temperature from p, h, and X or Xi |
density_phX | Return density from p, h, and X or Xi |
setState_pT | Return thermodynamic state from p and T |
setState_ph | Return thermodynamic state from p and h |
setState_ps | Return thermodynamic state from p and s |
setState_dT | Return thermodynamic state from d and T |
density_ph | Return density from p and h |
temperature_ph | Return temperature from p and h |
pressure_dT | Return pressure from d and T |
specificEnthalpy_dT | Return specific enthalpy from d and T |
specificEnthalpy_ps | Return specific enthalpy from p and s |
temperature_ps | Return temperature from p and s |
density_ps | Return density from p and s |
specificEnthalpy_pT | Return specific enthalpy from p and T |
density_pT | Return density from p and T |
ThermoStates | Enumeration 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=101325 | Reference pressure of Medium: default 1 atmosphere |
reference_T=298.15 | Reference temperature of Medium: default 25 deg Celsius |
reference_X=fill(1/nX, nX) | Default mass fractions of medium |
p_default=101325 | Default value for pressure of medium (for initialization) |
T_default=Modelica.SIunits.Conversions.from_degC(20) | Default value for temperature of medium (for initialization) |
h_default=specificEnthalpy_pTX(p_default, T_default, X_default) | Default value for specific enthalpy of medium (for initialization) |
X_default=reference_X | Default value for mass fractions of medium (for initialization) |
nS=size(substanceNames, 1) | Number of substances |
nX=nS | Number of mass fractions |
nXi=if fixedX then 0 else if reducedX then nS - 1 else nS | Number of structurally independent mass fractions (see docu for details) |
nC=size(extraPropertiesNames, 1) | Number of extra (outside of standard mass-balance) transported properties |
setSmoothState | Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b |
prandtlNumber | Return the Prandtl number |
heatCapacity_cp | alias for deprecated name |
heatCapacity_cv | alias for deprecated name |
isentropicEnthalpy | Return isentropic enthalpy |
isobaricExpansionCoefficient | Return overall the isobaric expansion coefficient beta |
beta | alias for isobaricExpansionCoefficient for user convenience |
isothermalCompressibility | Return overall the isothermal compressibility factor |
kappa | alias of isothermalCompressibility for user convenience |
density_derp_h | Return density derivative wrt pressure at const specific enthalpy |
density_derh_p | Return density derivative wrt specific enthalpy at constant pressure |
density_derp_T | Return density derivative wrt pressure at const temperature |
density_derT_p | Return density derivative wrt temperature at constant pressure |
density_derX | Return density derivative wrt mass fraction |
molarMass | Return the molar mass of the medium |
density_pTX | Return density from p, T, and X or Xi |
temperature_psX | Return temperature from p,s, and X or Xi |
density_psX | Return density from p, s, and X or Xi |
specificEnthalpy_psX | Return specific enthalpy from p, s, and X or Xi |
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 |
Choices | Types, constants to define menu choices |
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";
Type | Name | Default | Description |
---|---|---|---|
Boolean | standardOrderComponents | true | if true, and reducedX = true, the last element of X will be computed from the other ones |
Advanced | |||
Boolean | preferredMediumStates | false | = true if StateSelect.prefer shall be used for the independent property variables of the medium |
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"; 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); u = h-R*T; R = R_gas; // d = p/(R*T); d/dStp = p/pStp; MM = MM_const; state.T = T; state.p = p; end BaseProperties;
Type | Name | Default | Description |
---|---|---|---|
Density | d | density [kg/m3] | |
Temperature | T | Temperature [K] | |
MassFraction | X[:] | fill(0, 0) | Mass fractions [kg/kg] |
Type | Name | Description |
---|---|---|
ThermodynamicState | state |
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;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | thermodynamic state record |
Type | Name | Description |
---|---|---|
Density | d | Density [kg/m3] |
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;
Type | Name | Default | Description |
---|---|---|---|
ThermodynamicState | state | thermodynamic state record |
Type | Name | Description |
---|---|---|
SpecificEntropy | s | Specific entropy [J/(kg.K)] |
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;
connector InputAbsolutePressure = input SI.AbsolutePressure "Pressure as input signal connector";
connector InputMassFraction = input SI.MassFraction "Mass fraction as input signal connector";
connector InputSpecificEnthalpy = input SI.SpecificEnthalpy "Specific enthalpy as input signal connector";