LBL logo

Buildings.Media.Interfaces.PartialSimpleMedium

Medium model with linear dependency of u, h from temperature. Most other quantities are constant.

Information

This medium model is identical to Modelica.Media.Interfaces.PartialSimpleMedium, but it allows to define a compressibility of the medium. This helps breaking algebraic loops, but the system gets stiff. The compressibility is defined by the constant kappa_const. If kappa_const=0, then the density is constant. Otherwise, the density is
  rho(p) = rho(p0) * ( 1 + kappa_const * (p-p0))

Extends from Modelica.Media.Interfaces.PartialPureSubstance (base class for pure substances of one chemical substance).

Package Content

NameDescription
cp_constConstant specific heat capacity at constant pressure
cv_constConstant specific heat capacity at constant volume
d_constConstant density
eta_constConstant dynamic viscosity
lambda_constConstant thermal conductivity
a_constConstant velocity of sound
T_minMinimum temperature valid for medium model
T_maxMaximum temperature valid for medium model
T0=reference_TZero enthalpy temperature
MM_constMolar mass
fluidConstantsfluid constants
Buildings.Media.Interfaces.PartialSimpleMedium.ThermodynamicState ThermodynamicState Thermodynamic state
kappa_const=0Compressibility factor at constant temperature
p0=3E5Reference pressure for compressibility and default medium pressure
constantDensity=(kappa_const <= 1E-20)Flag, true if density is modeled as a constant
Buildings.Media.Interfaces.PartialSimpleMedium.BaseProperties BaseProperties  
Buildings.Media.Interfaces.PartialSimpleMedium.setState_pTX setState_pTX Return thermodynamic state from p, T, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.setState_phX setState_phX Return thermodynamic state from p, h, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.setState_psX setState_psX Return thermodynamic state from p, s, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.setState_dTX setState_dTX Return thermodynamic state from d, T, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.setSmoothState setSmoothState Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b
Buildings.Media.Interfaces.PartialSimpleMedium.dynamicViscosity dynamicViscosity Return dynamic viscosity
Buildings.Media.Interfaces.PartialSimpleMedium.thermalConductivity thermalConductivity Return thermal conductivity
Buildings.Media.Interfaces.PartialSimpleMedium.pressure pressure Return pressure
Buildings.Media.Interfaces.PartialSimpleMedium.temperature temperature Return temperature
Buildings.Media.Interfaces.PartialSimpleMedium.density density Return density
Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy specificEnthalpy Return specific enthalpy
Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCp specificHeatCapacityCp Return specific heat capacity at constant pressure
Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCv specificHeatCapacityCv Return specific heat capacity at constant volume
Buildings.Media.Interfaces.PartialSimpleMedium.isentropicExponent isentropicExponent Return isentropic exponent
Buildings.Media.Interfaces.PartialSimpleMedium.velocityOfSound velocityOfSound Return velocity of sound
Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy_pTX specificEnthalpy_pTX Return specific enthalpy from p, T, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.temperature_phX temperature_phX Return temperature from p, h, and X or Xi
Buildings.Media.Interfaces.PartialSimpleMedium.density_phX density_phX Return density from p, h, and X or Xi
Inherited
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.FluidConstants FluidConstants critical, triple, molecular and other standard data of fluid
Modelica.Media.Interfaces.PartialMedium.prandtlNumber prandtlNumber Return the Prandtl number
Modelica.Media.Interfaces.PartialMedium.specificInternalEnergy specificInternalEnergy Return specific internal energy
Modelica.Media.Interfaces.PartialMedium.specificEntropy specificEntropy Return specific entropy
Modelica.Media.Interfaces.PartialMedium.specificGibbsEnergy specificGibbsEnergy Return specific Gibbs energy
Modelica.Media.Interfaces.PartialMedium.specificHelmholtzEnergy specificHelmholtzEnergy Return specific Helmholtz energy
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.isentropicEnthalpy isentropicEnthalpy Return isentropic enthalpy
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.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 SpecificHeatCapacity cp_const 
  "Constant specific heat capacity at constant pressure";

  constant SpecificHeatCapacity cv_const 
  "Constant specific heat capacity at constant volume";

  constant Density d_const "Constant density";

  constant DynamicViscosity eta_const "Constant dynamic viscosity";

  constant ThermalConductivity lambda_const "Constant thermal conductivity";

  constant VelocityOfSound a_const "Constant velocity of sound";

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

  constant MolarMass MM_const "Molar mass";

  constant FluidConstants[nS] fluidConstants "fluid constants";

  constant Real kappa_const(unit="1/Pa") = 0 
  "Compressibility factor at constant temperature";

  constant Modelica.SIunits.AbsolutePressure p0 = 3E5 
  "Reference pressure for compressibility and default medium pressure";

  constant Boolean constantDensity = (kappa_const <= 1E-20) 
  "Flag, true if density is modeled as a constant";


Buildings.Media.Interfaces.PartialSimpleMedium.ThermodynamicState Buildings.Media.Interfaces.PartialSimpleMedium.ThermodynamicState

Thermodynamic state

Information

Extends from (Minimal variable set that is available as input argument to every medium function).

Modelica definition

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

Buildings.Media.Interfaces.PartialSimpleMedium.BaseProperties Buildings.Media.Interfaces.PartialSimpleMedium.BaseProperties

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.

Extends from .

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
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 = cp_const*(T-T0);
  h = specificEnthalpy_pTX(p,T,X);
  u = cv_const*(T-T0);
  // original equation d = d_const;
  d = if constantDensity then d_const else d_const * (1+kappa_const*(p-p0));
 // d = d_const * (1+kT*(T-T0)/T0); "this gives large coupled equations"
  R = 0;
  MM = MM_const;
  state.T = T;
  state.p = p;
end BaseProperties;

Buildings.Media.Interfaces.PartialSimpleMedium.setState_pTX Buildings.Media.Interfaces.PartialSimpleMedium.setState_pTX

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

Information

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

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestatethermodynamic state record

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[:]=reference_X "Mass fractions";
  output ThermodynamicState state "thermodynamic state record";
algorithm 
  state := ThermodynamicState(p=p,T=T);
end setState_pTX;

Buildings.Media.Interfaces.PartialSimpleMedium.setState_phX Buildings.Media.Interfaces.PartialSimpleMedium.setState_phX

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

Information

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

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestatethermodynamic state record

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[:]=reference_X "Mass fractions";
  output ThermodynamicState state "thermodynamic state record";
algorithm 
  state := ThermodynamicState(p=p,T=temperature_phX(p, h, X));
end setState_phX;

Buildings.Media.Interfaces.PartialSimpleMedium.setState_psX Buildings.Media.Interfaces.PartialSimpleMedium.setState_psX

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

Information

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

Inputs

TypeNameDefaultDescription
AbsolutePressurep Pressure [Pa]
SpecificEntropys Specific entropy [J/(kg.K)]
MassFractionX[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestatethermodynamic state record

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[:]=reference_X "Mass fractions";
  output ThermodynamicState state "thermodynamic state record";
algorithm 
  state := ThermodynamicState(p=p,T=Modelica.Math.exp(s/cp_const + Modelica.Math.log(T0))) 
    "here the incompressible limit is used, with cp as heat capacity";
end setState_psX;

Buildings.Media.Interfaces.PartialSimpleMedium.setState_dTX Buildings.Media.Interfaces.PartialSimpleMedium.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[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestatethermodynamic state record

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[:]=reference_X "Mass fractions";
  output ThermodynamicState state "thermodynamic state record";
algorithm 
  assert(false,"pressure can not be computed from temperature and density for an incompressible fluid!");
end setState_dTX;

Buildings.Media.Interfaces.PartialSimpleMedium.setSmoothState Buildings.Media.Interfaces.PartialSimpleMedium.setSmoothState

Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b

Information

Extends from (Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b).

Inputs

TypeNameDefaultDescription
Realx m_flow or dp
ThermodynamicStatestate_a Thermodynamic state if x > 0
ThermodynamicStatestate_b Thermodynamic state if x < 0
Realx_small Smooth transition in the region -x_small < x < x_small

Outputs

TypeNameDescription
ThermodynamicStatestateSmooth thermodynamic state for all x (continuous and differentiable)

Modelica definition

redeclare function extends setSmoothState 
  "Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b"
algorithm 
  state := ThermodynamicState(p=Modelica.Media.Common.smoothStep(
      x,
      state_a.p,
      state_b.p,
      x_small), T=Modelica.Media.Common.smoothStep(
      x,
      state_a.T,
      state_b.T,
      x_small));
end setSmoothState;

Buildings.Media.Interfaces.PartialSimpleMedium.dynamicViscosity Buildings.Media.Interfaces.PartialSimpleMedium.dynamicViscosity

Return dynamic viscosity

Information

Extends from (Return dynamic viscosity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

redeclare function extends dynamicViscosity 
  "Return dynamic viscosity"

algorithm 
  eta := eta_const;
end dynamicViscosity;

Buildings.Media.Interfaces.PartialSimpleMedium.thermalConductivity Buildings.Media.Interfaces.PartialSimpleMedium.thermalConductivity

Return thermal conductivity

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 
  "Return thermal conductivity"

algorithm 
  lambda := lambda_const;
end thermalConductivity;

Buildings.Media.Interfaces.PartialSimpleMedium.pressure Buildings.Media.Interfaces.PartialSimpleMedium.pressure

Return pressure

Information

Extends from (Return pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function extends pressure "Return pressure"

algorithm 
  p := state.p;
end pressure;

Buildings.Media.Interfaces.PartialSimpleMedium.temperature Buildings.Media.Interfaces.PartialSimpleMedium.temperature

Return temperature

Information

Extends from (Return temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function extends temperature "Return temperature"

algorithm 
  T := state.T;
end temperature;

Buildings.Media.Interfaces.PartialSimpleMedium.density Buildings.Media.Interfaces.PartialSimpleMedium.density

Return density

Information

Extends from (Return density).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function extends density "Return density"

algorithm 
  d := if constantDensity then d_const else d_const * (1+kappa_const*(state.p-p0));
end density;

Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy

Return 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 
  "Return specific enthalpy"

algorithm 
  h := cp_const*(state.T-T0);
end specificEnthalpy;

Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCp Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCp

Return specific heat capacity 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 function extends specificHeatCapacityCp 
  "Return specific heat capacity at constant pressure"

algorithm 
  cp := cp_const;
end specificHeatCapacityCp;

Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCv Buildings.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCv

Return specific heat capacity 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 function extends specificHeatCapacityCv 
  "Return specific heat capacity at constant volume"

algorithm 
  cv := cv_const;
end specificHeatCapacityCv;

Buildings.Media.Interfaces.PartialSimpleMedium.isentropicExponent Buildings.Media.Interfaces.PartialSimpleMedium.isentropicExponent

Return isentropic exponent

Information

Extends from (Return isentropic exponent).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

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.PartialSimpleMedium.velocityOfSound Buildings.Media.Interfaces.PartialSimpleMedium.velocityOfSound

Return velocity of sound

Information

Extends from (Return velocity of sound).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
VelocityOfSoundaVelocity of sound [m/s]

Modelica definition

redeclare function extends velocityOfSound 
  "Return velocity of sound "

algorithm 
  a := a_const;
end velocityOfSound;

Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy_pTX Buildings.Media.Interfaces.PartialSimpleMedium.specificEnthalpy_pTX

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

Information

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

Inputs

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

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [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";
algorithm 
  h := cp_const*(T-T0);
end specificEnthalpy_pTX;

Buildings.Media.Interfaces.PartialSimpleMedium.temperature_phX Buildings.Media.Interfaces.PartialSimpleMedium.temperature_phX

Return temperature from p, h, and X or Xi

Information

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

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 := T0 + h/cp_const;
end temperature_phX;

Buildings.Media.Interfaces.PartialSimpleMedium.density_phX Buildings.Media.Interfaces.PartialSimpleMedium.density_phX

Return density from p, h, and X or Xi

Information

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

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;

Automatically generated Thu Jul 26 10:22:33 2012.