Buildings.Media.Antifreeze.PropyleneGlycolWater

Package with model for propylene glycol - water with constant properties

Information

This medium package models propylene glycol - water mixtures.

The mass density, specific heat capacity, thermal conductivity and viscosity are assumed constant and evaluated at a set temperature and mass fraction of propylene glycol within the mixture. The dependence of the four properties are shown on the figure below.

Relative variation of specific heat capacity with temperature

The accuracy of the thermophysical properties is dependent on the temperature variations encountered during simulations. The figure below shows the relative error of the the four properties over a 10 °C range around the temperature used to evaluate the constant properties. The maximum errors are 0.8 % for mass density, 1.5 % for specific heat capacity, 3.2 % for thermal conductivity and 250 % for dynamic viscosity.

Relative variation of specific heat capacity with temperature

The figure below shows the relative error of the the four properties over a 20 °C range around the temperature used to evaluate the constant proepties. The maximum errors are 1.6 % for mass density, 3.0 % for specific heat capacity, 6.2 % for thermal conductivity and 950 % for dynamic viscosity.

Relative variation of specific heat capacity with temperature

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

Limitations

Density, specific heat capacity, thermal conductivity and viscosity are constant. The propylene glycol/water mixture is modeled as an incompressible liquid. There are no phase changes. The medium is limited to temperatures below 100 °C and mass fractions below 0.60. As is the case for Buildings.Media.Water, this medium package should not be used if the simulation relies on the dynamic viscosity.

Typical use and important parameters

The temperature and mass fraction must be specified for the evaluation of the constant thermophysical properties. A typical use of the package is (e.g. for a temperature of 20 °C and a mass fraction of 0.40):

Medium = Buildings.Media.Antifreeze.PropyleneGlycolWater(property_T=293.15, X_a=0.40)

Extends from Modelica.Media.Interfaces.PartialSimpleMedium (Medium model with linear dependency of u, h from temperature. All other quantities, especially density, are constant.).

Package Content

Name Description
property_T Temperature for evaluation of constant fluid properties
X_a Mass fraction of propylene glycol in water
Buildings.Media.Antifreeze.PropyleneGlycolWater.BaseProperties BaseProperties Base properties
X_a_min=0. Minimum allowed mass fraction of propylene glycol in water
X_a_max=0.6 Maximum allowed mass fraction of propylene glycol in water
simplePropyleneGlycolWaterConstants  
proCoe Coefficients for evaluation of thermo-physical properties
Buildings.Media.Antifreeze.PropyleneGlycolWater.density_TX_a density_TX_a Evaluate density of antifreeze-water mixture
Buildings.Media.Antifreeze.PropyleneGlycolWater.dynamicViscosity_TX_a dynamicViscosity_TX_a Evaluate dynamic viscosity of antifreeze-water mixture
Buildings.Media.Antifreeze.PropyleneGlycolWater.fusionTemperature_TX_a fusionTemperature_TX_a Evaluate temperature of fusion of antifreeze-water mixture
Buildings.Media.Antifreeze.PropyleneGlycolWater.polynomialProperty polynomialProperty Evaluates thermophysical property from 2-variable polynomial
Buildings.Media.Antifreeze.PropyleneGlycolWater.specificHeatCapacityCp_TX_a specificHeatCapacityCp_TX_a Evaluate specific heat capacity of antifreeze-water mixture
Buildings.Media.Antifreeze.PropyleneGlycolWater.thermalConductivity_TX_a thermalConductivity_TX_a Evaluate thermal conductivity of antifreeze-water mixture
Inherited
cp_const=specificHeatCapacityCp_TX_a(T=property_T, X_a=X_a) Constant specific heat capacity at constant pressure
cv_const=cp_const Constant specific heat capacity at constant volume
d_const=density_TX_a(T=property_T, X_a=X_a) Constant density
eta_const=dynamicViscosity_TX_a(T=property_T, X_a=X_a) Constant dynamic viscosity
lambda_const=thermalConductivity_TX_a(T=property_T, X_a=X_a) Constant thermal conductivity
a_const=1484 Constant velocity of sound
T_min=fusionTemperature_TX_a(T=property_T, X_a=X_a) Minimum temperature valid for medium model
T_max=Modelica.Units.Conversions.from_degC(100) Maximum temperature valid for medium model
T0=273.15 Zero enthalpy temperature
MM_const=(X_a/simplePropyleneGlycolWaterConstants[1].molarMass + (1 - X_a)/0.018015268)^(-1) Molar mass
fluidConstants=simplePropyleneGlycolWaterConstants Fluid constants
Modelica.Media.Interfaces.PartialSimpleMedium.ThermodynamicState ThermodynamicState Thermodynamic state
Modelica.Media.Interfaces.PartialSimpleMedium.setState_pTX setState_pTX Return thermodynamic state from p, T, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.setState_phX setState_phX Return thermodynamic state from p, h, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.setState_psX setState_psX Return thermodynamic state from p, s, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.setState_dTX setState_dTX Return thermodynamic state from d, T, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.setSmoothState setSmoothState Return thermodynamic state so that it smoothly approximates: if x > 0 then state_a else state_b
Modelica.Media.Interfaces.PartialSimpleMedium.dynamicViscosity dynamicViscosity Return dynamic viscosity
Modelica.Media.Interfaces.PartialSimpleMedium.thermalConductivity thermalConductivity Return thermal conductivity
Modelica.Media.Interfaces.PartialSimpleMedium.pressure pressure Return pressure
Modelica.Media.Interfaces.PartialSimpleMedium.temperature temperature Return temperature
Modelica.Media.Interfaces.PartialSimpleMedium.density density Return density
Modelica.Media.Interfaces.PartialSimpleMedium.specificEnthalpy specificEnthalpy Return specific enthalpy
Modelica.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCp specificHeatCapacityCp Return specific heat capacity at constant pressure
Modelica.Media.Interfaces.PartialSimpleMedium.specificHeatCapacityCv specificHeatCapacityCv Return specific heat capacity at constant volume
Modelica.Media.Interfaces.PartialSimpleMedium.isentropicExponent isentropicExponent Return isentropic exponent
Modelica.Media.Interfaces.PartialSimpleMedium.velocityOfSound velocityOfSound Return velocity of sound
Modelica.Media.Interfaces.PartialSimpleMedium.specificEnthalpy_pTX specificEnthalpy_pTX Return specific enthalpy from p, T, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.temperature_phX temperature_phX Return temperature from p, h, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.density_phX density_phX Return density from p, h, and X or Xi
Modelica.Media.Interfaces.PartialSimpleMedium.specificInternalEnergy specificInternalEnergy Return specific internal energy
Modelica.Media.Interfaces.PartialSimpleMedium.specificEntropy specificEntropy Return specific entropy
Modelica.Media.Interfaces.PartialSimpleMedium.specificGibbsEnergy specificGibbsEnergy Return specific Gibbs energy
Modelica.Media.Interfaces.PartialSimpleMedium.specificHelmholtzEnergy specificHelmholtzEnergy Return specific Helmholtz energy
Modelica.Media.Interfaces.PartialSimpleMedium.isentropicEnthalpy isentropicEnthalpy Return isentropic enthalpy
Modelica.Media.Interfaces.PartialSimpleMedium.isobaricExpansionCoefficient isobaricExpansionCoefficient Returns overall the isobaric expansion coefficient beta
Modelica.Media.Interfaces.PartialSimpleMedium.isothermalCompressibility isothermalCompressibility Returns overall the isothermal compressibility factor
Modelica.Media.Interfaces.PartialSimpleMedium.density_derp_T density_derp_T Returns the partial derivative of density with respect to pressure at constant temperature
Modelica.Media.Interfaces.PartialSimpleMedium.density_derT_p density_derT_p Returns the partial derivative of density with respect to temperature at constant pressure
Modelica.Media.Interfaces.PartialSimpleMedium.density_derX density_derX Returns the partial derivative of density with respect to mass fractions at constant pressure and temperature
Modelica.Media.Interfaces.PartialSimpleMedium.molarMass molarMass Return the molar mass of the medium
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.pT Enumeration type for independent variables
mediumName="PropyleneGlycolWater(X_a = " + String(X_a) + ", property_T = " + String(property_T) + ")" 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.Units.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.FluidConstants FluidConstants Critical, triple, molecular and other standard data of fluid
Modelica.Media.Interfaces.PartialMedium.prandtlNumber prandtlNumber Return the Prandtl number
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cp heatCapacity_cp Alias for deprecated name
Modelica.Media.Interfaces.PartialMedium.heatCapacity_cv heatCapacity_cv Alias for deprecated name
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.specificEntropy_pTX specificEntropy_pTX Return specific enthalpy from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_pTX density_pTX Return density from p, T, and X or Xi
Modelica.Media.Interfaces.PartialMedium.temperature_psX temperature_psX Return temperature from p,s, and X or Xi
Modelica.Media.Interfaces.PartialMedium.density_psX density_psX Return density from p, s, and X or Xi
Modelica.Media.Interfaces.PartialMedium.specificEnthalpy_psX specificEnthalpy_psX Return specific enthalpy from p, s, and X or Xi
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.Units.SI.Temperature property_T
    "Temperature for evaluation of constant fluid properties";
  constant Modelica.Units.SI.MassFraction X_a
    "Mass fraction of propylene glycol in water";
  constant Modelica.Units.SI.MassFraction X_a_min=0.
    "Minimum allowed mass fraction of propylene glycol in water";
  constant Modelica.Units.SI.MassFraction X_a_max=0.6
    "Maximum allowed mass fraction of propylene glycol in water";
  constant Modelica.Media.Interfaces.Types.Basic.FluidConstants[1]
    simplePropyleneGlycolWaterConstants(
    each chemicalFormula="C3H8O2",
    each structureFormula="CH3CH(OH)CH2OH",
    each casRegistryNumber="57-55-6",
    each iupacName="1,2-Propylene glycol",
    each molarMass=0.07609);
  constant Buildings.Media.Antifreeze.BaseClasses.PropertyCoefficients proCoe(
    X_a_ref=0.307031,
    T_ref=Modelica.Units.Conversions.from_degC(32.7083),
    nX_a=6,
    nT={4,4,4,3,2,1},
    nTot=18,
    a_d={1.018e3,-5.406e-1,-2.666e-3,1.347e-5,7.604e-1,-9.450e-3,5.541e-5,-1.343e-7,
        -2.498e-3,2.700e-5,-4.018e-7,3.376e-9,-1.550e-4,2.829e-6,-7.175e-9,-1.131e-6,
        -2.221e-8,2.342e-8},
    a_eta={6.837e-1,-3.045e-2,2.525e-4,-1.399e-6,3.328e-2,-3.984e-4,4.332e-6,-1.860e-8,
        5.453e-5,-8.600e-8,-1.593e-8,-4.465e-11,-3.900e-6,1.054e-7,-1.589e-9,-1.587e-8,
        4.475e-10,3.564e-9},
    a_Tf={-1.325e1,-3.820e-5,7.865e-7,-1.733e-9,-6.631e-1,6.774e-6,-6.242e-8,-7.819e-10,
        -1.094e-2,5.332e-8,-4.169e-9,3.288e-11,-2.283e-4,-1.131e-8,1.918e-10,-3.409e-6,
        8.035e-11,1.465e-8},
    a_cp={3.882e3,2.699e0,-1.659e-3,-1.032e-5,-1.304e1,5.070e-2,-4.752e-5,
        1.522e-6,-1.598e-1,9.534e-5,1.167e-5,-4.870e-8,3.539e-4,3.102e-5,-2.950e-7,
        5.000e-5,-7.135e-7,-4.959e-7},
    a_lambda={4.513e-1,7.955e-4,3.482e-8,-5.966e-9,-4.795e-3,-1.678e-5,8.941e-8,
        1.493e-10,2.076e-5,1.563e-7,-4.615e-9,9.897e-12,-9.083e-8,-2.518e-9,
        6.543e-11,-5.952e-10,-3.605e-11,2.104e-11})
    "Coefficients for evaluation of thermo-physical properties";

Buildings.Media.Antifreeze.PropyleneGlycolWater.BaseProperties

Base properties

Information

This base properties model is identical to Modelica.Media.Water.ConstantPropertyLiquidWater, except that the equation u = cv_const*(T - reference_T) has been replaced by u=h because cp_const=cv_const. Also, the model checks if the mass fraction of the mixture is within the allowed limits.

Parameters

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

Modelica definition

redeclare model BaseProperties "Base properties" Temperature T(stateSelect= if preferredMediumStates then StateSelect.prefer else StateSelect.default) "Temperature of medium"; InputAbsolutePressure p "Absolute pressure of medium"; InputMassFraction[nXi] Xi=fill(0, 0) "Structurally independent mass fractions"; InputSpecificEnthalpy h "Specific enthalpy of medium"; Modelica.Units.SI.SpecificInternalEnergy u "Specific internal energy of medium"; Modelica.Units.SI.Density d=d_const "Density of medium"; Modelica.Units.SI.MassFraction[nX] X={1} "Mass fractions (= (component mass)/total mass m_i/m)"; final Modelica.Units.SI.SpecificHeatCapacity R_s=0 "Gas constant (of mixture if applicable)"; final Modelica.Units.SI.MolarMass MM=MM_const "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"; final parameter Boolean standardOrderComponents=true "If true, and reducedX = true, the last element of X will be computed from the other ones"; Modelica.Units.NonSI.Temperature_degC T_degC= Modelica.Units.Conversions.to_degC(T) "Temperature of medium in [degC]"; Modelica.Units.NonSI.Pressure_bar p_bar=Modelica.Units.Conversions.to_bar(p) "Absolute pressure of medium in [bar]"; // Local connector definition, used for equation balancing check connector InputAbsolutePressure = input Modelica.Units.SI.AbsolutePressure "Pressure as input signal connector"; connector InputSpecificEnthalpy = input Modelica.Units.SI.SpecificEnthalpy "Specific enthalpy as input signal connector"; connector InputMassFraction = input Modelica.Units.SI.MassFraction "Mass fraction as input signal connector"; equation assert(T >= T_min, " In " + getInstanceName() + ": Temperature T exceeded its minimum allowed value of " + String(T_min-273.15) + " degC (" + String(T_min) + " Kelvin) as required from medium model \"" + mediumName + "\"."); assert(T <= T_max, " In " + getInstanceName() + ": Temperature T exceeded its maximum allowed value of " + String(T_max-273.15) + " degC (" + String(T_max) + " Kelvin) as required from medium model \"" + mediumName + "\"."); assert(X_a >= X_a_min, " In " + getInstanceName() + ": Mass fraction x_a exceeded its minimum allowed value of " + String(X_a_min) + " as required from medium model \"" + mediumName + "\"."); assert(X_a <= X_a_max, " In " + getInstanceName() + ": Mass fraction x_a exceeded its maximum allowed value of " + String(X_a_max) + " as required from medium model \"" + mediumName + "\"."); h = cp_const*(T-reference_T); u = h; state.T = T; state.p = p; end BaseProperties;

Buildings.Media.Antifreeze.PropyleneGlycolWater.density_TX_a Buildings.Media.Antifreeze.PropyleneGlycolWater.density_TX_a

Evaluate density of antifreeze-water mixture

Information

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature of antifreeze-water mixture [K]
MassFractionX_a Mass fraction of antifreeze [1]

Outputs

TypeNameDescription
DensitydDensity of antifreeze-water mixture [kg/m3]

Modelica definition

replaceable function density_TX_a "Evaluate density of antifreeze-water mixture" extends Modelica.Icons.Function; input Modelica.Units.SI.Temperature T "Temperature of antifreeze-water mixture"; input Modelica.Units.SI.MassFraction X_a "Mass fraction of antifreeze"; output Modelica.Units.SI.Density d "Density of antifreeze-water mixture"; algorithm d :=polynomialProperty( X_a, T, proCoe.a_d); end density_TX_a;

Buildings.Media.Antifreeze.PropyleneGlycolWater.dynamicViscosity_TX_a Buildings.Media.Antifreeze.PropyleneGlycolWater.dynamicViscosity_TX_a

Evaluate dynamic viscosity of antifreeze-water mixture

Information

Dynamic viscosity of antifreeze-water mixture at specified mass fraction and temperature, based on Melinder (2010).

References

Melinder, Åke. 2010. Properties of Secondary Working Fluids (Secondary Refrigerants or Coolants, Heat Transfer Fluids) for Indirect Systems. Paris: IIR/IIF.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature of antifreeze-water mixture [K]
MassFractionX_a Mass fraction of antifreeze [1]

Outputs

TypeNameDescription
DynamicViscosityetaDynamic Viscosity of antifreeze-water mixture [Pa.s]

Modelica definition

replaceable function dynamicViscosity_TX_a "Evaluate dynamic viscosity of antifreeze-water mixture" extends Modelica.Icons.Function; input Modelica.Units.SI.Temperature T "Temperature of antifreeze-water mixture"; input Modelica.Units.SI.MassFraction X_a "Mass fraction of antifreeze"; output Modelica.Units.SI.DynamicViscosity eta "Dynamic Viscosity of antifreeze-water mixture"; algorithm eta :=1e-3*exp(polynomialProperty( X_a, T, proCoe.a_eta)); end dynamicViscosity_TX_a;

Buildings.Media.Antifreeze.PropyleneGlycolWater.fusionTemperature_TX_a Buildings.Media.Antifreeze.PropyleneGlycolWater.fusionTemperature_TX_a

Evaluate temperature of fusion of antifreeze-water mixture

Information

Fusion temperature of antifreeze-water mixture at specified mass fraction and temperature, based on Melinder (2010).

References

Melinder, Åke. 2010. Properties of Secondary Working Fluids (Secondary Refrigerants or Coolants, Heat Transfer Fluids) for Indirect Systems. Paris: IIR/IIF.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature of antifreeze-water mixture [K]
MassFractionX_a Mass fraction of antifreeze [1]

Outputs

TypeNameDescription
TemperatureTfTemperature of fusion of antifreeze-water mixture [K]

Modelica definition

replaceable function fusionTemperature_TX_a "Evaluate temperature of fusion of antifreeze-water mixture" extends Modelica.Icons.Function; input Modelica.Units.SI.Temperature T "Temperature of antifreeze-water mixture"; input Modelica.Units.SI.MassFraction X_a "Mass fraction of antifreeze"; output Modelica.Units.SI.Temperature Tf "Temperature of fusion of antifreeze-water mixture"; algorithm Tf :=Modelica.Units.Conversions.from_degC(polynomialProperty( X_a, T, proCoe.a_Tf)); end fusionTemperature_TX_a;

Buildings.Media.Antifreeze.PropyleneGlycolWater.polynomialProperty Buildings.Media.Antifreeze.PropyleneGlycolWater.polynomialProperty

Evaluates thermophysical property from 2-variable polynomial

Information

Evaluates a thermophysical property of a mixture, based on correlations proposed by Melinder (2010).

The polynomial has the form

f = a1 (x-xm)0(y-ym)0 + a2 (x-xm)0(y-ym)1 + ... + any[1] (x-xm)0(y-ym)ny[1]-1 + ... + any[1])+1 (x-xm)1(y-ym)0 + ... + any[1]+ny[2] (x-xm)1(y-ym)ny[2]-1 + ...

References

Melinder, Åke. 2010. Properties of Secondary Working Fluids (Secondary Refrigerants or Coolants, Heat Transfer Fluids) for Indirect Systems. Paris: IIR/IIF.

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

Inputs

TypeNameDefaultDescription
Realx First independent variable
Realy Second independent variable
Reala[sum(proCoe.nT)] Polynomial coefficients

Outputs

TypeNameDescription
RealfValue of thermophysical property

Modelica definition

replaceable function polynomialProperty "Evaluates thermophysical property from 2-variable polynomial" extends Modelica.Icons.Function; input Real x "First independent variable"; input Real y "Second independent variable"; input Real a[sum(proCoe.nT)] "Polynomial coefficients"; output Real f "Value of thermophysical property"; protected Real dx; Real dy; Integer n; algorithm dx := 100*(x - proCoe.X_a_ref); dy := y - proCoe.T_ref; f := 0; n := 0; for i in 0:proCoe.nX_a - 1 loop for j in 0:proCoe.nT[i+1] - 1 loop n := n + 1; f := f + a[n]*dx^i*dy^j; end for; end for; end polynomialProperty;

Buildings.Media.Antifreeze.PropyleneGlycolWater.specificHeatCapacityCp_TX_a Buildings.Media.Antifreeze.PropyleneGlycolWater.specificHeatCapacityCp_TX_a

Evaluate specific heat capacity of antifreeze-water mixture

Information

Specific heat capacity of antifreeze-water mixture at specified mass fraction and temperature, based on Melinder (2010).

References

Melinder, Åke. 2010. Properties of Secondary Working Fluids (Secondary Refrigerants or Coolants, Heat Transfer Fluids) for Indirect Systems. Paris: IIR/IIF.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature of antifreeze-water mixture [K]
MassFractionX_a Mass fraction of antifreeze [1]

Outputs

TypeNameDescription
SpecificHeatCapacitycpSpecific heat capacity of antifreeze-water mixture [J/(kg.K)]

Modelica definition

replaceable function specificHeatCapacityCp_TX_a "Evaluate specific heat capacity of antifreeze-water mixture" extends Modelica.Icons.Function; input Modelica.Units.SI.Temperature T "Temperature of antifreeze-water mixture"; input Modelica.Units.SI.MassFraction X_a "Mass fraction of antifreeze"; output Modelica.Units.SI.SpecificHeatCapacity cp "Specific heat capacity of antifreeze-water mixture"; algorithm cp :=polynomialProperty( X_a, T, proCoe.a_cp); end specificHeatCapacityCp_TX_a;

Buildings.Media.Antifreeze.PropyleneGlycolWater.thermalConductivity_TX_a Buildings.Media.Antifreeze.PropyleneGlycolWater.thermalConductivity_TX_a

Evaluate thermal conductivity of antifreeze-water mixture

Information

Thermal conductivity of antifreeze-water mixture at specified mass fraction and temperature, based on Melinder (2010).

References

Melinder, Åke. 2010. Properties of Secondary Working Fluids (Secondary Refrigerants or Coolants, Heat Transfer Fluids) for Indirect Systems. Paris: IIR/IIF.

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

Inputs

TypeNameDefaultDescription
TemperatureT Temperature of antifreeze-water mixture [K]
MassFractionX_a Mass fraction of antifreeze [1]

Outputs

TypeNameDescription
ThermalConductivitylambdaThermal conductivity of antifreeze-water mixture [W/(m.K)]

Modelica definition

replaceable function thermalConductivity_TX_a "Evaluate thermal conductivity of antifreeze-water mixture" extends Modelica.Icons.Function; input Modelica.Units.SI.Temperature T "Temperature of antifreeze-water mixture"; input Modelica.Units.SI.MassFraction X_a "Mass fraction of antifreeze"; output Modelica.Units.SI.ThermalConductivity lambda "Thermal conductivity of antifreeze-water mixture"; algorithm lambda :=polynomialProperty( X_a, T, proCoe.a_lambda); end thermalConductivity_TX_a;

Buildings.Media.Antifreeze.PropyleneGlycolWater.BaseProperties.InputAbsolutePressure

Pressure as input signal connector

Modelica definition

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

Buildings.Media.Antifreeze.PropyleneGlycolWater.BaseProperties.InputMassFraction

Mass fraction as input signal connector

Modelica definition

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

Buildings.Media.Antifreeze.PropyleneGlycolWater.BaseProperties.InputSpecificEnthalpy

Specific enthalpy as input signal connector

Modelica definition

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