FCSys.Characteristics.BaseClasses.Characteristic

Package of thermodynamic and diffusive properties

Information

This package is compatible with NASA CEA thermodynamic data [McBride2002] and the virial equation of state [Dymond2002].

Notes regarding the constants:

  • Currently, formula may not contain parentheses or brackets.
  • d is the Van der Waals diameter or the diameter for the rigid-sphere ("billiard-ball") approximation of the kinetic theory of gases [Present1958].
  • bc: The rows give the coefficients for the temperature intervals bounded by the values in Tlim c. The powers of T increase by column. By default, the powers of T for the first column are each -2, which corresponds to [McBride2002]. In that case, the dimensionalities of the coefficients are {L4.M2/(N2.T4), L2.M/(N.T2), 1, …} for each row, where L is length, M is mass, N is particle number, and T is time. (In FCSys, temperature is a potential with dimension L2.M/(N.T2); see the Units package.)
  • Bc: As in bc, the rows correspond to different temperature intervals. The first column is for specific enthalpy and has dimensionality L2.M/(N.T2). The second is for specific entropy and is dimensionless. The integration constants for enthalpy are defined such that the enthalpy at 25 °C is the specific enthalpy of formation at that temperature and reference pressure [McBride2002, p. 2]. The integration constants for specific entropy are defined such that specific entropy is absolute.
  • Tlim c: The first and last entries are the minimum and maximum valid temperatures. The intermediate entries are the thresholds between rows of bc (and Bc). Therefore, if there are n temperature intervals (and rows in bc and Bc), then Tlim c must have n + 1 entries.
  • The reference pressure is po. In the NASA CEA data [McBride2002], it is 1 bar for gases and 1 atm for condensed species. For gases, the reference state is the ideal gas at po. For example, the enthalpy of a non-ideal (real) gas at 25 °C and po with ReferenceEnthalpy.zeroAt25degC is not exactly zero.
  • If the material is gaseous (phase == Phase.gas), then the first virial coefficient must be independent of temperature. Otherwise, the function for specific enthalpy (h) will be ill-posed. Typically, the first virial coefficient is one (or equivalently U.R), which satisfies this requirement.
Extends from CharacteristicEOS (Base thermodynamic package with only the p-v-T relations).

Package Content

NameDescription
formulaChemical formula
phaseMaterial phase
mSpecific mass
dSpecific diameter
z=charge(formula)Charge number
referenceEnthalpy=ReferenceEnthalpy.enthalpyOfFormationAt25degCChoice of enthalpy reference
Deltah0_fEnthalpy of formation at 298.15 K, pohof)
Deltah0ho(298.15 K) - ho(0 K) (Δho)
h_offset=0Additional enthalpy offset (hoffset)
n_c=-2Power of T for 1st column of bc (nc)
T_lim_c={0,Modelica.Constants.inf}Temperature limits for the rows of bc and Bc (Tlim c)
b_cCoefficients of isobaric specific heat capacity at po as a polynomial in T (bc)
B_cIntegration constants for specific enthalpy and entropy (Bc)
alpha=1.5*U.pi^1.5*d^2*U.qScaled specific intercept area
FCSys.Characteristics.BaseClasses.Characteristic.omega omega Root mean square of thermal velocity in one dimension as a function of temperature (ω = √ T/m )
FCSys.Characteristics.BaseClasses.Characteristic.c_p c_p Isobaric specific heat capacity (cp) as a function of temperature and pressure
FCSys.Characteristics.BaseClasses.Characteristic.c_v c_v Isochoric specific heat capacity (cv) as a function of temperature and pressure
FCSys.Characteristics.BaseClasses.Characteristic.D D Diffusivity as a function of temperature and specific volume
FCSys.Characteristics.BaseClasses.Characteristic.g g Gibbs potential as a function of temperature and pressure
FCSys.Characteristics.BaseClasses.Characteristic.h h Specific enthalpy as a function of temperature and pressure
FCSys.Characteristics.BaseClasses.Characteristic.s s Specific entropy as a function of temperature and pressure
FCSys.Characteristics.BaseClasses.Characteristic.zeta zeta ** (ζ) as a function of temperature
FCSys.Characteristics.BaseClasses.Characteristic.eta eta Fluidity (η) as a function of temperature
FCSys.Characteristics.BaseClasses.Characteristic.theta theta Thermal resistivity (θ) as a function of temperature and specific volume
FCSys.Characteristics.BaseClasses.Characteristic.tauprime tauprime Phase change interval (τ′) as a function of temperature and specific volume
FCSys.Characteristics.BaseClasses.Characteristic.mu mu Mobility (μ) as a function of temperature and specific volume
FCSys.Characteristics.BaseClasses.Characteristic.nu nu Thermal independity (ν) as a function of temperature and specific volume
Inherited
p0=U.barReference pressure (po)
n_v={-1,0}Powers of p/T and T for 1st row and column of bv (nv)
b_v=[1]Coefficients for specific volume as a polynomial in p/T and T (bv)
isCompressible=anyTrue({anyTrue({abs(b_v[i, j]) > Modelica.Constants.small and n_v[1] + i - 1 <> 0 for i in 1:size(b_v, 1)}) for j in 1:size(b_v, 2)})true, if density depends on pressure
hasThermalExpansion=anyTrue({anyTrue({abs(b_v[i, j]) > Modelica.Constants.small and n_v[2] + j - n_v[1] - i <> 0 for i in 1:size(b_v, 1)}) for j in 1:size(b_v, 2)})true, if density depends on temperature
n_p={n_v[1] - size(b_v, 1) + 1,n_v[2] + 1}Powers of v and T for 1st row and column of bp
b_p=if size(b_v, 1) == 1 then b_v .^ (-n_p[1]) else {(if n_v[1] + i == 0 or n_v[1] + i == 1 or size(b_v, 1) == 1 then b_v[i, :] else (if n_v[1] + i == 2 and n_v[1] <= 0 then b_v[i, :] + b_v[i - 1, :] .^ 2 else (if n_v[1] + i == 3 and n_v[1] <= 0 then b_v[i, :] + b_v[i - 2, :] .* (b_v[i - 2, :] .^ 2 + 3*b_v[i - 1, :]) else zeros(size(b_v, 2))))) for i in size(b_v, 1):-1:1}Coefficients of p as a polynomial in v and T
FCSys.Characteristics.BaseClasses.CharacteristicEOS.dp_Tv dp_Tv Derivative of pressure as defined by pT v()
FCSys.Characteristics.BaseClasses.CharacteristicEOS.dv_Tp dv_Tp Derivative of specific volume as defined by vT p()
FCSys.Characteristics.BaseClasses.CharacteristicEOS.p_Tv p_Tv Pressure as a function of temperature and specific volume (pT v())
FCSys.Characteristics.BaseClasses.CharacteristicEOS.v_Tp v_Tp Specific volume as a function of temperature and pressure (vT p())
FCSys.Characteristics.BaseClasses.CharacteristicEOS.beta beta Isothermal compressibility as a function of temperature and pressure (β)

Types and constants

constant String formula "Chemical formula";
constant Phase phase "Material phase";
constant Q.MassSpecific m "Specific mass";
constant Q.LengthSpecific d "Specific diameter";
constant Integer z=charge(formula) "Charge number";
constant ReferenceEnthalpy referenceEnthalpy=ReferenceEnthalpy.enthalpyOfFormationAt25degC
  "Choice of enthalpy reference";
constant Q.PotentialChemical Deltah0_f 
  "Enthalpy of formation at 298.15 K, pohof)";
constant Q.PotentialChemical Deltah0 
  "ho(298.15 K) - ho(0 K) (Δho)";
constant Q.PotentialChemical h_offset=0 
  "Additional enthalpy offset (hoffset)";
constant Integer n_c=-2 
  "Power of T for 1st column of bc (nc)";
constant Q.TemperatureAbsolute T_lim_c[:]={0,Modelica.Constants.inf} 
  "Temperature limits for the rows of bc and Bc (Tlim c)";
constant Real b_c[size(T_lim_c, 1) - 1, :] 
  "Coefficients of isobaric specific heat capacity at po as a polynomial in T (bc)";
constant Real B_c[size(T_lim_c, 1) - 1, 2] 
  "Integration constants for specific enthalpy and entropy (Bc)";
constant Q.AreaSpecific alpha=1.5*U.pi^1.5*d^2*U.q 
  "Scaled specific intercept area";

FCSys.Characteristics.BaseClasses.Characteristic.omega FCSys.Characteristics.BaseClasses.Characteristic.omega

Root mean square of thermal velocity in one dimension as a function of temperature (ω = √ T/m )

Information

This function outputs the root mean square of the thermal velocity in any one dimension, assuming the speeds of the particles follow the Maxwell-Boltzmann distribution. In combination with α, this refactorization is convenient for calculating τ′, μ, ν, ζ, η, and θ.

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]

Outputs

TypeNameDescription
VelocityomegaRoot mean square of thermal velocity in one dimension [L/T]

Modelica definition

function omega 
  "Root mean square of thermal velocity in one dimension as a function of temperature (ω = √ T/m )"
  extends Modelica.Icons.Function;

  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  output Q.Velocity omega 
    "Root mean square of thermal velocity in one dimension";

algorithm 
  omega := sqrt(T/m);
end omega;

FCSys.Characteristics.BaseClasses.Characteristic.c_p FCSys.Characteristics.BaseClasses.Characteristic.c_p

Isobaric specific heat capacity (cp) as a function of temperature and pressure

Information

For an ideal gas, this function is independent of pressure (although pressure remains as a valid input).

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
PressureAbsolutepp0Pressure [M/(L.T2)]

Outputs

TypeNameDescription
CapacityThermalSpecificc_pIsobaric specific heat capacity [1]

Modelica definition

function c_p 
  "Isobaric specific heat capacity (cp) as a function of temperature and pressure"
  import FCSys.Utilities.Polynomial;
  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.PressureAbsolute p=p0 "Pressure";
  output Q.CapacityThermalSpecific c_p "Isobaric specific heat capacity";

protected 

    "Isobaric specific heat capacity at reference pressure as a function of temperature"
    import FCSys.Utilities.Polynomial;
    input Q.TemperatureAbsolute T "Temperature";
    output Q.CapacityThermalSpecific c0_p 
      "Isobaric specific heat capacity at reference pressure";

  algorithm 
    /*
  assert(T_lim_c[1] <= T and T <= T_lim_c[end], "Temperature " + String(T/U.K) + " K is out of bounds for "
     + formula + " ([" + String(T_lim_c[1]/U.K) + ", " + String(T_lim_c[size(
    T_lim_c, 1)]/U.K) + "] K).");
  */
    // Note:  This is commented out so that the function can be inlined.
    c0_p := smooth(0, sum(if (T_lim_c[i] <= T or i == 1) and (T < T_lim_c[i + 1]
       or i == size(T_lim_c, 1) - 1) then Polynomial.f(
        T,
        b_c[i, :],
        n_c) else 0 for i in 1:size(T_lim_c, 1) - 1));
  end c0_p;

"Residual isobaric specific heat capacity for pressure adjustment"
    import FCSys.Utilities.Polynomial;
    input Q.TemperatureAbsolute T "Temperature";
    input Q.PressureAbsolute p "Pressure";
    input Integer rowLimits[2]={1,size(b_v, 1)} 
      "Beginning and ending indices of rows of b_v to be included";
    output Q.CapacityThermalSpecific c_p_resid 
      "Temperature times the partial derivative of the integral of (dels/delp)_T*dp up to p w.r.t. T";

  algorithm 
    c_p_resid := Polynomial.F(
        p,
        {Polynomial.f(
          T,
          b_v[i, :] .* {(n_v[2] - n_v[1] + j - i)*(n_v[1] - n_v[2] + i - j + 1)
          for j in 1:size(b_v, 2)},
          n_v[2] - n_v[1] - i) for i in rowLimits[1]:rowLimits[2]},
        n_v[1]);
    // See s_resid() in Characteristic.s for the integral of (dels/delp)_T*dp.
    // This is temperature times the isobaric partial derivative of that
    // function with respect to temperature.  It is zero for an ideal gas.

  end c_p_resid;

algorithm 
  c_p := c0_p(T) + c_p_resid(T, p) - (if phase <> Phase.gas then c_p_resid(T,
    p0) else c_p_resid(
    T,
    p0,
    {1,-n_v[1]}));
  // See the notes in the algorithm of Characteristic.s.
  // Note:  [Dymond2002, p.17, eqs. 1.45 & 1.46] may be incorrect.
end c_p;

FCSys.Characteristics.BaseClasses.Characteristic.c_v FCSys.Characteristics.BaseClasses.Characteristic.c_v

Isochoric specific heat capacity (cv) as a function of temperature and pressure

Information

For an ideal gas, this function is independent of pressure (although pressure remains as a valid input).

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
PressureAbsolutepp0Pressure [M/(L.T2)]

Outputs

TypeNameDescription
CapacityThermalSpecificc_vIsochoric specific heat capacity [1]

Modelica definition

function c_v 
  "Isochoric specific heat capacity (cv) as a function of temperature and pressure"
  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.PressureAbsolute p=p0 "Pressure";
  output Q.CapacityThermalSpecific c_v "Isochoric specific heat capacity";

algorithm 
  c_v := c_p(T, p) - T*dp_Tv(
    T,
    v_Tp(T, p),
    dT=1,
    dv=0)*dv_Tp(
    T,
    p,
    dT=1,
    dp=0) "[Moran2004, p.&nbsp;546, eq. 11.66]";
  // Note 1:  This reduces to c_v = c_p - 1 for an ideal gas (where in
  // FCSys 1 = U.R).
  // Note 2:  [Dymond2002, p.17, eqs. 1.43 & 1.44] may be incorrect.
end c_v;

FCSys.Characteristics.BaseClasses.Characteristic.D FCSys.Characteristics.BaseClasses.Characteristic.D

Diffusivity as a function of temperature and specific volume

Information

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.
Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
DiffusivityDDiffusivity [L2.M/(N.T)]

Modelica definition

replaceable function D 
  "Diffusivity as a function of temperature and specific volume"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.Diffusivity D "Diffusivity";

algorithm 
  D := v*omega(T)/alpha;
end D;

FCSys.Characteristics.BaseClasses.Characteristic.g FCSys.Characteristics.BaseClasses.Characteristic.g

Gibbs potential as a function of temperature and pressure

Information

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
PressureAbsolutepp0Pressure [M/(L.T2)]

Outputs

TypeNameDescription
PotentialgGibbs potential [L2.M/(N.T2)]

Modelica definition

function g 
  "Gibbs potential as a function of temperature and pressure"
  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.PressureAbsolute p=p0 "Pressure";
  output Q.Potential g "Gibbs potential";

algorithm 
  g := h(T, p) - T*s(T, p);
end g;

FCSys.Characteristics.BaseClasses.Characteristic.h FCSys.Characteristics.BaseClasses.Characteristic.h

Specific enthalpy as a function of temperature and pressure

Information

For an ideal gas, this function is independent of pressure (although pressure remains as a valid input).

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
PressureAbsolutepp0Pressure [M/(L.T2)]

Outputs

TypeNameDescription
PotentialhSpecific enthalpy [L2.M/(N.T2)]

Modelica definition

function h 
  "Specific enthalpy as a function of temperature and pressure"
  import FCSys.Utilities.Polynomial;
  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.PressureAbsolute p=p0 "Pressure";
  output Q.Potential h "Specific enthalpy";

protected 

    "Return h0 as a function of T using one of the temperature intervals"
  annotation(derivative=dh0_i);
    import FCSys.Utilities.Polynomial;
    input Q.TemperatureAbsolute T "Temperature";
    input Integer i "Index of the temperature interval";
    output Q.Potential h0 
      "Specific enthalpy at given temperature relative to enthalpy of formation at 25 degC, both at reference pressure";

  algorithm 
    h0 := Polynomial.F(
        T,
        b_c[i, :],
        n_c) + B_c[i, 1];
    // This is the integral of c0_p*dT up to T at p0.  The lower bound is the
    // enthalpy of formation (of ideal gas, if the material is gaseous) at
    // 25 degC [McBride2002, p. 2].

  end h0_i;

  "Derivative of h0_i"
    // Note:  This function is necessary for Dymola 7.4 to differentiate h().
    import FCSys.Utilities.Polynomial;
    input Q.TemperatureAbsolute T "Temperature";
    input Integer i "Index of the temperature interval";
    input Q.Temperature dT "Derivative of temperature";
    output Q.Potential dh0 
      "Derivative of specific enthalpy at reference pressure";

  algorithm 
    dh0 := Polynomial.f(
        T,
        b_c[i, :],
        n_c)*dT;

  end dh0_i;

"Residual specific enthalpy for pressure adjustment for selected rows of b_v"
    import FCSys.Utilities.Polynomial;
    input Q.TemperatureAbsolute T "Temperature";
    input Q.PressureAbsolute p "Pressure";
    input Integer rowLimits[2]={1,size(b_v, 1)} 
      "Beginning and ending indices of rows of b_v to be included";
    output Q.Potential h_resid 
      "Integral of (delh/delp)_T*dp up to p with zero integration constant (for selected rows)";

  algorithm 
    h_resid := Polynomial.F(
        p,
        {Polynomial.f(
          T,
          b_v[i, :] .* {n_v[1] - n_v[2] + i - j + 1 for j in 1:size(b_v, 2)},
          n_v[2] - n_v[1] - i + 1) for i in rowLimits[1]:rowLimits[2]},
        n_v[1] + rowLimits[1] - 1);
    // Note:  The partial derivative (delh/delp)_T is equal to v +
    // T*(dels/delp)_T by definition of enthalpy change (dh = T*ds + v*dp)
    // and then to v - T*(delv/delT)_p by applying the appropriate Maxwell
    // relation, (dels/delp)_T = -(delv/delT)_p.
    // Note:  This is zero for an ideal gas.

  end h_resid;

algorithm 
  /*
    assert(T_lim_c[1] <= T and T <= T_lim_c[end], "Temperature " +
    String(T/U.K) + " K is out of bounds for " + name + " ([" + String(T_lim_c[1]
    /U.K) + ", " + String(T_lim_c[size(T_lim_c, 1)]/U.K) + "] K).");
  */
  // Note:  This is commented out so that the function can be inlined.
  h := smooth(1, sum(if (T_lim_c[i] <= T or i == 1) and (T < T_lim_c[i + 1] or 
    i == size(b_c, 1)) then h0_i(T, i) else 0 for i in 1:size(b_c, 1))) + (if 
    referenceEnthalpy == ReferenceEnthalpy.zeroAt0K then Deltah0 else 0) - (if 
    referenceEnthalpy == ReferenceEnthalpy.enthalpyOfFormationAt25degC then 0
     else Deltah0_f) + h_offset + h_resid(T, p) - (if phase <> Phase.gas then 
    h_resid(T, p0) else h_resid(
    T,
    p0,
    {1,-n_v[1]}));
  // The last two terms adjust for the actual pressure relative to the
  // reference.  In general, the lower limit of the integral of
  // (delh/delp)_T*dp is the reference pressure (p0).  However, if the
  // material is gaseous, then the reference is the corresponding ideal gas.
  // In that case, the lower limit of the real gas terms of the integral is
  // p=0, where a real gas behaves as an ideal gas.  See [Rao1997, p. 271].
end h;

FCSys.Characteristics.BaseClasses.Characteristic.s FCSys.Characteristics.BaseClasses.Characteristic.s

Specific entropy as a function of temperature and pressure

Information

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
PressureAbsolutepp0Pressure [M/(L.T2)]

Outputs

TypeNameDescription
NumberAbsolutesSpecific entropy [1]

Modelica definition

function s 
  "Specific entropy as a function of temperature and pressure"
  import FCSys.Utilities.Polynomial;
  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.PressureAbsolute p=p0 "Pressure";
  output Q.NumberAbsolute s "Specific entropy";

protected 

    "Return s0 as a function of T using one of the temperature intervals"
  annotation(derivative=ds0_i);
    import FCSys.Utilities.Polynomial.F;
    input Q.TemperatureAbsolute T "Temperature";
    input Integer i "Index of the temperature interval";
    output Q.NumberAbsolute s0 
      "Specific entropy at given temperature and reference pressure";

  algorithm 
    s0 := F(
        T,
        b_c[i, :],
        n_c - 1) + B_c[i, 2];
    // This is the integral of c0_p/T*dT up to T at p0 with the absolute
    // entropy at the lower bound [McBride2002, p. 2].

  end s0_i;

  "Derivative of s0_i"
    // Note:  This function is necessary for Dymola 7.4 to differentiate s().
    import FCSys.Utilities.Polynomial.f;
    input Q.TemperatureAbsolute T "Temperature";
    input Integer i "Index of the temperature interval";
    input Q.Temperature dT "Derivative of temperature";
    output Q.Number ds0 
      "Derivative of specific entropy at given temperature and reference pressure";

  algorithm 
    ds0 := f(
        T,
        b_c[i, :],
        n_c - 1)*dT;
  end ds0_i;

"Residual specific entropy for pressure adjustment for selected rows of b_v"
  annotation(derivative=ds_resid);
    input Q.TemperatureAbsolute T "Temperature";
    input Q.PressureAbsolute p "Pressure";
    input Integer rowLimits[2]={1,size(b_v, 1)} 
      "Beginning and ending indices of rows of b_v to be included";
    output Q.NumberAbsolute s_resid 
      "Integral of (dels/delp)_T*dp up to p with zero integration constant (for selected rows)";

  algorithm 
    s_resid := Polynomial.F(
        p,
        {Polynomial.f(
          T,
          b_v[i, :] .* {n_v[1] - n_v[2] + i - j for j in 1:size(b_v, 2)},
          n_v[2] - n_v[1] - i) for i in rowLimits[1]:rowLimits[2]},
        n_v[1] + rowLimits[1] - 1);
    // Note:  According to the Maxwell relations,
    // (dels/delp)_T = -(delv/delT)_p.

  end s_resid;

  "Derivative of s_resid"
    // Note:  This function is necessary for Dymola 7.4 to differentiate s().
    input Q.TemperatureAbsolute T "Temperature";
    input Q.PressureAbsolute p "Pressure";
    input Integer rowLimits[2]={1,size(b_v, 1)} 
      "Beginning and ending indices of rows of b_v to be included";
    input Q.Temperature dT "Derivative of temperature";
    input Q.Pressure dp "Derivative of pressure";
    output Q.Number ds_resid 
      "Derivative of integral of (dels/delp)_T*dp up to p with zero integration constant (for selected rows)";

  algorithm 
    ds_resid := Polynomial.dF(
        p,
        {Polynomial.df(
          T,
          b_v[i, :] .* {n_v[1] - n_v[2] + i - j for j in 1:size(b_v, 2)},
          n_v[2] - n_v[1] - i,
          dT) for i in rowLimits[1]:rowLimits[2]},
        n_v[1] + rowLimits[1] - 1,
        dp);
  end ds_resid;

algorithm 
  /*
  assert(T_lim_c[1] <= T and T <= T_lim_c[end], "Temperature " +
    String(T/U.K) + " K is out of bounds for " + name + " ([" + String(T_lim_c[1]
    /U.K) + ", " + String(T_lim_c[size(T_lim_c, 1)]/U.K) + "] K).");
  */
  // Note:  This is commented out so that the function can be inlined.
  s := smooth(1, sum(if (T_lim_c[i] <= T or i == 1) and (T < T_lim_c[i + 1] or 
    i == size(b_c, 1)) then s0_i(T, i) else 0 for i in 1:size(b_c, 1))) + 
    s_resid(T, p) - (if phase <> Phase.gas then s_resid(T, p0) else s_resid(
    T,
    p0,
    {1,-n_v[1]}));
  // The first term gives the specific entropy at the given temperature and
  // reference pressure (p0).  The following terms adjust for the actual
  // pressure (p) by integrating (dels/delp)_T*dp.  In general, the
  // integration is from p0 to p.  However, for gases the reference state
  // is the ideal gas at the reference pressure.  Therefore, the lower
  // integration limit for the higher-order real gas terms (i.e., terms with
  // power of p greater than -1) is p=0.  This is pressure limit at which a
  // real gas behaves as an ideal gas, and it implies that those terms are
  // zero.  See [Rao1997, p. 272].  Note that the first V_m inside the curly
  // brackets of the related eq. 1.47 in [Dymond2002, p. 17] should be a
  // subscript rather than a multiplicative factor.
  // Note:  If the phase is gas, the virial equation of state (as defined by
  // b_v and n_v) must include an ideal gas term (v = … + f(T)/p +
  // …).  Otherwise, an indexing error will occur.
end s;

FCSys.Characteristics.BaseClasses.Characteristic.zeta FCSys.Characteristics.BaseClasses.Characteristic.zeta

** (ζ) as a function of temperature

Information

**Continuity is a property is defined in FCSys resistivity to axial compression or material storage during transport.

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
DiffusivityMassSpecificzeta** [L2.M/(N.T)]

Modelica definition

replaceable function zeta 
  "** (ζ) as a function of temperature"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.DiffusivityMassSpecific zeta "**";

algorithm 
  zeta := m*D(T, v);

end zeta;

FCSys.Characteristics.BaseClasses.Characteristic.eta FCSys.Characteristics.BaseClasses.Characteristic.eta

Fluidity (η) as a function of temperature

Information

Fluidity is defined as the reciprocal of dynamic viscosity (see http://en.wikipedia.org/wiki/Viscosity#Fluidity).

Although specific volume is an input to this function, the result is independent of specific volume. According to Present [Present1958], this independence very accurately matches the measured fluidity of gases. However, the fluidity varies by species and generally falls more rapidly with temperature than indicated [Present1958, p. 41].

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.
Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
FluidityetaFluidity [L.T/M]

Modelica definition

replaceable function eta 
  "Fluidity (η) as a function of temperature"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  // Note:  Specific volume isn't used here but is included for generality.
  output Q.Fluidity eta "Fluidity";

algorithm 
  eta := alpha/(m*omega(T));
end eta;

FCSys.Characteristics.BaseClasses.Characteristic.theta FCSys.Characteristics.BaseClasses.Characteristic.theta

Thermal resistivity (θ) as a function of temperature and specific volume

Information

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.
Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
ResistivityThermalthetaThermal resistivity [L.T/N]

Modelica definition

replaceable function theta 
  "Thermal resistivity (θ) as a function of temperature and specific volume"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.ResistivityThermal theta "Thermal resistivity";

algorithm 
  theta := alpha/(c_v(T, p_Tv(T, v))*omega(T));
end theta;

FCSys.Characteristics.BaseClasses.Characteristic.tauprime FCSys.Characteristics.BaseClasses.Characteristic.tauprime

Phase change interval (τ′) as a function of temperature and specific volume

Information

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.

Also, it is assumed that the Einstein relation applies.

Please see [Davies2013, Ch. 3] for a derivation of the rate of phase change from kinetic theory.

Although specific volume is an input to this function, the result is independent of specific volume.

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
TimeAbsolutetauprimePhase change interval [T]

Modelica definition

replaceable function tauprime 
  "Phase change interval (τ′) as a function of temperature and specific volume"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.TimeAbsolute tauprime "Phase change interval";

algorithm 
  tauprime := v/(alpha*omega(T));
end tauprime;

FCSys.Characteristics.BaseClasses.Characteristic.mu FCSys.Characteristics.BaseClasses.Characteristic.mu

Mobility (μ) as a function of temperature and specific volume

Information

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.

Also, it is assumed that the Einstein relation applies.

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
MobilitymuMobility [N.T/M]

Modelica definition

replaceable function mu 
  "Mobility (μ) as a function of temperature and specific volume"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.Mobility mu "Mobility";

algorithm 
  mu := v/(m*alpha*omega(T));
end mu;

FCSys.Characteristics.BaseClasses.Characteristic.nu FCSys.Characteristics.BaseClasses.Characteristic.nu

Thermal independity (ν) as a function of temperature and specific volume

Information

Thermal independity describes the extent to which an exchange of thermal energy between species causes or requires a temperature difference.

This function is based on the kinetic theory of gases under the following assumptions [Present1958]:

  1. The particles are smooth and rigid but elastic spheres with identical radii. This is the "billiard-ball" assumption, and it implies that the collisions are instantaneous and conserve kinetic energy.
  2. Between collisions, particles have no influence on one another.
  3. The mean free path, or average distance a particle travels between collisions, is much larger than the diameter of a particle.
  4. The properties carried by a particle depend only on those of the last particle with which it collided.
  5. The speeds of the particles follow the Maxwell-Boltzmann distribution.

Also, it is assumed that the Einstein relation applies.

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

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT298.15*U.KTemperature [L2.M/(N.T2)]
VolumeSpecificv298.15*U.K/p0Specific volume [L3/N]

Outputs

TypeNameDescription
TimeAbsolutenuThermal independity [T]

Modelica definition

replaceable function nu 
  "Thermal independity (ν) as a function of temperature and specific volume"

  extends Modelica.Icons.Function;
  input Q.TemperatureAbsolute T=298.15*U.K "Temperature";
  input Q.VolumeSpecific v=298.15*U.K/p0 "Specific volume";
  output Q.TimeAbsolute nu "Thermal independity";

algorithm 
  nu := v/(c_p(T, p_Tv(T, v))*alpha*omega(T));
end nu;

FCSys.Characteristics.BaseClasses.Characteristic.h.dh0_i

Derivative of h0_i

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT Temperature [L2.M/(N.T2)]
Integeri Index of the temperature interval
TemperaturedT Derivative of temperature [L2.M/(N.T2)]

Outputs

TypeNameDescription
Potentialdh0Derivative of specific enthalpy at reference pressure [L2.M/(N.T2)]

Modelica definition

function dh0_i "Derivative of h0_i"
  // Note:  This function is necessary for Dymola 7.4 to differentiate h().
  import FCSys.Utilities.Polynomial;
  input Q.TemperatureAbsolute T "Temperature";
  input Integer i "Index of the temperature interval";
  input Q.Temperature dT "Derivative of temperature";
  output Q.Potential dh0 
    "Derivative of specific enthalpy at reference pressure";

algorithm 
  dh0 := Polynomial.f(
    T,
    b_c[i, :],
    n_c)*dT;

end dh0_i;

FCSys.Characteristics.BaseClasses.Characteristic.s.ds0_i

Derivative of s0_i

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT Temperature [L2.M/(N.T2)]
Integeri Index of the temperature interval
TemperaturedT Derivative of temperature [L2.M/(N.T2)]

Outputs

TypeNameDescription
Numberds0Derivative of specific entropy at given temperature and reference pressure [1]

Modelica definition

function ds0_i "Derivative of s0_i"
  // Note:  This function is necessary for Dymola 7.4 to differentiate s().
  import FCSys.Utilities.Polynomial.f;
  input Q.TemperatureAbsolute T "Temperature";
  input Integer i "Index of the temperature interval";
  input Q.Temperature dT "Derivative of temperature";
  output Q.Number ds0 
    "Derivative of specific entropy at given temperature and reference pressure";

algorithm 
  ds0 := f(
    T,
    b_c[i, :],
    n_c - 1)*dT;
end ds0_i;

FCSys.Characteristics.BaseClasses.Characteristic.s.ds_resid

Derivative of s_resid

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT Temperature [L2.M/(N.T2)]
PressureAbsolutep Pressure [M/(L.T2)]
IntegerrowLimits[2]{1,size(b_v, 1)}Beginning and ending indices of rows of b_v to be included
TemperaturedT Derivative of temperature [L2.M/(N.T2)]
Pressuredp Derivative of pressure [M/(L.T2)]

Outputs

TypeNameDescription
Numberds_residDerivative of integral of (dels/delp)_T*dp up to p with zero integration constant (for selected rows) [1]

Modelica definition

function ds_resid "Derivative of s_resid"
  // Note:  This function is necessary for Dymola 7.4 to differentiate s().
  input Q.TemperatureAbsolute T "Temperature";
  input Q.PressureAbsolute p "Pressure";
  input Integer rowLimits[2]={1,size(b_v, 1)} 
    "Beginning and ending indices of rows of b_v to be included";
  input Q.Temperature dT "Derivative of temperature";
  input Q.Pressure dp "Derivative of pressure";
  output Q.Number ds_resid 
    "Derivative of integral of (dels/delp)_T*dp up to p with zero integration constant (for selected rows)";

algorithm 
  ds_resid := Polynomial.dF(
    p,
    {Polynomial.df(
      T,
      b_v[i, :] .* {n_v[1] - n_v[2] + i - j for j in 1:size(b_v, 2)},
      n_v[2] - n_v[1] - i,
      dT) for i in rowLimits[1]:rowLimits[2]},
    n_v[1] + rowLimits[1] - 1,
    dp);
end ds_resid;