This model calculates medium properties for water in the liquid, gas and two phase regions according to the IAPWS/IF97 standard, i.e., the accepted industrial standard and best compromise between accuracy and computation time. For more details see Modelica.Media.Water.IF97_Utilities. Three variable pairs can be the independent variables of the model:
| Variable | Unit | Description | 
| T | K | temperature | 
| u | J/kg | specific internal energy | 
| d | kg/m^3 | density | 
| p | Pa | pressure | 
| h | J/kg | specific enthalpy | 
In some cases additional medium properties are needed. A component that needs these optional properties has to call one of the functions listed in Modelica.Media.UsersGuide.MediumUsage.OptionalProperties and in Modelica.Media.UsersGuide.MediumUsage.TwoPhase.
Many further properties can be computed. Using the well-known Bridgman's Tables, all first partial derivatives of the standard thermodynamic variables can be computed easily.
Extends from Interfaces.PartialTwoPhaseMedium (Base class for two phase medium of one substance).
| Name | Description | 
|---|---|
|  ThermodynamicState | thermodynamic state | 
| ph_explicit | true if explicit in pressure and specific enthalpy | 
| dT_explicit | true if explicit in density and temperature | 
| pT_explicit | true if explicit in pressure and temperature | 
|  BaseProperties | Base properties of water | 
|  density_ph | Computes density as a function of pressure and specific enthalpy | 
|  temperature_ph | Computes temperature as a function of pressure and specific enthalpy | 
|  temperature_ps | Compute temperature from pressure and specific enthalpy | 
|  density_ps | Computes density as a function of pressure and specific enthalpy | 
|  pressure_dT | Computes pressure as a function of density and temperature | 
|  specificEnthalpy_dT | Computes specific enthalpy as a function of density and temperature | 
|  specificEnthalpy_pT | Computes specific enthalpy as a function of pressure and temperature | 
|  specificEnthalpy_ps | Computes specific enthalpy as a function of pressure and temperature | 
|  density_pT | Computes density as a function of pressure and temperature | 
|  setDewState | set the thermodynamic state on the dew line | 
|  setBubbleState | set the thermodynamic state on the bubble line | 
|  dynamicViscosity | Dynamic viscosity of water | 
|  thermalConductivity | Thermal conductivity of water | 
|  surfaceTension | Surface tension in two phase region of water | 
|  pressure | return pressure of ideal gas | 
|  temperature | return temperature of ideal gas | 
|  density | return density of ideal gas | 
|  specificEnthalpy | Return specific enthalpy | 
|  specificInternalEnergy | Return specific internal energy | 
|  specificGibbsEnergy | Return specific Gibbs energy | 
|  specificHelmholtzEnergy | Return specific Helmholtz energy | 
|  specificEntropy | specific entropy of water | 
|  specificHeatCapacityCp | specific heat capacity at constant pressure of water | 
|  specificHeatCapacityCv | specific heat capacity at constant volume of water | 
|  isentropicExponent | Return isentropic exponent | 
|  isothermalCompressibility | Isothermal compressibility of water | 
|  isobaricExpansionCoefficient | isobaric expansion coefficient of water | 
|  velocityOfSound | Return velocity of sound as a function of the thermodynamic state record | 
|  isentropicEnthalpy | compute h(p,s) | 
|  density_derh_p | density derivative by specific enthalpy | 
|  density_derp_h | density derivative by pressure | 
|  bubbleEnthalpy | boiling curve specific enthalpy of water | 
|  dewEnthalpy | dew curve specific enthalpy of water | 
|  bubbleEntropy | boiling curve specific entropy of water | 
|  dewEntropy | dew curve specific entropy of water | 
|  bubbleDensity | boiling curve specific density of water | 
|  dewDensity | dew curve specific density of water | 
|  saturationTemperature | saturation temperature of water | 
|  saturationTemperature_derp | derivative of saturation temperature w.r.t. pressure | 
|  saturationPressure | saturation pressure of water | 
|  dBubbleDensity_dPressure | bubble point density derivative | 
|  dDewDensity_dPressure | dew point density derivative | 
|  dBubbleEnthalpy_dPressure | bubble point specific enthalpy derivative | 
|  dDewEnthalpy_dPressure | dew point specific enthalpy derivative | 
|  setState_dTX | Return thermodynamic state of water as function of d and T | 
|  setState_phX | Return thermodynamic state of water as function of p and h | 
|  setState_psX | Return thermodynamic state of water as function of p and s | 
|  setState_pTX | Return thermodynamic state of water as function of p and T | 
|  setSmoothState | Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b | 
| Inherited | |
| smoothModel | true if the (derived) model should not generate state events | 
| onePhase | true if the (derived) model should never be called with two-phase inputs | 
|  FluidLimits | validity limits for fluid model | 
|  FluidConstants | extended fluid constants | 
| fluidConstants | constant data for the fluid | 
|  SaturationProperties | Saturation properties of two phase medium | 
| FixedPhase | phase of the fluid: 1 for 1-phase, 2 for two-phase, 0 for not known, e.g. interactive use | 
|  setSat_T | Return saturation property record from temperature | 
|  setSat_p | Return saturation property record from pressure | 
|  saturationPressure_sat | Return saturation temperature | 
|  saturationTemperature_sat | Return saturation temperature | 
|  saturationTemperature_derp_sat | Return derivative of saturation temperature w.r.t. pressure | 
|  molarMass | Return the molar mass of the medium | 
|  specificEnthalpy_pTX | Return specific enthalpy from pressure, temperature and mass fraction | 
|  temperature_phX | Return temperature from p, h, and X or Xi | 
|  density_phX | Return density from p, h, 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 | 
|  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 | 
|  setState_px | Return thermodynamic state from pressure and vapour quality | 
|  setState_Tx | Return thermodynamic state from temperature and vapour quality | 
|  vapourQuality | Return vapour quality | 
| 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 | 
|  prandtlNumber | Return the Prandtl number | 
|  heatCapacity_cp | alias for deprecated name | 
|  heatCapacity_cv | alias for deprecated name | 
|  beta | alias for isobaricExpansionCoefficient for user convenience | 
|  kappa | alias of isothermalCompressibility for user convenience | 
|  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 | 
|  density_pTX | Return density from p, T, 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 Boolean ph_explicit "true if explicit in pressure and specific enthalpy";
constant Boolean dT_explicit "true if explicit in density and temperature";
constant Boolean pT_explicit "true if explicit in pressure and temperature";
 Modelica.Media.Water.WaterIF97_base.ThermodynamicState
Modelica.Media.Water.WaterIF97_base.ThermodynamicState
redeclare record extends ThermodynamicState "thermodynamic state" SpecificEnthalpy h "specific enthalpy"; Density d "density"; Temperature T "temperature"; AbsolutePressure p "pressure"; end ThermodynamicState;
 Modelica.Media.Water.WaterIF97_base.BaseProperties
Modelica.Media.Water.WaterIF97_base.BaseProperties
| Type | Name | Default | Description | 
|---|---|---|---|
| Advanced | |||
| Boolean | preferredMediumStates | false | = true if StateSelect.prefer shall be used for the independent property variables of the medium | 
redeclare replaceable model extends BaseProperties(
  h(stateSelect=if ph_explicit and preferredMediumStates then StateSelect.prefer else StateSelect.default),
  d(stateSelect=if dT_explicit and preferredMediumStates then StateSelect.prefer else StateSelect.default),
  T(stateSelect=if (pT_explicit or dT_explicit) and preferredMediumStates then StateSelect.prefer else StateSelect.default),
  p(stateSelect=if (pT_explicit or ph_explicit) and preferredMediumStates then StateSelect.prefer else StateSelect.default)) 
  "Base properties of water"
  Integer phase(min=0, max=2, start=1,fixed=false) 
    "2 for two-phase, 1 for one-phase, 0 if not known";
  SaturationProperties sat(Tsat(start=300.0), psat(start=1.0e5)) 
    "saturation temperature and pressure";
equation 
  MM = fluidConstants[1].molarMass;
  if smoothModel then
    if onePhase then
      phase = 1;
      if ph_explicit then
        assert(((h < bubbleEnthalpy(sat) or h > dewEnthalpy(sat)) or p >
  fluidConstants[1].criticalPressure),
 "With onePhase=true this model may only be called with one-phase states h < hl or h > hv!"
 + "(p = " + String(p) + ", h = " + String(h) + ")");
      else
 if dT_explicit then
   assert(not ((d < bubbleDensity(sat) and d > dewDensity(sat)) and T <
 fluidConstants[1].criticalTemperature),
   "With onePhase=true this model may only be called with one-phase states d > dl or d < dv!"
   + "(d = " + String(d) + ", T = " + String(T) + ")");
 else
   assert(true,"no events for pT-model");
 end if;
      end if;
    else
      phase = 0;
    end if;
  else
    if ph_explicit then
      phase = if ((h < bubbleEnthalpy(sat) or h > dewEnthalpy(sat)) or p >
        fluidConstants[1].criticalPressure) then 1 else 2;
    elseif dT_explicit then
      phase = if not ((d < bubbleDensity(sat) and d > dewDensity(sat)) and T
         < fluidConstants[1].criticalTemperature) then 1 else 2;
    else
      phase = 1;
      //this is for the one-phase only case pT
    end if;
  end if;
  if dT_explicit then
    p = pressure_dT(d, T, phase);
    h = specificEnthalpy_dT(d, T, phase);
    sat.Tsat = T;
    sat.psat = saturationPressure(T);
  elseif ph_explicit then
    d = density_ph(p, h, phase);
    T = temperature_ph(p, h, phase);
    sat.Tsat = saturationTemperature(p);
    sat.psat = p;
  else
    h = specificEnthalpy_pT(p, T);
    d = density_pT(p, T);
    sat.psat = p;
    sat.Tsat = saturationTemperature(p);
  end if;
  u = h - p/d;
  R = Modelica.Constants.R/fluidConstants[1].molarMass;
  h = state.h;
  p = state.p;
  T = state.T;
  d = state.d;
  phase = state.phase;
end BaseProperties;
 Modelica.Media.Water.WaterIF97_base.density_ph
Modelica.Media.Water.WaterIF97_base.density_ph
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| Density | d | Density [kg/m3] | 
redeclare function density_ph "Computes density as a function of pressure and specific enthalpy" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input SpecificEnthalpy h "Specific enthalpy"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output Density d "Density"; algorithm d := IF97_Utilities.rho_ph(p, h, phase); end density_ph;
 Modelica.Media.Water.WaterIF97_base.temperature_ph
Modelica.Media.Water.WaterIF97_base.temperature_ph
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | Temperature [K] | 
redeclare function temperature_ph "Computes temperature as a function of pressure and specific enthalpy" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input SpecificEnthalpy h "Specific enthalpy"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output Temperature T "Temperature"; algorithm T := IF97_Utilities.T_ph(p, h, phase); end temperature_ph;
 Modelica.Media.Water.WaterIF97_base.temperature_ps
Modelica.Media.Water.WaterIF97_base.temperature_ps
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | Temperature [K] | 
redeclare function temperature_ps "Compute temperature from pressure and specific enthalpy" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input SpecificEntropy s "Specific entropy"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output Temperature T "Temperature"; algorithm T := IF97_Utilities.T_ps(p, s, phase); end temperature_ps;
 Modelica.Media.Water.WaterIF97_base.density_ps
Modelica.Media.Water.WaterIF97_base.density_ps
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| Density | d | density [kg/m3] | 
redeclare function density_ps 
  "Computes density as a function of pressure and specific enthalpy"
    extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEntropy s "Specific entropy";
  input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known";
  output Density d "density";
algorithm 
  d := IF97_Utilities.rho_ps(p, s, phase);
end density_ps;
 Modelica.Media.Water.WaterIF97_base.pressure_dT
Modelica.Media.Water.WaterIF97_base.pressure_dT
| Type | Name | Default | Description | 
|---|---|---|---|
| Density | d | Density [kg/m3] | |
| Temperature | T | Temperature [K] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | 
redeclare function pressure_dT "Computes pressure as a function of density and temperature" extends Modelica.Icons.Function; input Density d "Density"; input Temperature T "Temperature"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output AbsolutePressure p "Pressure"; algorithm p := IF97_Utilities.p_dT(d, T, phase); end pressure_dT;
 Modelica.Media.Water.WaterIF97_base.specificEnthalpy_dT
Modelica.Media.Water.WaterIF97_base.specificEnthalpy_dT
| Type | Name | Default | Description | 
|---|---|---|---|
| Density | d | Density [kg/m3] | |
| Temperature | T | Temperature [K] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | specific enthalpy [J/kg] | 
redeclare function specificEnthalpy_dT "Computes specific enthalpy as a function of density and temperature" extends Modelica.Icons.Function; input Density d "Density"; input Temperature T "Temperature"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output SpecificEnthalpy h "specific enthalpy"; algorithm h := IF97_Utilities.h_dT(d, T, phase); end specificEnthalpy_dT;
 Modelica.Media.Water.WaterIF97_base.specificEnthalpy_pT
Modelica.Media.Water.WaterIF97_base.specificEnthalpy_pT
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| Temperature | T | Temperature [K] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | specific enthalpy [J/kg] | 
redeclare function specificEnthalpy_pT "Computes specific enthalpy as a function of pressure and temperature" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input Temperature T "Temperature"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output SpecificEnthalpy h "specific enthalpy"; algorithm h := IF97_Utilities.h_pT(p, T); end specificEnthalpy_pT;
 Modelica.Media.Water.WaterIF97_base.specificEnthalpy_ps
Modelica.Media.Water.WaterIF97_base.specificEnthalpy_ps
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | specific enthalpy [J/kg] | 
redeclare function specificEnthalpy_ps 
  "Computes specific enthalpy as a function of pressure and temperature"
    extends Modelica.Icons.Function;
  input AbsolutePressure p "Pressure";
  input SpecificEntropy s "Specific entropy";
  input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known";
  output SpecificEnthalpy h "specific enthalpy";
algorithm 
  h := IF97_Utilities.h_ps(p, s, phase);
end specificEnthalpy_ps;
 Modelica.Media.Water.WaterIF97_base.density_pT
Modelica.Media.Water.WaterIF97_base.density_pT
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | |
| Temperature | T | Temperature [K] | |
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Type | Name | Description | 
|---|---|---|
| Density | d | Density [kg/m3] | 
redeclare function density_pT "Computes density as a function of pressure and temperature" extends Modelica.Icons.Function; input AbsolutePressure p "Pressure"; input Temperature T "Temperature"; input FixedPhase phase=0 "2 for two-phase, 1 for one-phase, 0 if not known"; output Density d "Density"; algorithm d := IF97_Utilities.rho_pT(p, T); end density_pT;
 Modelica.Media.Water.WaterIF97_base.setDewState
Modelica.Media.Water.WaterIF97_base.setDewState
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation point | |
| FixedPhase | phase | 1 | phase: default is one phase | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | complete thermodynamic state info | 
redeclare function extends setDewState 
  "set the thermodynamic state on the dew line"
algorithm 
  state := ThermodynamicState(
     phase=  phase,
     p=  sat.psat,
     T=  sat.Tsat,
     h=  dewEnthalpy(sat),
     d=  dewDensity(sat));
end setDewState;
 Modelica.Media.Water.WaterIF97_base.setBubbleState
Modelica.Media.Water.WaterIF97_base.setBubbleState
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation point | |
| FixedPhase | phase | 1 | phase: default is one phase | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | complete thermodynamic state info | 
redeclare function extends setBubbleState 
  "set the thermodynamic state on the bubble line"
algorithm 
  state := ThermodynamicState(
     phase=  phase,
     p=  sat.psat,
     T=  sat.Tsat,
     h=  bubbleEnthalpy(sat),
     d=  bubbleDensity(sat));
end setBubbleState;
 Modelica.Media.Water.WaterIF97_base.dynamicViscosity
Modelica.Media.Water.WaterIF97_base.dynamicViscosity
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| DynamicViscosity | eta | Dynamic viscosity [Pa.s] | 
redeclare function extends dynamicViscosity 
  "Dynamic viscosity of water"
algorithm 
  eta := IF97_Utilities.dynamicViscosity(state.d, state.T, state.p, state.
    phase);
end dynamicViscosity;
 Modelica.Media.Water.WaterIF97_base.thermalConductivity
Modelica.Media.Water.WaterIF97_base.thermalConductivity
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| ThermalConductivity | lambda | Thermal conductivity [W/(m.K)] | 
redeclare function extends thermalConductivity 
  "Thermal conductivity of water"
algorithm 
  lambda := IF97_Utilities.thermalConductivity(state.d, state.T, state.p,
    state.phase);
end thermalConductivity;
 Modelica.Media.Water.WaterIF97_base.surfaceTension
Modelica.Media.Water.WaterIF97_base.surfaceTension
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| SurfaceTension | sigma | Surface tension sigma in the two phase region [N/m] | 
redeclare function extends surfaceTension "Surface tension in two phase region of water" algorithm sigma := IF97_Utilities.surfaceTension(sat.Tsat); end surfaceTension;
 Modelica.Media.Water.WaterIF97_base.pressure
Modelica.Media.Water.WaterIF97_base.pressure
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | p | Pressure [Pa] | 
redeclare function extends pressure "return pressure of ideal gas" algorithm p := state.p; end pressure;
 Modelica.Media.Water.WaterIF97_base.temperature
Modelica.Media.Water.WaterIF97_base.temperature
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | Temperature [K] | 
redeclare function extends temperature "return temperature of ideal gas" algorithm T := state.T; end temperature;
 Modelica.Media.Water.WaterIF97_base.density
Modelica.Media.Water.WaterIF97_base.density
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| Density | d | Density [kg/m3] | 
redeclare function extends density "return density of ideal gas" algorithm d := state.d; end density;
 Modelica.Media.Water.WaterIF97_base.specificEnthalpy
Modelica.Media.Water.WaterIF97_base.specificEnthalpy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | 
redeclare function extends specificEnthalpy "Return specific enthalpy" extends Modelica.Icons.Function; algorithm h := state.h; end specificEnthalpy;
 Modelica.Media.Water.WaterIF97_base.specificInternalEnergy
Modelica.Media.Water.WaterIF97_base.specificInternalEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | u | Specific internal energy [J/kg] | 
redeclare function extends specificInternalEnergy "Return specific internal energy" extends Modelica.Icons.Function; algorithm u := state.h - state.p/state.d; end specificInternalEnergy;
 Modelica.Media.Water.WaterIF97_base.specificGibbsEnergy
Modelica.Media.Water.WaterIF97_base.specificGibbsEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | g | Specific Gibbs energy [J/kg] | 
redeclare function extends specificGibbsEnergy "Return specific Gibbs energy" extends Modelica.Icons.Function; algorithm g := state.h - state.T*specificEntropy(state); end specificGibbsEnergy;
 Modelica.Media.Water.WaterIF97_base.specificHelmholtzEnergy
Modelica.Media.Water.WaterIF97_base.specificHelmholtzEnergy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnergy | f | Specific Helmholtz energy [J/kg] | 
redeclare function extends specificHelmholtzEnergy "Return specific Helmholtz energy" extends Modelica.Icons.Function; algorithm f := state.h - state.p/state.d - state.T*specificEntropy(state); end specificHelmholtzEnergy;
 Modelica.Media.Water.WaterIF97_base.specificEntropy
Modelica.Media.Water.WaterIF97_base.specificEntropy
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | 
redeclare function extends specificEntropy 
  "specific entropy of water"
algorithm 
  if dT_explicit then
    s := IF97_Utilities.s_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    s := IF97_Utilities.s_pT(state.p, state.T);
  else
    s := IF97_Utilities.s_ph(state.p, state.h, state.phase);
  end if;
end specificEntropy;
 Modelica.Media.Water.WaterIF97_base.specificHeatCapacityCp
Modelica.Media.Water.WaterIF97_base.specificHeatCapacityCp
In the two phase region this function returns the interpolated heat capacity between the liquid and vapour state heat capacities.
Error:Found no end-tag in HTML-documentationExtends from (Return specific heat capacity at constant pressure).
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | cp | Specific heat capacity at constant pressure [J/(kg.K)] | 
redeclare function extends specificHeatCapacityCp 
  "specific heat capacity at constant pressure of water"
algorithm 
  if dT_explicit then
    cp := IF97_Utilities.cp_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    cp := IF97_Utilities.cp_pT(state.p, state.T);
  else
    cp := IF97_Utilities.cp_ph(state.p, state.h, state.phase);
  end if;
end specificHeatCapacityCp;
 
 Modelica.Media.Water.WaterIF97_base.specificHeatCapacityCv
Modelica.Media.Water.WaterIF97_base.specificHeatCapacityCv
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| SpecificHeatCapacity | cv | Specific heat capacity at constant volume [J/(kg.K)] | 
redeclare function extends specificHeatCapacityCv 
  "specific heat capacity at constant volume of water"
algorithm 
  if dT_explicit then
    cv := IF97_Utilities.cv_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    cv := IF97_Utilities.cv_pT(state.p, state.T);
  else
    cv := IF97_Utilities.cv_ph(state.p, state.h, state.phase);
  end if;
end specificHeatCapacityCv;
 Modelica.Media.Water.WaterIF97_base.isentropicExponent
Modelica.Media.Water.WaterIF97_base.isentropicExponent
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| IsentropicExponent | gamma | Isentropic exponent [1] | 
redeclare function extends isentropicExponent 
  "Return isentropic exponent"
algorithm 
  if dT_explicit then
    gamma := IF97_Utilities.isentropicExponent_dT(state.d, state.T, state.
      phase);
  elseif pT_explicit then
    gamma := IF97_Utilities.isentropicExponent_pT(state.p, state.T);
  else
    gamma := IF97_Utilities.isentropicExponent_ph(state.p, state.h, state.
      phase);
  end if;
end isentropicExponent;
 Modelica.Media.Water.WaterIF97_base.isothermalCompressibility
Modelica.Media.Water.WaterIF97_base.isothermalCompressibility
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| IsothermalCompressibility | kappa | Isothermal compressibility [1/Pa] | 
redeclare function extends isothermalCompressibility 
  "Isothermal compressibility of water"
algorithm 
  if dT_explicit then
    kappa := IF97_Utilities.kappa_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    kappa := IF97_Utilities.kappa_pT(state.p, state.T);
  else
    kappa := IF97_Utilities.kappa_ph(state.p, state.h, state.phase);
  end if;
end isothermalCompressibility;
 Modelica.Media.Water.WaterIF97_base.isobaricExpansionCoefficient
Modelica.Media.Water.WaterIF97_base.isobaricExpansionCoefficient
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| IsobaricExpansionCoefficient | beta | Isobaric expansion coefficient [1/K] | 
redeclare function extends isobaricExpansionCoefficient 
  "isobaric expansion coefficient of water"
algorithm 
  if dT_explicit then
    beta := IF97_Utilities.beta_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    beta := IF97_Utilities.beta_pT(state.p, state.T);
  else
    beta := IF97_Utilities.beta_ph(state.p, state.h, state.phase);
  end if;
end isobaricExpansionCoefficient;
 Modelica.Media.Water.WaterIF97_base.velocityOfSound
Modelica.Media.Water.WaterIF97_base.velocityOfSound
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| VelocityOfSound | a | Velocity of sound [m/s] | 
redeclare function extends velocityOfSound 
  "Return velocity of sound as a function of the thermodynamic state record"
algorithm 
  if dT_explicit then
    a := IF97_Utilities.velocityOfSound_dT(state.d, state.T, state.phase);
  elseif pT_explicit then
    a := IF97_Utilities.velocityOfSound_pT(state.p, state.T);
  else
    a := IF97_Utilities.velocityOfSound_ph(state.p, state.h, state.phase);
  end if;
end velocityOfSound;
 Modelica.Media.Water.WaterIF97_base.isentropicEnthalpy
Modelica.Media.Water.WaterIF97_base.isentropicEnthalpy
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p_downstream | downstream pressure [Pa] | |
| ThermodynamicState | refState | reference state for entropy | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | h_is | Isentropic enthalpy [J/kg] | 
redeclare function extends isentropicEnthalpy "compute h(p,s)"
algorithm 
  h_is := IF97_Utilities.isentropicEnthalpy(p_downstream, specificEntropy(
    refState), 0);
end isentropicEnthalpy;
 Modelica.Media.Water.WaterIF97_base.density_derh_p
Modelica.Media.Water.WaterIF97_base.density_derh_p
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| DerDensityByEnthalpy | ddhp | Density derivative wrt specific enthalpy [kg.s2/m5] | 
redeclare function extends density_derh_p "density derivative by specific enthalpy" algorithm ddhp := IF97_Utilities.ddhp(state.p, state.h, state.phase); end density_derh_p;
 Modelica.Media.Water.WaterIF97_base.density_derp_h
Modelica.Media.Water.WaterIF97_base.density_derp_h
| Type | Name | Default | Description | 
|---|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
| Type | Name | Description | 
|---|---|---|
| DerDensityByPressure | ddph | Density derivative wrt pressure [s2/m2] | 
redeclare function extends density_derp_h "density derivative by pressure" algorithm ddph := IF97_Utilities.ddph(state.p, state.h, state.phase); end density_derp_h;
 Modelica.Media.Water.WaterIF97_base.bubbleEnthalpy
Modelica.Media.Water.WaterIF97_base.bubbleEnthalpy
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | hl | boiling curve specific enthalpy [J/kg] | 
redeclare function extends bubbleEnthalpy "boiling curve specific enthalpy of water" algorithm hl := IF97_Utilities.BaseIF97.Regions.hl_p(sat.psat); end bubbleEnthalpy;
 Modelica.Media.Water.WaterIF97_base.dewEnthalpy
Modelica.Media.Water.WaterIF97_base.dewEnthalpy
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEnthalpy | hv | dew curve specific enthalpy [J/kg] | 
redeclare function extends dewEnthalpy "dew curve specific enthalpy of water" algorithm hv := IF97_Utilities.BaseIF97.Regions.hv_p(sat.psat); end dewEnthalpy;
 Modelica.Media.Water.WaterIF97_base.bubbleEntropy
Modelica.Media.Water.WaterIF97_base.bubbleEntropy
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEntropy | sl | boiling curve specific entropy [J/(kg.K)] | 
redeclare function extends bubbleEntropy "boiling curve specific entropy of water" algorithm sl := IF97_Utilities.BaseIF97.Regions.sl_p(sat.psat); end bubbleEntropy;
 Modelica.Media.Water.WaterIF97_base.dewEntropy
Modelica.Media.Water.WaterIF97_base.dewEntropy
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| SpecificEntropy | sv | dew curve specific entropy [J/(kg.K)] | 
redeclare function extends dewEntropy "dew curve specific entropy of water" algorithm sv := IF97_Utilities.BaseIF97.Regions.sv_p(sat.psat); end dewEntropy;
 Modelica.Media.Water.WaterIF97_base.bubbleDensity
Modelica.Media.Water.WaterIF97_base.bubbleDensity
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| Density | dl | boiling curve density [kg/m3] | 
redeclare function extends bubbleDensity 
  "boiling curve specific density of water"
algorithm 
  if ph_explicit or pT_explicit then
    dl := IF97_Utilities.BaseIF97.Regions.rhol_p(sat.psat);
  else
    dl := IF97_Utilities.BaseIF97.Regions.rhol_T(sat.Tsat);
  end if;
end bubbleDensity;
 Modelica.Media.Water.WaterIF97_base.dewDensity
Modelica.Media.Water.WaterIF97_base.dewDensity
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| Density | dv | dew curve density [kg/m3] | 
redeclare function extends dewDensity 
  "dew curve specific density of water"
algorithm 
  if ph_explicit or pT_explicit then
    dv := IF97_Utilities.BaseIF97.Regions.rhov_p(sat.psat);
  else
    dv := IF97_Utilities.BaseIF97.Regions.rhov_T(sat.Tsat);
  end if;
end dewDensity;
 Modelica.Media.Water.WaterIF97_base.saturationTemperature
Modelica.Media.Water.WaterIF97_base.saturationTemperature
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | pressure [Pa] | 
| Type | Name | Description | 
|---|---|---|
| Temperature | T | saturation temperature [K] | 
redeclare function extends saturationTemperature "saturation temperature of water" algorithm T := IF97_Utilities.BaseIF97.Basic.tsat(p); end saturationTemperature;
 Modelica.Media.Water.WaterIF97_base.saturationTemperature_derp
Modelica.Media.Water.WaterIF97_base.saturationTemperature_derp
| Type | Name | Default | Description | 
|---|---|---|---|
| AbsolutePressure | p | pressure [Pa] | 
| Type | Name | Description | 
|---|---|---|
| Real | dTp | derivative of saturation temperature w.r.t. pressure | 
redeclare function extends saturationTemperature_derp "derivative of saturation temperature w.r.t. pressure" algorithm dTp := IF97_Utilities.BaseIF97.Basic.dtsatofp(p); end saturationTemperature_derp;
 Modelica.Media.Water.WaterIF97_base.saturationPressure
Modelica.Media.Water.WaterIF97_base.saturationPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| Temperature | T | temperature [K] | 
| Type | Name | Description | 
|---|---|---|
| AbsolutePressure | p | saturation pressure [Pa] | 
redeclare function extends saturationPressure "saturation pressure of water" algorithm p := IF97_Utilities.BaseIF97.Basic.psat(T); end saturationPressure;
 Modelica.Media.Water.WaterIF97_base.dBubbleDensity_dPressure
Modelica.Media.Water.WaterIF97_base.dBubbleDensity_dPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| DerDensityByPressure | ddldp | boiling curve density derivative [s2/m2] | 
redeclare function extends dBubbleDensity_dPressure "bubble point density derivative" algorithm ddldp := IF97_Utilities.BaseIF97.Regions.drhol_dp(sat.psat); end dBubbleDensity_dPressure;
 Modelica.Media.Water.WaterIF97_base.dDewDensity_dPressure
Modelica.Media.Water.WaterIF97_base.dDewDensity_dPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| DerDensityByPressure | ddvdp | saturated steam density derivative [s2/m2] | 
redeclare function extends dDewDensity_dPressure "dew point density derivative" algorithm ddvdp := IF97_Utilities.BaseIF97.Regions.drhov_dp(sat.psat); end dDewDensity_dPressure;
 Modelica.Media.Water.WaterIF97_base.dBubbleEnthalpy_dPressure
Modelica.Media.Water.WaterIF97_base.dBubbleEnthalpy_dPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| DerEnthalpyByPressure | dhldp | boiling curve specific enthalpy derivative [J.m.s2/kg2] | 
redeclare function extends dBubbleEnthalpy_dPressure "bubble point specific enthalpy derivative" algorithm dhldp := IF97_Utilities.BaseIF97.Regions.dhl_dp(sat.psat); end dBubbleEnthalpy_dPressure;
 Modelica.Media.Water.WaterIF97_base.dDewEnthalpy_dPressure
Modelica.Media.Water.WaterIF97_base.dDewEnthalpy_dPressure
| Type | Name | Default | Description | 
|---|---|---|---|
| SaturationProperties | sat | saturation property record | 
| Type | Name | Description | 
|---|---|---|
| DerEnthalpyByPressure | dhvdp | saturated steam specific enthalpy derivative [J.m.s2/kg2] | 
redeclare function extends dDewEnthalpy_dPressure "dew point specific enthalpy derivative" algorithm dhvdp := IF97_Utilities.BaseIF97.Regions.dhv_dp(sat.psat); end dDewEnthalpy_dPressure;
 Modelica.Media.Water.WaterIF97_base.setState_dTX
Modelica.Media.Water.WaterIF97_base.setState_dTX
| Type | Name | Default | Description | 
|---|---|---|---|
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| Density | d | density [kg/m3] | |
| Temperature | T | Temperature [K] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
redeclare function extends setState_dTX 
  "Return thermodynamic state of water as function of d and T"
algorithm 
  state := ThermodynamicState(
    d=d,
    T=T,
    phase=0,
    h=specificEnthalpy_dT(d,T),
    p=pressure_dT(d,T));
end setState_dTX;
 Modelica.Media.Water.WaterIF97_base.setState_phX
Modelica.Media.Water.WaterIF97_base.setState_phX
| Type | Name | Default | Description | 
|---|---|---|---|
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEnthalpy | h | Specific enthalpy [J/kg] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
redeclare function extends setState_phX 
  "Return thermodynamic state of water as function of p and h"
algorithm 
  state := ThermodynamicState(
    d=density_ph(p,h),
    T=temperature_ph(p,h),
    phase=0,
    h=h,
    p=p);
end setState_phX;
 Modelica.Media.Water.WaterIF97_base.setState_psX
Modelica.Media.Water.WaterIF97_base.setState_psX
| Type | Name | Default | Description | 
|---|---|---|---|
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| AbsolutePressure | p | Pressure [Pa] | |
| SpecificEntropy | s | Specific entropy [J/(kg.K)] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
redeclare function extends setState_psX 
  "Return thermodynamic state of water as function of p and s"
algorithm 
  state := ThermodynamicState(
    d=density_ps(p,s),
    T=temperature_ps(p,s),
    phase=0,
    h=specificEnthalpy_ps(p,s),
    p=p);
end setState_psX;
 Modelica.Media.Water.WaterIF97_base.setState_pTX
Modelica.Media.Water.WaterIF97_base.setState_pTX
| Type | Name | Default | Description | 
|---|---|---|---|
| FixedPhase | phase | 0 | 2 for two-phase, 1 for one-phase, 0 if not known | 
| AbsolutePressure | p | Pressure [Pa] | |
| Temperature | T | Temperature [K] | |
| MassFraction | X[:] | reference_X | Mass fractions [kg/kg] | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | thermodynamic state record | 
redeclare function extends setState_pTX 
  "Return thermodynamic state of water as function of p and T"
algorithm 
  state := ThermodynamicState(
    d=density_pT(p,T),
    T=T,
    phase=1,
    h=specificEnthalpy_pT(p,T),
    p=p);
end setState_pTX;
 Modelica.Media.Water.WaterIF97_base.setSmoothState
Modelica.Media.Water.WaterIF97_base.setSmoothState
| Type | Name | Default | Description | 
|---|---|---|---|
| Real | x | m_flow or dp | |
| ThermodynamicState | state_a | Thermodynamic state if x > 0 | |
| ThermodynamicState | state_b | Thermodynamic state if x < 0 | |
| Real | x_small | Smooth transition in the region -x_small < x < x_small | 
| Type | Name | Description | 
|---|---|---|
| ThermodynamicState | state | Smooth thermodynamic state for all x (continuous and differentiable) | 
redeclare function extends setSmoothState 
  "Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b"
  import Modelica.Media.Common.smoothStep;
algorithm 
  state :=ThermodynamicState(
    p=smoothStep(x, state_a.p, state_b.p, x_small),
    h=smoothStep(x, state_a.h, state_b.h, x_small),
    d=density_ph(smoothStep(x, state_a.p, state_b.p, x_small),
                 smoothStep(x, state_a.h, state_b.h, x_small)),
    T=temperature_ph(smoothStep(x, state_a.p, state_b.p, x_small),
                     smoothStep(x, state_a.h, state_b.h, x_small)),
    phase=0);
end setSmoothState;