Buildings.Media.Specialized.Water.TemperatureDependentDensity

Package with model for liquid water with temperature-dependent density

Information

This medium package models liquid water.

The mass density is computed using a 3rd order polynomial, which yields the density as a function of temperature as shown in the figure below. Note, however, that computing density as a function of temperature can lead to considerably slower computing time compared to using Buildings.Media.Water in which the density is a constant. We therefore recommend to use Buildings.Media.Water for typical building energy simulations.

Mass density as a function of temperature

For the specific heat capacities at constant pressure and at constant volume, a constant value of 4184 J/(kg K), which corresponds to 20°C is used. The figure below shows the relative error of the specific heat capacity that is introduced by this simplification. Using a constant value for the specific heat capacity allows to compute temperature from enthalpy without having to solve an implicit equation, and therefore leads to faster simulation.

Relative variation of specific heat capacity with temperature

Thermal conductivity is calculated as a function of temperature as shown in the figure below. The correlation used to calculate the thermal conductivity is

λ(T) = λ(298.15 K) ⋅ (-1.48445+4.12292⋅(T/298.15)-1.63866⋅(T/298.15)2),

where λ(298.15 K) = 0.6065 W/(m ⋅ K) is the adopted standard value of the thermal conductivity of water at 298.15 K and 0.1 MPa.

Thermal conductivity as a function of temperature

Dynamic viscosity is calculated as the product of density and kinematic viscosity, both temperature dependent. However, the kinematic viscosity has its own temperature dependent correlation, implemented at Buildings.Media.Specialized.Water.TemperatureDependentDensity.kinematicViscosity. Results of the kinematic viscosity as a function of temperature are shown in the figure below.

Kinematic viscosity as a function of temperature

The enthalpy is computed using the convention that h=0 if T=0 °C.

Limitations

Phase changes are not modeled.

Extends from Modelica.Media.Interfaces.PartialPureSubstance (Base class for pure substances of one chemical substance), Modelica.Icons.Package (Icon for standard packages).

Package Content

Name Description
Buildings.Media.Specialized.Water.TemperatureDependentDensity.FluidConstants FluidConstants  
Buildings.Media.Specialized.Water.TemperatureDependentDensity.ThermodynamicState ThermodynamicState Thermodynamic state variables
cp_const=4184 Specific heat capacity at constant pressure
Buildings.Media.Specialized.Water.TemperatureDependentDensity.BaseProperties BaseProperties Base properties
Buildings.Media.Specialized.Water.TemperatureDependentDensity.density density Return the density
Buildings.Media.Specialized.Water.TemperatureDependentDensity.dynamicViscosity dynamicViscosity Return the dynamic viscosity
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEnthalpy specificEnthalpy Return the specific enthalpy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.enthalpyOfLiquid enthalpyOfLiquid Return the specific enthalpy of liquid
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificInternalEnergy specificInternalEnergy Return the specific enthalpy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEntropy specificEntropy Return the specific entropy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificGibbsEnergy specificGibbsEnergy Return the specific Gibbs energy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHelmholtzEnergy specificHelmholtzEnergy Return the specific Helmholtz energy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.isentropicEnthalpy isentropicEnthalpy Return the isentropic enthalpy
Buildings.Media.Specialized.Water.TemperatureDependentDensity.isobaricExpansionCoefficient isobaricExpansionCoefficient Return the isobaric expansion coefficient
Buildings.Media.Specialized.Water.TemperatureDependentDensity.isothermalCompressibility isothermalCompressibility Return the isothermal compressibility factor
Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derp_T density_derp_T Return the partial derivative of density with respect to pressure at constant temperature
Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derT_p density_derT_p Return the partial derivative of density with respect to temperature at constant pressure
Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derX density_derX Return the partial derivative of density with respect to mass fractions at constant pressure and temperature
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCp specificHeatCapacityCp Return the specific heat capacity at constant pressure
Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCv specificHeatCapacityCv Return the specific heat capacity at constant volume
Buildings.Media.Specialized.Water.TemperatureDependentDensity.thermalConductivity thermalConductivity Return the thermal conductivity
Buildings.Media.Specialized.Water.TemperatureDependentDensity.pressure pressure Return the pressure
Buildings.Media.Specialized.Water.TemperatureDependentDensity.temperature temperature Return the temperature
Buildings.Media.Specialized.Water.TemperatureDependentDensity.molarMass molarMass Return the molar mass
Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_dTX setState_dTX Return thermodynamic state from d, T, and X or Xi
Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_phX setState_phX Return the thermodynamic state as function of pressure p, specific enthalpy h and composition X or Xi
Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_pTX setState_pTX Return the thermodynamic state as function of p, T and composition X or Xi
Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_psX setState_psX Return the thermodynamic state as function of p, s and composition X or Xi
cv_const=cp_const Specific heat capacity at constant volume
a_const=1484 Constant velocity of sound
MM_const=0.018015268 Molar mass
Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_specificHeatCapacityCp der_specificHeatCapacityCp Return the derivative of the specific heat capacity at constant pressure
Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_enthalpyOfLiquid der_enthalpyOfLiquid Temperature derivative of enthalpy of liquid per unit mass of liquid
Buildings.Media.Specialized.Water.TemperatureDependentDensity.kinematicViscosity kinematicViscosity Return the kinematic viscosity
Inherited
Modelica.Media.Interfaces.PartialPureSubstance.setState_pT setState_pT Return thermodynamic state from p and T
Modelica.Media.Interfaces.PartialPureSubstance.setState_ph setState_ph Return thermodynamic state from p and h
Modelica.Media.Interfaces.PartialPureSubstance.setState_ps setState_ps Return thermodynamic state from p and s
Modelica.Media.Interfaces.PartialPureSubstance.setState_dT setState_dT Return thermodynamic state from d and T
Modelica.Media.Interfaces.PartialPureSubstance.density_ph density_ph Return density from p and h
Modelica.Media.Interfaces.PartialPureSubstance.temperature_ph temperature_ph Return temperature from p and h
Modelica.Media.Interfaces.PartialPureSubstance.pressure_dT pressure_dT Return pressure from d and T
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_dT specificEnthalpy_dT Return specific enthalpy from d and T
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_ps specificEnthalpy_ps Return specific enthalpy from p and s
Modelica.Media.Interfaces.PartialPureSubstance.temperature_ps temperature_ps Return temperature from p and s
Modelica.Media.Interfaces.PartialPureSubstance.density_ps density_ps Return density from p and s
Modelica.Media.Interfaces.PartialPureSubstance.specificEnthalpy_pT specificEnthalpy_pT Return specific enthalpy from p and T
Modelica.Media.Interfaces.PartialPureSubstance.density_pT density_pT Return density from p and T
ThermoStates=Modelica.Media.Interfaces.Choices.IndependentVariables.T Enumeration type for independent variables
mediumName="WaterDetailed" 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 = 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=true = true if medium contains the equation X = reference_X
reference_p=300000 Reference pressure of Medium: default 1 atmosphere
reference_T=273.15 Reference temperature of Medium: default 25 deg Celsius
reference_X={1} Default mass fractions of medium
p_default=300000 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)
C_default=fill(0, nC) Default value for trace substances 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
C_nominal=1.0e-6*ones(nC) Default for the nominal values for the extra properties
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.prandtlNumber prandtlNumber Return the Prandtl number
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp heatCapacity_cp Alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cv heatCapacity_cv Alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.isentropicExponent isentropicExponent Return isentropic exponent
Modelica.Media.Interfaces.PartialMedium.velocityOfSound velocityOfSound Return velocity of sound
Modelica.Media.Interfaces.PartialMedium.beta beta Alias for isobaricExpansionCoefficient for user convenience
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.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
Modelica.Media.Interfaces.PartialMedium.MassFlowRate MassFlowRate Type for mass flow rate with medium specific attributes
Modelica.Media.Interfaces.Types.AbsolutePressure AbsolutePressure Type for absolute pressure with medium specific attributes
Modelica.Media.Interfaces.Types.Density Density Type for density with medium specific attributes
Modelica.Media.Interfaces.Types.DynamicViscosity DynamicViscosity Type for dynamic viscosity with medium specific attributes
Modelica.Media.Interfaces.Types.EnthalpyFlowRate EnthalpyFlowRate Type for enthalpy flow rate with medium specific attributes
Modelica.Media.Interfaces.Types.MassFraction MassFraction Type for mass fraction with medium specific attributes
Modelica.Media.Interfaces.Types.MoleFraction MoleFraction Type for mole fraction with medium specific attributes
Modelica.Media.Interfaces.Types.MolarMass MolarMass Type for molar mass with medium specific attributes
Modelica.Media.Interfaces.Types.MolarVolume MolarVolume Type for molar volume with medium specific attributes
Modelica.Media.Interfaces.Types.IsentropicExponent IsentropicExponent Type for isentropic exponent with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEnergy SpecificEnergy Type for specific energy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificInternalEnergy SpecificInternalEnergy Type for specific internal energy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEnthalpy SpecificEnthalpy Type for specific enthalpy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificEntropy SpecificEntropy Type for specific entropy with medium specific attributes
Modelica.Media.Interfaces.Types.SpecificHeatCapacity SpecificHeatCapacity Type for specific heat capacity with medium specific attributes
Modelica.Media.Interfaces.Types.SurfaceTension SurfaceTension Type for surface tension with medium specific attributes
Modelica.Media.Interfaces.Types.Temperature Temperature Type for temperature with medium specific attributes
Modelica.Media.Interfaces.Types.ThermalConductivity ThermalConductivity Type for thermal conductivity with medium specific attributes
Modelica.Media.Interfaces.Types.PrandtlNumber PrandtlNumber Type for Prandtl number with medium specific attributes
Modelica.Media.Interfaces.Types.VelocityOfSound VelocityOfSound Type for velocity of sound with medium specific attributes
Modelica.Media.Interfaces.Types.ExtraProperty ExtraProperty Type for unspecified, mass-specific property transported by flow
Modelica.Media.Interfaces.Types.CumulativeExtraProperty CumulativeExtraProperty Type for conserved integral of unspecified, mass specific property
Modelica.Media.Interfaces.Types.ExtraPropertyFlowRate ExtraPropertyFlowRate Type for flow rate of unspecified, mass-specific property
Modelica.Media.Interfaces.Types.IsobaricExpansionCoefficient IsobaricExpansionCoefficient Type for isobaric expansion coefficient with medium specific attributes
Modelica.Media.Interfaces.Types.DipoleMoment DipoleMoment Type for dipole moment with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByPressure DerDensityByPressure Type for partial derivative of density with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByEnthalpy DerDensityByEnthalpy Type for partial derivative of density with respect to enthalpy with medium specific attributes
Modelica.Media.Interfaces.Types.DerEnthalpyByPressure DerEnthalpyByPressure Type for partial derivative of enthalpy with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.DerDensityByTemperature DerDensityByTemperature Type for partial derivative of density with respect to temperature with medium specific attributes
Modelica.Media.Interfaces.Types.DerTemperatureByPressure DerTemperatureByPressure Type for partial derivative of temperature with respect to pressure with medium specific attributes
Modelica.Media.Interfaces.Types.SaturationProperties SaturationProperties Saturation properties of two phase medium
Modelica.Media.Interfaces.Types.FluidLimits FluidLimits Validity limits for fluid model
Modelica.Media.Interfaces.Types.FixedPhase FixedPhase Phase of the fluid: 1 for 1-phase, 2 for two-phase, 0 for not known, e.g., interactive use
Modelica.Media.Interfaces.Types.Basic Basic The most basic version of a record used in several degrees of detail
Modelica.Media.Interfaces.Types.IdealGas IdealGas The ideal gas version of a record used in several degrees of detail
Modelica.Media.Interfaces.Types.TwoPhase TwoPhase The two phase fluid version of a record used in several degrees of detail

Types and constants

  constant Modelica.SIunits.SpecificHeatCapacity cp_const = 4184
    "Specific heat capacity at constant pressure";
  final constant Modelica.SIunits.SpecificHeatCapacity cv_const = cp_const
    "Specific heat capacity at constant volume";
  constant Modelica.SIunits.VelocityOfSound a_const=1484
    "Constant velocity of sound";
  constant Modelica.SIunits.MolarMass MM_const=0.018015268 "Molar mass";

Buildings.Media.Specialized.Water.TemperatureDependentDensity.FluidConstants Buildings.Media.Specialized.Water.TemperatureDependentDensity.FluidConstants

Parameters

TypeNameDefaultDescription
Custom Parameters
StringiupacName"oxidane"Complete IUPAC name (or common name, if non-existent)
StringcasRegistryNumber"7732-18-5"Chemical abstracts sequencing number (if it exists)
StringchemicalFormula"H2O"Chemical formula, (brutto, nomenclature according to Hill
StringstructureFormula"H2O"Chemical structure formula
MolarMassmolarMassMM_constMolar mass [kg/mol]

Modelica definition

redeclare record FluidConstants = Modelica.Media.Interfaces.Types.Basic.FluidConstants ( each chemicalFormula="H2O", each structureFormula="H2O", each casRegistryNumber="7732-18-5", each iupacName="oxidane", each molarMass=MM_const);

Buildings.Media.Specialized.Water.TemperatureDependentDensity.ThermodynamicState Buildings.Media.Specialized.Water.TemperatureDependentDensity.ThermodynamicState

Thermodynamic state variables

Information

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

Modelica definition

redeclare record extends ThermodynamicState "Thermodynamic state variables" Temperature T(start=T_default) "Temperature of medium"; AbsolutePressure p(start=p_default) "Pressure of medium"; end ThermodynamicState;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.BaseProperties Buildings.Media.Specialized.Water.TemperatureDependentDensity.BaseProperties

Base properties

Information

Base properties of the medium.

Extends from .

Parameters

TypeNameDefaultDescription
Advanced
BooleanpreferredMediumStatestrue= true if StateSelect.prefer shall be used for the independent property variables of the medium

Modelica definition

redeclare model extends BaseProperties( preferredMediumStates=true) "Base properties" equation h = (T - reference_T)*cp_const; u = h-reference_p/d; d = density(state); state.T = T; state.p = p; R=Modelica.Constants.R; MM=MM_const; end BaseProperties;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.density Buildings.Media.Specialized.Water.TemperatureDependentDensity.density

Return the density

Information

This function computes the density as a function of temperature.

Implementation

The function is based on the IDA implementation in therpro.nmf, which implements

d := 1000.12 + 1.43711e-2*T_degC -
 5.83576e-3*T_degC^2 + 1.5009e-5*T_degC^3;
 

This has been converted to Kelvin, which resulted in the above expression. In addition, below 5 °C and above 100 °C, the density is replaced by a linear function to avoid inflection points. This linear extension is such that the density is once continuously differentiable.

Extends from (Return density).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DensitydDensity [kg/m3]

Modelica definition

redeclare function extends density "Return the density" algorithm d := smooth(1, if state.T < 278.15 then -0.042860825*state.T + 1011.9695761 elseif state.T < 373.15 then 0.000015009*state.T^3 - 0.01813488505*state.T^2 + 6.5619527954075*state.T + 254.900074971947 else -0.7025109*state.T + 1220.35045233); end density;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.dynamicViscosity Buildings.Media.Specialized.Water.TemperatureDependentDensity.dynamicViscosity

Return the dynamic viscosity

Information

This function computes the dynamic viscosity.

Extends from (Return dynamic viscosity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
DynamicViscosityetaDynamic viscosity [Pa.s]

Modelica definition

redeclare function extends dynamicViscosity "Return the dynamic viscosity" algorithm eta := density(state)*kinematicViscosity(state.T); end dynamicViscosity;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEnthalpy Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEnthalpy

Return the specific enthalpy

Information

This function computes the specific enthalpy.

Extends from (Return specific enthalpy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

redeclare function extends specificEnthalpy "Return the specific enthalpy" algorithm h := (state.T - reference_T)*cp_const; end specificEnthalpy;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.enthalpyOfLiquid Buildings.Media.Specialized.Water.TemperatureDependentDensity.enthalpyOfLiquid

Return the specific enthalpy of liquid

Information

This function computes the specific enthalpy of liquid water.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]

Outputs

TypeNameDescription
SpecificEnthalpyhSpecific enthalpy [J/kg]

Modelica definition

function enthalpyOfLiquid "Return the specific enthalpy of liquid" annotation(derivative=der_enthalpyOfLiquid); extends Modelica.Icons.Function; input Modelica.SIunits.Temperature T "Temperature"; output Modelica.SIunits.SpecificEnthalpy h "Specific enthalpy"; algorithm h := (T - reference_T)*cp_const; end enthalpyOfLiquid;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificInternalEnergy Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificInternalEnergy

Return the specific enthalpy

Information

This function computes the specific internal energy.

Extends from (Return specific internal energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyuSpecific internal energy [J/kg]

Modelica definition

redeclare function extends specificInternalEnergy "Return the specific enthalpy" algorithm u := specificEnthalpy(state) - reference_p/density(state); end specificInternalEnergy;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEntropy Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEntropy

Return the specific entropy

Information

This function computes the specific entropy.

To obtain the state for a given pressure, entropy and mass fraction, use Buildings.Media.Air.setState_psX.

Extends from Modelica.Icons.Function (Icon for functions), (Return specific entropy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends specificEntropy "Return the specific entropy" extends Modelica.Icons.Function; algorithm s := cv_const*Modelica.Math.log(state.T/reference_T); end specificEntropy;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificGibbsEnergy Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificGibbsEnergy

Return the specific Gibbs energy

Information

This function computes the specific Gibbs energy.

Extends from Modelica.Icons.Function (Icon for functions), (Return specific Gibbs energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergygSpecific Gibbs energy [J/kg]

Modelica definition

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

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHelmholtzEnergy Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHelmholtzEnergy

Return the specific Helmholtz energy

Information

This function computes the specific Helmholtz energy.

Extends from Modelica.Icons.Function (Icon for functions), (Return specific Helmholtz energy).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
SpecificEnergyfSpecific Helmholtz energy [J/kg]

Modelica definition

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

Buildings.Media.Specialized.Water.TemperatureDependentDensity.isentropicEnthalpy Buildings.Media.Specialized.Water.TemperatureDependentDensity.isentropicEnthalpy

Return the isentropic enthalpy

Information

This function computes the specific enthalpy for an isentropic state change from the temperature that corresponds to the state refState to reference_T.

Extends from (Return isentropic enthalpy).

Inputs

TypeNameDefaultDescription
AbsolutePressurep_downstream Downstream pressure [Pa]
ThermodynamicStaterefState Reference state for entropy

Outputs

TypeNameDescription
SpecificEnthalpyh_isIsentropic enthalpy [J/kg]

Modelica definition

redeclare function extends isentropicEnthalpy "Return the isentropic enthalpy" algorithm h_is := specificEnthalpy(setState_psX( p=p_downstream, s=specificEntropy(refState), X={1})); end isentropicEnthalpy;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.isobaricExpansionCoefficient Buildings.Media.Specialized.Water.TemperatureDependentDensity.isobaricExpansionCoefficient

Return the isobaric expansion coefficient

Information

This function returns the isobaric expansion coefficient,

βp = - 1 ⁄ v   (∂ v ⁄ ∂ T)p,

where v is the specific volume, T is the temperature and p is the pressure.

Extends from (Return overall the isobaric expansion coefficient beta).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
IsobaricExpansionCoefficientbetaIsobaric expansion coefficient [1/K]

Modelica definition

redeclare function extends isobaricExpansionCoefficient "Return the isobaric expansion coefficient" algorithm beta := -smooth(0, if state.T < 278.15 then 0.042860825*(0.042860825*state.T - 1011.9695761)/(-0.042860825*state.T + 1011.9695761)^2 elseif state.T < 373.15 then (4.5027e-5*state.T^2 - 0.0362697701*state.T + 6.5619527954075)/ (1.5009e-5*state.T^3 - 0.01813488505*state.T^2 + 6.5619527954075*state.T + 254.900074971947) else 0.7025109*(0.7025109*state.T - 1220.35045233)/(-0.7025109*state.T + 1220.35045233)^2); // Symbolic conversion of degC to Kelvin // ((4.5027e-05)*T_degC^2 - 0.01167152*state.T + // 3.202446788)/((1.5009e-05)*T_degC^3 - 0.00583576*T_degC^2 + // 0.0143711*state.T + 996.194534035) end isobaricExpansionCoefficient;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.isothermalCompressibility Buildings.Media.Specialized.Water.TemperatureDependentDensity.isothermalCompressibility

Return the isothermal compressibility factor

Information

This function returns the isothermal compressibility coefficient, which is zero as this medium is incompressible. The isothermal compressibility is defined as

κT = - 1 ⁄ v   (∂ v ⁄ ∂ p)T,

where v is the specific volume, T is the temperature and p is the pressure.

Extends from (Return overall the isothermal compressibility factor).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
IsothermalCompressibilitykappaIsothermal compressibility [1/Pa]

Modelica definition

redeclare function extends isothermalCompressibility "Return the isothermal compressibility factor" algorithm kappa := 0; end isothermalCompressibility;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derp_T Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derp_T

Return the partial derivative of density with respect to pressure at constant temperature

Information

This function returns the partial derivative of density with respect to pressure at constant temperature, which is zero as the medium is incompressible.

Extends from (Return density derivative w.r.t. pressure at const temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends density_derp_T "Return the partial derivative of density with respect to pressure at constant temperature" algorithm ddpT := 0; end density_derp_T;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derT_p Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derT_p

Return the partial derivative of density with respect to temperature at constant pressure

Information

This function computes the derivative of density with respect to temperature at constant pressure.

Extends from (Return density derivative w.r.t. temperature at constant pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends density_derT_p "Return the partial derivative of density with respect to temperature at constant pressure" algorithm ddTp := if state.T < 278.15 then -0.042860825 elseif state.T < 373.15 then (0.0000450270000000000*state.T^2 - 0.0362697701000000*state.T + 6.56195279540750) else -0.7025109; end density_derT_p;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derX Buildings.Media.Specialized.Water.TemperatureDependentDensity.density_derX

Return the partial derivative of density with respect to mass fractions at constant pressure and temperature

Information

This function returns the partial derivative of density with respect to mass fraction, which is zero as the medium is a single substance.

Extends from (Return density derivative w.r.t. mass fraction).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends density_derX "Return the partial derivative of density with respect to mass fractions at constant pressure and temperature" algorithm dddX := fill(0, nX); end density_derX;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCp Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCp

Return the specific heat capacity at constant pressure

Information

This function returns the specific heat capacity at constant pressure.

Extends from (Return specific heat capacity at constant pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare replaceable function extends specificHeatCapacityCp "Return the specific heat capacity at constant pressure" annotation(derivative=der_specificHeatCapacityCp); algorithm cp := cp_const; end specificHeatCapacityCp;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCv Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificHeatCapacityCv

Return the specific heat capacity at constant volume

Information

This function computes the specific heat capacity at constant volume.

Extends from (Return specific heat capacity at constant volume).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare replaceable function extends specificHeatCapacityCv "Return the specific heat capacity at constant volume" annotation(derivative=der_specificHeatCapacityCp); algorithm cv := cv_const; end specificHeatCapacityCv;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.thermalConductivity Buildings.Media.Specialized.Water.TemperatureDependentDensity.thermalConductivity

Return the thermal conductivity

Information

This function returns the thermal conductivity. The expression is obtained from Ramires et al. (1995).

References

Ramires, Maria L. V. and Nieto de Castro, Carlos A. and Nagasaka, Yuchi and Nagashima, Akira and Assael, Marc J. and Wakeham, William A. Standard Reference Data for the Thermal Conductivity of Water. Journal of Physical and Chemical Reference Data, 24, p. 1377-1381, 1995. DOI:10.1063/1.555963.

Extends from (Return thermal conductivity).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

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

Modelica definition

redeclare function extends thermalConductivity "Return the thermal conductivity" algorithm lambda :=0.6065*(-1.48445 + 4.12292*(state.T/298.15) - 1.63866*(state.T/298.15)^2); end thermalConductivity;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.pressure Buildings.Media.Specialized.Water.TemperatureDependentDensity.pressure

Return the pressure

Information

This function returns the pressure.

Extends from (Return pressure).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
AbsolutePressurepPressure [Pa]

Modelica definition

redeclare function extends pressure "Return the pressure" algorithm p := state.p; end pressure;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.temperature Buildings.Media.Specialized.Water.TemperatureDependentDensity.temperature

Return the temperature

Information

This function returns the temperature.

Extends from (Return temperature).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
TemperatureTTemperature [K]

Modelica definition

redeclare function extends temperature "Return the temperature" algorithm T := state.T; end temperature;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.molarMass Buildings.Media.Specialized.Water.TemperatureDependentDensity.molarMass

Return the molar mass

Information

This function returns the molar mass, which is assumed to be constant.

Extends from (Return the molar mass of the medium).

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state record

Outputs

TypeNameDescription
MolarMassMMMixture molar mass [kg/mol]

Modelica definition

redeclare function extends molarMass "Return the molar mass" algorithm MM := MM_const; end molarMass;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_dTX Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_dTX

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

Information

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

Inputs

TypeNameDefaultDescription
Densityd Density [kg/m3]
TemperatureT Temperature [K]
MassFractionX[:]reference_XMass fractions [kg/kg]

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function setState_dTX "Return thermodynamic state from d, T, and X or Xi" extends Modelica.Icons.Function; input Density d "Density"; input Temperature T "Temperature"; input MassFraction X[:]=reference_X "Mass fractions"; output ThermodynamicState state "Thermodynamic state record"; algorithm assert(false, "Pressure can not be computed from temperature and density for an incompressible fluid!"); end setState_dTX;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_phX Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_phX

Return the thermodynamic state as function of pressure p, specific enthalpy h and composition X or Xi

Information

This function returns the thermodynamic state for a given pressure, specific enthalpy and composition.

Extends from (Return thermodynamic state as function of p, h and composition X or Xi).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_phX "Return the thermodynamic state as function of pressure p, specific enthalpy h and composition X or Xi" algorithm state := ThermodynamicState(p=p, T=reference_T + h/cp_const); end setState_phX;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_pTX Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_pTX

Return the thermodynamic state as function of p, T and composition X or Xi

Information

This function returns the thermodynamic state for a given pressure, temperature and composition.

Extends from (Return thermodynamic state as function of p, T and composition X or Xi).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_pTX "Return the thermodynamic state as function of p, T and composition X or Xi" algorithm state := ThermodynamicState(p=p, T=T); end setState_pTX;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_psX Buildings.Media.Specialized.Water.TemperatureDependentDensity.setState_psX

Return the thermodynamic state as function of p, s and composition X or Xi

Information

This function returns the thermodynamic state based on pressure, specific entropy and mass fraction.

The state is computed by symbolically solving Buildings.Media.Specialized.Water.TemperatureDependentDensity.specificEntropy for temperature.

Extends from (Return thermodynamic state as function of p, s and composition X or Xi).

Inputs

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

Outputs

TypeNameDescription
ThermodynamicStatestateThermodynamic state record

Modelica definition

redeclare function extends setState_psX "Return the thermodynamic state as function of p, s and composition X or Xi" algorithm // The temperature is obtained from symbolic solving the // specificEntropy function for T, i.e., // s := cv_const*Modelica.Math.log(state.T/reference_T) state := ThermodynamicState(p=p, T=reference_T * Modelica.Math.exp(s/cv_const)); end setState_psX;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_specificHeatCapacityCp Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_specificHeatCapacityCp

Return the derivative of the specific heat capacity at constant pressure

Information

This function computes the derivative of the specific heat capacity at constant pressure with respect to the state.

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

Inputs

TypeNameDefaultDescription
ThermodynamicStatestate Thermodynamic state
ThermodynamicStateder_state Derivative of thermodynamic state

Outputs

TypeNameDescription
Realder_cpDerivative of specific heat capacity [J/(kg.K.s)]

Modelica definition

replaceable function der_specificHeatCapacityCp "Return the derivative of the specific heat capacity at constant pressure" extends Modelica.Icons.Function; input ThermodynamicState state "Thermodynamic state"; input ThermodynamicState der_state "Derivative of thermodynamic state"; output Real der_cp(unit="J/(kg.K.s)") "Derivative of specific heat capacity"; algorithm der_cp := 0; end der_specificHeatCapacityCp;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_enthalpyOfLiquid Buildings.Media.Specialized.Water.TemperatureDependentDensity.der_enthalpyOfLiquid

Temperature derivative of enthalpy of liquid per unit mass of liquid

Information

This function computes the temperature derivative of the enthalpy of liquid water per unit mass.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]
Realder_T Temperature derivative

Outputs

TypeNameDescription
Realder_hDerivative of liquid enthalpy

Modelica definition

replaceable function der_enthalpyOfLiquid "Temperature derivative of enthalpy of liquid per unit mass of liquid" extends Modelica.Icons.Function; input Temperature T "Temperature"; input Real der_T "Temperature derivative"; output Real der_h "Derivative of liquid enthalpy"; algorithm der_h := cp_const*der_T; end der_enthalpyOfLiquid;

Buildings.Media.Specialized.Water.TemperatureDependentDensity.kinematicViscosity Buildings.Media.Specialized.Water.TemperatureDependentDensity.kinematicViscosity

Return the kinematic viscosity

Information

This function computes the kinematic viscosity as a function of temperature.

Implementation

The function is based on the IDA implementation in therpro.nmf. The original equation is

kinVis :=1E-6*Modelica.Math.exp(0.577449 - 3.253945e-2*T_degC + 2.17369e-4*
      T_degC^2 - 7.22111e-7*T_degC^3);
      

This has been converted to Kelvin, which resulted in the above expression. In addition, at 5 °C the kinematic viscosity is linearly extrapolated to avoid a large gradient at very low temperatures. We selected the same point for the linearization as we used for the density, as the density and the kinematic viscosity are combined in Buildings.Media.Specialized.Water.TemperatureDependentDensity.dynamicViscosity.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature [K]

Outputs

TypeNameDescription
KinematicViscositykinVisKinematic viscosity [m2/s]

Modelica definition

function kinematicViscosity "Return the kinematic viscosity" extends Modelica.Icons.Function; input Modelica.SIunits.Temperature T "Temperature"; output Modelica.SIunits.KinematicViscosity kinVis "Kinematic viscosity"; algorithm kinVis := smooth(1, if T < 278.15 then -(4.63023776563e-08)*T + 1.44011135763e-05 else 1.0e-6*Modelica.Math.exp( -(7.22111000000000e-7)*T^3 + 0.000809102858950000*T^2 - 0.312920238272193*T + 40.4003044106506)); end kinematicViscosity;