Modelica.Media.Interfaces.PartialMedium

Partial medium properties (base package of all media packages)

Information


PartialMedium is a package and contains all declarations for a medium. This means that constants, models, and functions are defined that every medium is supposed to support (some of them are optional). A medium package inherits from PartialMedium and provides the equations for the medium. The details of this package are described in Modelica.Media.UsersGuide.

Extends from Modelica.Icons.MaterialPropertiesPackage (Icon for package containing property classes).

Package Content

NameDescription
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.ThermodynamicState ThermodynamicState Minimal variable set that is available as input argument to every medium function
Modelica.Media.Interfaces.PartialMedium.BaseProperties BaseProperties Base properties (p, d, T, h, u, R, MM and, if applicable, X and Xi) of a medium
Modelica.Media.Interfaces.PartialMedium.setState_pTX setState_pTX Return thermodynamic state as function of p, T and composition X or Xi
Modelica.Media.Interfaces.PartialMedium.setState_phX setState_phX Return thermodynamic state as function of p, h and composition X or Xi
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.setState_dTX setState_dTX Return thermodynamic state as function of d, T 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.dynamicViscosity dynamicViscosity Return dynamic viscosity
Modelica.Media.Interfaces.PartialMedium.thermalConductivity thermalConductivity Return thermal conductivity
Modelica.Media.Interfaces.PartialMedium.prandtlNumber prandtlNumber Return the Prandtl number
Modelica.Media.Interfaces.PartialMedium.pressure pressure Return pressure
Modelica.Media.Interfaces.PartialMedium.temperature temperature Return temperature
Modelica.Media.Interfaces.PartialMedium.density density Return density
Modelica.Media.Interfaces.PartialMedium.specificEnthalpy specificEnthalpy Return specific enthalpy
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.specificHeatCapacityCp specificHeatCapacityCp Return specific heat capacity at constant pressure
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp heatCapacity_cp alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.specificHeatCapacityCv specificHeatCapacityCv Return specific heat capacity at constant volume
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 Modelica.Media.Interfaces.PartialMedium.Choices.IndependentVariables
    ThermoStates "Enumeration type for independent variables";

  constant String mediumName = "unusablePartialMedium" "Name of the medium";

  constant String substanceNames[:]={mediumName} 
  "Names of the mixture substances. Set substanceNames={mediumName} if only one substance.";

  constant String extraPropertiesNames[:]=fill("", 0) 
  "Names of the additional (extra) transported properties. Set extraPropertiesNames=fill(\"\",0) if unused";

  constant Boolean singleState 
  "= true, if u and d are not a function of pressure";

  constant Boolean reducedX=true 
  "= true if medium contains the equation sum(X) = 1.0; set reducedX=true if only one substance (see docu for details)";

  constant Boolean fixedX=false 
  "= true if medium contains the equation X = reference_X";

  constant AbsolutePressure reference_p=101325 
  "Reference pressure of Medium: default 1 atmosphere";

  constant Temperature reference_T=298.15 
  "Reference temperature of Medium: default 25 deg Celsius";

  constant MassFraction reference_X[nX]= fill(1/nX, nX) 
  "Default mass fractions of medium";

  constant AbsolutePressure p_default=101325 
  "Default value for pressure of medium (for initialization)";

  constant Temperature T_default = Modelica.SIunits.Conversions.from_degC(20) 
  "Default value for temperature of medium (for initialization)";

  constant SpecificEnthalpy h_default = specificEnthalpy_pTX(p_default, T_default, X_default) 
  "Default value for specific enthalpy of medium (for initialization)";

  constant MassFraction X_default[nX]=reference_X 
  "Default value for mass fractions of medium (for initialization)";

  final constant Integer nS=size(substanceNames, 1) "Number of substances";

  constant Integer nX = nS "Number of mass fractions";

  constant Integer nXi=if fixedX then 0 else if reducedX then nS - 1 else nS 
  "Number of structurally independent mass fractions (see docu for details)";

  final constant Integer nC=size(extraPropertiesNames, 1) 
  "Number of extra (outside of standard mass-balance) transported properties";

  constant Real C_nominal[nC](min=fill(Modelica.Constants.eps, nC)) = 1.0e-6*ones(nC) 
  "Default for the nominal values for the extra properties";

  type AbsolutePressure = SI.AbsolutePressure (
      min=0,
      max=1.e8,
      nominal=1.e5,
      start=1.e5) "Type for absolute pressure with medium specific attributes";

  type Density = SI.Density (
      min=0,
      max=1.e5,
      nominal=1,
      start=1) "Type for density with medium specific attributes";

  type DynamicViscosity = SI.DynamicViscosity (
      min=0,
      max=1.e8,
      nominal=1.e-3,
      start=1.e-3) "Type for dynamic viscosity with medium specific attributes";

  type EnthalpyFlowRate = SI.EnthalpyFlowRate (
      nominal=1000.0,
      min=-1.0e8,
      max=1.e8) "Type for enthalpy flow rate with medium specific attributes";

  type MassFlowRate = SI.MassFlowRate (
      quantity="MassFlowRate." + mediumName,
      min=-1.0e5,
      max=1.e5) "Type for mass flow rate with medium specific attributes";

  type MassFraction = Real (
      quantity="MassFraction",
      final unit="kg/kg",
      min=0,
      max=1,
      nominal=0.1) "Type for mass fraction with medium specific attributes";

  type MoleFraction = Real (
      quantity="MoleFraction",
      final unit="mol/mol",
      min=0,
      max=1,
      nominal=0.1) "Type for mole fraction with medium specific attributes";

  type MolarMass = SI.MolarMass (
      min=0.001,
      max=0.25,
      nominal=0.032) "Type for molar mass with medium specific attributes";

  type MolarVolume = SI.MolarVolume (
      min=1e-6,
      max=1.0e6,
      nominal=1.0) "Type for molar volume with medium specific attributes";

  type IsentropicExponent = SI.RatioOfSpecificHeatCapacities (
      min=1,
      max=500000,
      nominal=1.2,
      start=1.2) "Type for isentropic exponent with medium specific attributes";

  type SpecificEnergy = SI.SpecificEnergy (
      min=-1.0e8,
      max=1.e8,
      nominal=1.e6) "Type for specific energy with medium specific attributes";

  type SpecificInternalEnergy = SpecificEnergy 
  "Type for specific internal energy with medium specific attributes";

  type SpecificEnthalpy = SI.SpecificEnthalpy (
      min=-1.0e10,
      max=1.e10,
      nominal=1.e6) 
  "Type for specific enthalpy with medium specific attributes";

  type SpecificEntropy = SI.SpecificEntropy (
      min=-1.e7,
      max=1.e7,
      nominal=1.e3) "Type for specific entropy with medium specific attributes";

  type SpecificHeatCapacity = SI.SpecificHeatCapacity (
      min=0,
      max=1.e7,
      nominal=1.e3,
      start=1.e3) 
  "Type for specific heat capacity with medium specific attributes";

  type SurfaceTension = SI.SurfaceTension 
  "Type for surface tension with medium specific attributes";

  type Temperature = SI.Temperature (
      min=1,
      max=1.e4,
      nominal=300,
      start=300) "Type for temperature with medium specific attributes";

  type ThermalConductivity = SI.ThermalConductivity (
      min=0,
      max=500,
      nominal=1,
      start=1) "Type for thermal conductivity with medium specific attributes";

  type PrandtlNumber = SI.PrandtlNumber (
      min=1e-3,
      max=1e5,
      nominal=1.0) "Type for Prandtl number with medium specific attributes";

  type VelocityOfSound = SI.Velocity (
      min=0,
      max=1.e5,
      nominal=1000,
      start=1000) "Type for velocity of sound with medium specific attributes";

  type ExtraProperty = Real (min=0.0, start=1.0) 
  "Type for unspecified, mass-specific property transported by flow";

  type CumulativeExtraProperty = Real (min=0.0, start=1.0) 
  "Type for conserved integral of unspecified, mass specific property";

  type ExtraPropertyFlowRate = Real(unit="kg/s") 
  "Type for flow rate of unspecified, mass-specific property";

  type IsobaricExpansionCoefficient = Real (
      min=0,
      max=1.0e8,
      unit="1/K") 
  "Type for isobaric expansion coefficient with medium specific attributes";

  type DipoleMoment = Real (
      min=0.0,
      max=2.0,
      unit="debye",
      quantity="ElectricDipoleMoment") 
  "Type for dipole moment with medium specific attributes";

  type DerDensityByPressure = SI.DerDensityByPressure 
  "Type for partial derivative of density with resect to pressure with medium specific attributes";

  type DerDensityByEnthalpy = SI.DerDensityByEnthalpy 
  "Type for partial derivative of density with resect to enthalpy with medium specific attributes";

  type DerEnthalpyByPressure = SI.DerEnthalpyByPressure 
  "Type for partial derivative of enthalpy with resect to pressure with medium specific attributes";

  type DerDensityByTemperature = SI.DerDensityByTemperature 
  "Type for partial derivative of density with resect to temperature with medium specific attributes";


Modelica.Media.Interfaces.PartialMedium.FluidConstants Modelica.Media.Interfaces.PartialMedium.FluidConstants

critical, triple, molecular and other standard data of fluid

Information

Extends from Modelica.Icons.Record (Icon for records).

Modelica definition

replaceable record FluidConstants 
  "critical, triple, molecular and other standard data of fluid"
  extends Modelica.Icons.Record;
  String iupacName "complete IUPAC name (or common name, if non-existent)";
  String casRegistryNumber 
    "chemical abstracts sequencing number (if it exists)";
  String chemicalFormula 
    "Chemical formula, (brutto, nomenclature according to Hill";
  String structureFormula "Chemical structure formula";
  MolarMass molarMass "molar mass";
end FluidConstants;

Modelica.Media.Interfaces.PartialMedium.ThermodynamicState Modelica.Media.Interfaces.PartialMedium.ThermodynamicState

Minimal variable set that is available as input argument to every medium function

Information

Extends from Modelica.Icons.Record (Icon for records).

Modelica definition

replaceable record ThermodynamicState 
  "Minimal variable set that is available as input argument to every medium function"
  extends Modelica.Icons.Record;
end ThermodynamicState;

Modelica.Media.Interfaces.PartialMedium.BaseProperties Modelica.Media.Interfaces.PartialMedium.BaseProperties

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

Information


Model BaseProperties is a model within package PartialMedium and contains the declarations of the minimum number of variables that every medium model is supposed to support. A specific medium inherits from model BaseProperties and provides the equations for the basic properties.

The BaseProperties model contains the following 7+nXi variables (nXi is the number of independent mass fractions defined in package PartialMedium):

Variable Unit Description
T K temperature
p Pa absolute pressure
d kg/m3 density
h J/kg specific enthalpy
u J/kg specific internal energy
Xi[nXi] kg/kg independent mass fractions m_i/m
R J/kg.K gas constant
M kg/mol molar mass

In order to implement an actual medium model, one can extend from this base model and add 5 equations that provide relations among these variables. Equations will also have to be added in order to set all the variables within the ThermodynamicState record state.

If standardOrderComponents=true, the full composition vector X[nX] is determined by the equations contained in this base class, depending on the independent mass fraction vector Xi[nXi].

Additional 2 + nXi equations will have to be provided when using the BaseProperties model, in order to fully specify the thermodynamic conditions. The input connector qualifier applied to p, h, and nXi indirectly declares the number of missing equations, permitting advanced equation balance checking by Modelica tools. Please note that this doesn't mean that the additional equations should be connection equations, nor that exactly those variables should be supplied, in order to complete the model. For further information, see the Modelica.Media User's guide, and Section 4.7 (Balanced Models) of the Modelica 3.0 specification.

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

replaceable partial model BaseProperties 
  "Base properties (p, d, T, h, u, R, MM and, if applicable, X and Xi) of a medium"
  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";

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)");
end BaseProperties;

Modelica.Media.Interfaces.PartialMedium.setState_pTX Modelica.Media.Interfaces.PartialMedium.setState_pTX

Return thermodynamic state as function of p, T and composition 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

replaceable partial function setState_pTX 
  "Return thermodynamic state as function of p, T and composition 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";
end setState_pTX;

Modelica.Media.Interfaces.PartialMedium.setState_phX Modelica.Media.Interfaces.PartialMedium.setState_phX

Return thermodynamic state as function of p, h and composition 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

replaceable partial function setState_phX 
  "Return thermodynamic state as function of p, h and composition 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";
end setState_phX;

Modelica.Media.Interfaces.PartialMedium.setState_psX Modelica.Media.Interfaces.PartialMedium.setState_psX

Return thermodynamic state as function of p, s and composition 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

replaceable partial function setState_psX 
  "Return thermodynamic state as function of p, s and composition 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";
end setState_psX;

Modelica.Media.Interfaces.PartialMedium.setState_dTX Modelica.Media.Interfaces.PartialMedium.setState_dTX

Return thermodynamic state as function of d, T and composition 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

replaceable partial function setState_dTX 
  "Return thermodynamic state as function of d, T and composition 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";
end setState_dTX;

Modelica.Media.Interfaces.PartialMedium.setSmoothState Modelica.Media.Interfaces.PartialMedium.setSmoothState

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

Information


This function is used to approximate the equation

    state = if x > 0 then state_a else state_b;

by a smooth characteristic, so that the expression is continuous and differentiable:

   state := smooth(1, if x >  x_small then state_a else
                      if x < -x_small then state_b else f(state_a, state_b));

This is performed by applying function Media.Common.smoothStep(..) on every element of the thermodynamic state record.

If mass fractions X[:] are approximated with this function then this can be performed for all nX mass fractions, instead of applying it for nX-1 mass fractions and computing the last one by the mass fraction constraint sum(X)=1. The reason is that the approximating function has the property that sum(state.X) = 1, provided sum(state_a.X) = sum(state_b.X) = 1. This can be shown by evaluating the approximating function in the abs(x) < x_small region (otherwise state.X is either state_a.X or state_b.X):

    X[1]  = smoothStep(x, X_a[1] , X_b[1] , x_small);
    X[2]  = smoothStep(x, X_a[2] , X_b[2] , x_small);
       ...
    X[nX] = smoothStep(x, X_a[nX], X_b[nX], x_small);

or

    X[1]  = c*(X_a[1]  - X_b[1])  + (X_a[1]  + X_b[1])/2
    X[2]  = c*(X_a[2]  - X_b[2])  + (X_a[2]  + X_b[2])/2;
       ...
    X[nX] = c*(X_a[nX] - X_b[nX]) + (X_a[nX] + X_b[nX])/2;
    c     = (x/x_small)*((x/x_small)^2 - 3)/4

Summing all mass fractions together results in

    sum(X) = c*(sum(X_a) - sum(X_b)) + (sum(X_a) + sum(X_b))/2
           = c*(1 - 1) + (1 + 1)/2
           = 1

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

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

replaceable partial function setSmoothState 
  "Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b"
  extends Modelica.Icons.Function;
  input Real x "m_flow or dp";
  input ThermodynamicState state_a "Thermodynamic state if x > 0";
  input ThermodynamicState state_b "Thermodynamic state if x < 0";
  input Real x_small(min=0) 
    "Smooth transition in the region -x_small < x < x_small";
  output ThermodynamicState state 
    "Smooth thermodynamic state for all x (continuous and differentiable)";
end setSmoothState;

Modelica.Media.Interfaces.PartialMedium.dynamicViscosity Modelica.Media.Interfaces.PartialMedium.dynamicViscosity

Return dynamic viscosity

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

replaceable partial function dynamicViscosity 
  "Return dynamic viscosity"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output DynamicViscosity eta "Dynamic viscosity";
end dynamicViscosity;

Modelica.Media.Interfaces.PartialMedium.thermalConductivity Modelica.Media.Interfaces.PartialMedium.thermalConductivity

Return thermal conductivity

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

replaceable partial function thermalConductivity 
  "Return thermal conductivity"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output ThermalConductivity lambda "Thermal conductivity";
end thermalConductivity;

Modelica.Media.Interfaces.PartialMedium.prandtlNumber Modelica.Media.Interfaces.PartialMedium.prandtlNumber

Return the Prandtl number

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
PrandtlNumberPrPrandtl number [1]

Modelica definition

replaceable function prandtlNumber "Return the Prandtl number"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output PrandtlNumber Pr "Prandtl number";
algorithm 
  Pr := dynamicViscosity(state)*specificHeatCapacityCp(state)/thermalConductivity(
    state);
end prandtlNumber;

Modelica.Media.Interfaces.PartialMedium.pressure Modelica.Media.Interfaces.PartialMedium.pressure

Return pressure

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

replaceable partial function pressure "Return pressure"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output AbsolutePressure p "Pressure";
end pressure;

Modelica.Media.Interfaces.PartialMedium.temperature Modelica.Media.Interfaces.PartialMedium.temperature

Return temperature

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

replaceable partial function temperature "Return temperature"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output Temperature T "Temperature";
end temperature;

Modelica.Media.Interfaces.PartialMedium.density Modelica.Media.Interfaces.PartialMedium.density

Return density

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

replaceable partial function density "Return density"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output Density d "Density";
end density;

Modelica.Media.Interfaces.PartialMedium.specificEnthalpy Modelica.Media.Interfaces.PartialMedium.specificEnthalpy

Return specific enthalpy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

replaceable partial function specificEnthalpy 
  "Return specific enthalpy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEnthalpy h "Specific enthalpy";
end specificEnthalpy;

Modelica.Media.Interfaces.PartialMedium.specificInternalEnergy Modelica.Media.Interfaces.PartialMedium.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

replaceable partial function specificInternalEnergy 
  "Return specific internal energy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEnergy u "Specific internal energy";
end specificInternalEnergy;

Modelica.Media.Interfaces.PartialMedium.specificEntropy Modelica.Media.Interfaces.PartialMedium.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

replaceable partial function specificEntropy 
  "Return specific entropy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEntropy s "Specific entropy";
end specificEntropy;

Modelica.Media.Interfaces.PartialMedium.specificGibbsEnergy Modelica.Media.Interfaces.PartialMedium.specificGibbsEnergy

Return specific Gibbs energy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

replaceable partial function specificGibbsEnergy 
  "Return specific Gibbs energy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEnergy g "Specific Gibbs energy";
end specificGibbsEnergy;

Modelica.Media.Interfaces.PartialMedium.specificHelmholtzEnergy Modelica.Media.Interfaces.PartialMedium.specificHelmholtzEnergy

Return specific Helmholtz energy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

replaceable partial function specificHelmholtzEnergy 
  "Return specific Helmholtz energy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificEnergy f "Specific Helmholtz energy";
end specificHelmholtzEnergy;

Modelica.Media.Interfaces.PartialMedium.specificHeatCapacityCp Modelica.Media.Interfaces.PartialMedium.specificHeatCapacityCp

Return specific heat capacity at constant pressure

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

replaceable partial function specificHeatCapacityCp 
  "Return specific heat capacity at constant pressure"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificHeatCapacity cp "Specific heat capacity at constant pressure";
end specificHeatCapacityCp;

Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp

alias for deprecated name

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

function heatCapacity_cp = specificHeatCapacityCp "alias for deprecated name";

Modelica.Media.Interfaces.PartialMedium.specificHeatCapacityCv Modelica.Media.Interfaces.PartialMedium.specificHeatCapacityCv

Return specific heat capacity at constant volume

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

replaceable partial function specificHeatCapacityCv 
  "Return specific heat capacity at constant volume"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SpecificHeatCapacity cv "Specific heat capacity at constant volume";
end specificHeatCapacityCv;

Modelica.Media.Interfaces.PartialMedium.heatCapacity_cv Modelica.Media.Interfaces.PartialMedium.heatCapacity_cv

alias for deprecated name

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

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

Modelica definition

function heatCapacity_cv = specificHeatCapacityCv "alias for deprecated name";

Modelica.Media.Interfaces.PartialMedium.isentropicExponent Modelica.Media.Interfaces.PartialMedium.isentropicExponent

Return isentropic exponent

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
IsentropicExponentgammaIsentropic exponent [1]

Modelica definition

replaceable partial function isentropicExponent 
  "Return isentropic exponent"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output IsentropicExponent gamma "Isentropic exponent";
end isentropicExponent;

Modelica.Media.Interfaces.PartialMedium.isentropicEnthalpy Modelica.Media.Interfaces.PartialMedium.isentropicEnthalpy

Return isentropic enthalpy

Information


This function computes an isentropic state transformation:

  1. A medium is in a particular state, refState.
  2. The enhalpy at another state (h_is) shall be computed under the assumption that the state transformation from refState to h_is is performed with a change of specific entropy ds = 0 and the pressure of state h_is is p_downstream and the composition X upstream and downstream is assumed to be the same.

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

Inputs

TypeNameDefaultDescription
AbsolutePressurep_downstream downstream pressure [Pa]
ThermodynamicStaterefState reference state for entropy

Outputs

TypeNameDescription
SpecificEnthalpyh_isIsentropic enthalpy [J/kg]

Modelica definition

replaceable partial function isentropicEnthalpy 
  "Return isentropic enthalpy"
  extends Modelica.Icons.Function;
  input AbsolutePressure p_downstream "downstream pressure";
  input ThermodynamicState refState "reference state for entropy";
  output SpecificEnthalpy h_is "Isentropic enthalpy";
end isentropicEnthalpy;

Modelica.Media.Interfaces.PartialMedium.velocityOfSound Modelica.Media.Interfaces.PartialMedium.velocityOfSound

Return velocity of sound

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
VelocityOfSoundaVelocity of sound [m/s]

Modelica definition

replaceable partial function velocityOfSound 
  "Return velocity of sound"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output VelocityOfSound a "Velocity of sound";
end velocityOfSound;

Modelica.Media.Interfaces.PartialMedium.isobaricExpansionCoefficient Modelica.Media.Interfaces.PartialMedium.isobaricExpansionCoefficient

Return overall the isobaric expansion coefficient beta

Information


beta is defined as  1/v * der(v,T), with v = 1/d, at constant pressure p.

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
IsobaricExpansionCoefficientbetaIsobaric expansion coefficient [1/K]

Modelica definition

replaceable partial function isobaricExpansionCoefficient 
  "Return overall the isobaric expansion coefficient beta"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output IsobaricExpansionCoefficient beta "Isobaric expansion coefficient";
end isobaricExpansionCoefficient;

Modelica.Media.Interfaces.PartialMedium.beta Modelica.Media.Interfaces.PartialMedium.beta

alias for isobaricExpansionCoefficient for user convenience

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
IsobaricExpansionCoefficientbetaIsobaric expansion coefficient [1/K]

Modelica definition

function beta = isobaricExpansionCoefficient 
  "alias for isobaricExpansionCoefficient for user convenience";

Modelica.Media.Interfaces.PartialMedium.isothermalCompressibility Modelica.Media.Interfaces.PartialMedium.isothermalCompressibility

Return overall the isothermal compressibility factor

Information



kappa is defined as - 1/v * der(v,p), with v = 1/d at constant temperature T.

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
IsothermalCompressibilitykappaIsothermal compressibility [1/Pa]

Modelica definition

replaceable partial function isothermalCompressibility 
  "Return overall the isothermal compressibility factor"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output SI.IsothermalCompressibility kappa "Isothermal compressibility";
end isothermalCompressibility;

Modelica.Media.Interfaces.PartialMedium.kappa Modelica.Media.Interfaces.PartialMedium.kappa

alias of isothermalCompressibility for user convenience

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
IsothermalCompressibilitykappaIsothermal compressibility [1/Pa]

Modelica definition

function kappa = isothermalCompressibility 
  "alias of isothermalCompressibility for user convenience";

Modelica.Media.Interfaces.PartialMedium.density_derp_h Modelica.Media.Interfaces.PartialMedium.density_derp_h

Return density derivative w.r.t. pressure at const specific enthalpy

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DerDensityByPressureddphDensity derivative w.r.t. pressure [s2/m2]

Modelica definition

replaceable partial function density_derp_h 
  "Return density derivative w.r.t. pressure at const specific enthalpy"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output DerDensityByPressure ddph "Density derivative w.r.t. pressure";
end density_derp_h;

Modelica.Media.Interfaces.PartialMedium.density_derh_p Modelica.Media.Interfaces.PartialMedium.density_derh_p

Return density derivative w.r.t. specific enthalpy at constant pressure

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DerDensityByEnthalpyddhpDensity derivative w.r.t. specific enthalpy [kg.s2/m5]

Modelica definition

replaceable partial function density_derh_p 
  "Return density derivative w.r.t. specific enthalpy at constant pressure"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output DerDensityByEnthalpy ddhp 
    "Density derivative w.r.t. specific enthalpy";
end density_derh_p;

Modelica.Media.Interfaces.PartialMedium.density_derp_T Modelica.Media.Interfaces.PartialMedium.density_derp_T

Return density derivative w.r.t. pressure at const temperature

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DerDensityByPressureddpTDensity derivative w.r.t. pressure [s2/m2]

Modelica definition

replaceable partial function density_derp_T 
  "Return density derivative w.r.t. pressure at const temperature"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output DerDensityByPressure ddpT "Density derivative w.r.t. pressure";
end density_derp_T;

Modelica.Media.Interfaces.PartialMedium.density_derT_p Modelica.Media.Interfaces.PartialMedium.density_derT_p

Return density derivative w.r.t. temperature at constant pressure

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DerDensityByTemperatureddTpDensity derivative w.r.t. temperature [kg/(m3.K)]

Modelica definition

replaceable partial function density_derT_p 
  "Return density derivative w.r.t. temperature at constant pressure"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output DerDensityByTemperature ddTp "Density derivative w.r.t. temperature";
end density_derT_p;

Modelica.Media.Interfaces.PartialMedium.density_derX Modelica.Media.Interfaces.PartialMedium.density_derX

Return density derivative w.r.t. mass fraction

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
DensitydddX[nX]Derivative of density w.r.t. mass fraction [kg/m3]

Modelica definition

replaceable partial function density_derX 
  "Return density derivative w.r.t. mass fraction"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output Density[nX] dddX "Derivative of density w.r.t. mass fraction";
end density_derX;

Modelica.Media.Interfaces.PartialMedium.molarMass Modelica.Media.Interfaces.PartialMedium.molarMass

Return the molar mass of the medium

Information

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate thermodynamic state record

Outputs

TypeNameDescription
MolarMassMMMixture molar mass [kg/mol]

Modelica definition

replaceable partial function molarMass 
  "Return the molar mass of the medium"
  extends Modelica.Icons.Function;
  input ThermodynamicState state "thermodynamic state record";
  output MolarMass MM "Mixture molar mass";
end molarMass;

Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_pTX Modelica.Media.Interfaces.PartialMedium.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[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

replaceable 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[:]=reference_X "Mass fractions";
  output SpecificEnthalpy h "Specific enthalpy";
algorithm 
  h := specificEnthalpy(setState_pTX(p,T,X));
end specificEnthalpy_pTX;

Modelica.Media.Interfaces.PartialMedium.specificEntropy_pTX Modelica.Media.Interfaces.PartialMedium.specificEntropy_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[:]reference_XMass fractions [kg/kg]

Outputs

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

Modelica definition

replaceable function specificEntropy_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[:]=reference_X "Mass fractions";
  output SpecificEntropy s "Specific entropy";
algorithm 
  s := specificEntropy(setState_pTX(p,T,X));

end specificEntropy_pTX;

Modelica.Media.Interfaces.PartialMedium.density_pTX Modelica.Media.Interfaces.PartialMedium.density_pTX

Return density 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[:] Mass fractions [kg/kg]

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

replaceable function density_pTX 
  "Return density from p, T, and X or Xi"
  extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input Temperature T "Temperature";
  input MassFraction X[:] "Mass fractions";
  output Density d "Density";
algorithm 
  d := density(setState_pTX(p,T,X));
end density_pTX;

Modelica.Media.Interfaces.PartialMedium.temperature_phX Modelica.Media.Interfaces.PartialMedium.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[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

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

Modelica.Media.Interfaces.PartialMedium.density_phX Modelica.Media.Interfaces.PartialMedium.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[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

replaceable 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[:]=reference_X "Mass fractions";
  output Density d "Density";
algorithm 
  d := density(setState_phX(p,h,X));
end density_phX;

Modelica.Media.Interfaces.PartialMedium.temperature_psX Modelica.Media.Interfaces.PartialMedium.temperature_psX

Return temperature 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
TemperatureTTemperature [K]

Modelica definition

replaceable function temperature_psX 
  "Return temperature 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 Temperature T "Temperature";
algorithm 
  T := temperature(setState_psX(p,s,X));
end temperature_psX;

Modelica.Media.Interfaces.PartialMedium.density_psX Modelica.Media.Interfaces.PartialMedium.density_psX

Return density 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
DensitydDensity [kg/m3]

Modelica definition

replaceable function density_psX 
  "Return density 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 Density d "Density";
algorithm 
  d := density(setState_psX(p,s,X));
end density_psX;

Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_psX Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_psX

Return specific enthalpy 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
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

replaceable function specificEnthalpy_psX 
  "Return specific enthalpy 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 SpecificEnthalpy h "Specific enthalpy";
algorithm 
  h := specificEnthalpy(setState_psX(p,s,X));
end specificEnthalpy_psX;

Modelica.Media.Interfaces.PartialMedium.BaseProperties.InputAbsolutePressure

Pressure as input signal connector

Modelica definition

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

Modelica.Media.Interfaces.PartialMedium.BaseProperties.InputMassFraction

Mass fraction as input signal connector

Modelica definition

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

Modelica.Media.Interfaces.PartialMedium.BaseProperties.InputSpecificEnthalpy

Specific enthalpy as input signal connector

Modelica definition

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

Automatically generated Fri Nov 12 16:31:29 2010.