FCSys.Units

Constants and units of physical measure

Information

The information below has been updated and adapted from [Davies and Paredis, 2012]. That paper also offers suggestions as to how the approach might be better integrated in Modelica. For more information, please also see the documentation of the Quantities package.

Introduction and Overview:

Mathematical models of physical systems use variables to represent physical quantities. As stated by the Bureau International des Poids et Mesures (BIPM) [BIPM2006, p. 103]:

"The value of a quantity is generally expressed as the product of a number and a unit. The unit is simply a particular example of the quantity concerned which is used as a reference, and the number is the ratio of the value of the quantity to the unit."

In general, a unit may be the product of powers of other units, whether they are base units or units derived from the base units in the same manner.

In Modelica, a physical quantity is represented by a variable which is an instance of the Real type. Its value attribute is a number associated with the value of the quantity (not the value of the quantity itself, as will be shown). Usually the value attribute is not referenced explicitly because it is automatically returned when a variable is referenced. The unit attribute is a string that describes the unit by which the value of the quantity has been divided to arrive at the number. The displayUnit attribute (also a string) describes the unit by which the value of the quantity should be divided to arrive at the number as it is entered by or presented to the user. The Real type contains other attributes as well, including the quantity string.

The SIunits package of the Modelica Standard Library contains types that extend from the Real type. The type definitions modify the unit, displayUnit, and quantity attributes (among others) to represent various physical quantities. The unit and displayUnit attributes are based on the International System of Units (Système international d'unités, SI). The quantity string is the name of the physical quantity. For example, the Velocity type has a unit of "m/s" and a quantity of "Velocity". If an instance of Velocity has a value of one (v = 1), then it is meant that "the value of velocity is one meter per second." Again, the value attribute represents the number, or the value of the quantity divided by the unit, not the value of the quantity itself.

This apparent conflict is solved in FCSys by establishing units (including the meter and the second) as mathematical entities and writing v = 1 m/s (in code, v = 1*U.m/U.s or simply v = U.m/U.s). Here, the variable v directly represents the quantity. Its value attribute is truly the value of the quantity in the context of the statement by BIPM (above). One advantage is that unit conversion is handled naturally. The essence of unit conversion is the phrase "value of a quantity in a unit" typically means "value of a quantity divided by a unit." Continuing with the previous example, v is divided by m/s in order to display v in meters per second (as a number). If another unit of length like the foot is established by the appropriate relation (ft ≈ 0.3048 m) and v is divided by ft/s, the result is velocity in feet per second (∼3.2894). Some units such as °C, Pag, and dB involve offsets or nonlinear transformations between the value of the quantity and the number; these are described by functions besides simple division.

As another example, frequency is sometimes represented by a variable in hertz or cycles per second (e.g., ν) and other times by a variable in radians per second (e.g., ω). If the variable represents the quantity directly, then there is no need to specify which units it is in. The units are included; they have not been factored out by division (or another function). A common variable (e.g., f) can be used in both cases, which simplifies and standardizes the equations of a model. The forms are algebraically equivalent due to the relationships among units (e.g., 1 cycle = 2π rad).

Method:

In FCSys, each unit is a constant quantity. The values of the units, like other quantities, is the product of a number and a unit. Therefore, units may be derived from other units (e.g., cycle = 2π rad). This recursive definition leaves several units (in SI, 7) that are locally independent and must be established universally. These base units are established by the "particular example of the quantity concerned which is used as a reference" quoted previously [BIPM2006]. The choice of the base units is somewhat arbitrary [Fritzson2004, p. 375], but regardless, there are a number of units that must be defined by example.

If only SI will be used, then it is easiest to set each of the base units of SI equal to one—the meter (m), kilogram (kg), second (s), ampere (A), kelvin (K), mole (mol), and candela (cd). This is implicitly the case in the SIunits package, but again, it hardly captures the idea that the value of a quantity is the product of a number and a unit.

Instead, in FCSys, the base units are established by universal physical constants (except the candela, which is physically arbitrary). The "particular example of the quantity" [BIPM2006] is an experiment that yields precise and universally repeatable results in determining a constant rather than a prototype (e.g., the International Prototype of the Kilogram) which is carefully controlled and distributed via replicas. This method of defining the base units from measured physical quantities (rather than vice versa) is natural and reflects the way that standards organizations (e.g., NIST) define units. A system of units is considered to be natural if all of its base units are established by universal physical constants. Often, those universal constants are defined to be equal to one, but the values can be chosen to scale the numerical values of variables during simulation.

There are physical systems where typical quantities are many orders of magnitude larger or smaller than the related product of powers of base SI units (e.g., the domains of astrophysics and atomic physics). In modeling and simulating those systems, it may be advantageous to choose appropriately small or large values (respectively) for the corresponding base units so that the product of the number (large or small in magnitude) and the unit (small or large, respectively) is well-scaled. Products of this type are often involved in initial conditions or parameter expressions, which are not time-varying. Therefore, the number and the unit can be multiplied before the dynamic simulation. During the simulation, only the value is important. After the simulation, the trajectory of the value may be divided by the unit for display. This scaling is usually unnecessary due to the wide range and appropriate distribution of the real numbers that are representable in floating point. The Modelica language specification recommends that floating point numbers be represented in at least IEEE double precision, which covers magnitudes from ∼2.225×10-308 to ∼1.798×10308 [Modelica2010, p. 13]. However, in some cases it may be preferable to scale the units and use lower precision for the sake of computational performance. There are fields of research where, even today, simulations are sometimes performed in single precision [Brown2011, Hess2008] and where scaling is a concern [Rapaport2004, p. 29].

The method is neutral with regards to not only the values of the base units, but also the choice of the base units and even the number of base units. This is an advantage because many systems of units besides SI are used in science and engineering. As mentioned previously, the choice of base units is somewhat arbitrary, and different systems of units are based on different choices. Some systems of units have fewer base units (lower rank) than SI, since additional constraints are added that exchange base units for derived units. For example, the Planck, Stoney, Hartree, and Rydberg systems of units define the Boltzmann constant to be equal to one (kB = 1) [http://en.wikipedia.org/wiki/Natural_units]. The unit K is eliminated [Greiner1995, p. 386] or, more precisely, considered a derived unit instead of a base unit. In SI, the kelvin would be derived from the units kilogram, meter, and second (K ≈ 1.381×10-23 kg m2/s2).

There are six independent units and constants in the Units package (see Units.Bases), but SI has seven base units (m, kg, s, A, K, mol, and cd). In FCSys, two additional constraints are imposed in order to simplify the model equations and allow electrons and chemical species to be to represented by the same base Species model. First, the Faraday constant (kF or 96485.3399 C/mol) is normalized to one. This implies that the mole (mol) is proportional to the coulomb (C), which is considered a number of reference particles given a charge number of one. Also, the gas constant (R or 8.314472 J/(mol K)) is normalized to one. Therefore, the kelvin (K) is proportional to the volt (V or J/C). In addition, the radian (rad) is defined as a base constant. However, it must be set equal to one in the current specification of the International System of Units (SI) [BIPM2006].

Implementation:

The units and constants are defined as variables in this Units package. Each is a constant of the appropriate type from the Quantities package. The first section of this package establishes mathematical constants. The next section establishes the base constants and units, which grouped in a replaceable subpackage. The third section establishes the constants and units which may be derived from the base units and constants using accepted empirical relations. The rest of the code establishes the SI prefixes and the remaining derived units and constants. The SI prefixes are included in their unabbreviated form in order to avoid naming conflicts. All of the primary units of SI are included (Tables 1 and 3 of [BIPM2006]) except for °C, since it involves an offset. Other units such as the atmosphere (atm) are included for convenience. Some units that include prefixes are defined as well (e.g., kg, mm, and kPa). However, most prefixes must be given as explicit factors (e.g., U.kilo*U.m).

Besides the units and constants, this package also contains functions (e.g., to_degC) that convert quantities from the unit system defined in FCSys to quantities expressed in units. These functions are included for units that involve offsets. For conversions that require just a scaling factor, it is simpler to use the units directly. For example, to convert from potential in volts use v = v_V*U.V, where v is potential and v_V is potential expressed in volts.

This package (Units) is abbreviated as U for convenience throughout the rest of FCSys. For example, an initial pressure might be defined as pIC = U.atm.

An instance of the Environment model is usually included at the top level of a model. It records the base units and constants so that it is possible to re-derive all of the other units and constants. This is important in order to properly interpret simulation results if the base units and constants are later re-adjusted.

The Units.setup function establishes unit conversions using the values of the units, constants, and prefixes. These unit conversions may include offsets. The function also sets the default display units. It is automatically called when FCSys is loaded from the load.mos script. It can also be called manually from the "Re-initialize the units" command available in the Modelica development environment from the Units package or any subpackage. A spreadsheet (Resources/quantities.xls) is available to help maintain the quantities, default units, and the setup function.

The values of the units, constants, and prefixes can be evaluated by translating the Units.Examples.Evaluate model. This defines the values in the workspace of the Modelica development environment. For convenience, the load.mos script automatically does this and saves the result as "units.mos" in the working directory.

Where the der operator is used in models, it is explicitly divided by the unit second (e.g., der(x)/U.s). This is necessary because the global variable time is in seconds (i.e., time is a number, not a quantity).

Although it is not necessary due to the acausal nature of Modelica, the declarations in this package are sorted so that they can be easily ported to imperative or causal languages (e.g., Python and C). In fact, this has been done in the included FCRes module for plotting and analysis.

Extends from Modelica.Icons.Package (Icon for standard packages).

Package Content

NameDescription
FCSys.Units.setup setup Establish conversions to display quantities in units
FCSys.Units.Examples Examples Examples
FCSys.Units.Bases Bases Sets of base constants and units
FCSys.Units.from_degC from_degC Convert from temperature in degree Celsius to temperature as a quantity
FCSys.Units.to_degC to_degC Convert from temperature as a quantity to temperature in degree Celsius
FCSys.Units.from_kPag from_kPag Convert from gauge pressure in kilopascals to absolute pressure as a quantity
FCSys.Units.to_kPag to_kPag Convert from absolute pressure as a quantity to gauge pressure in kilopascals
pi=2*acos(0)pi (π)
baseScalable base constants and units
rad=base.radradian
R_inf=base.R_infRydberg constant (R)
c=base.cspeed of light in vacuum (c)
k_J=base.k_JJosephson constant (kJ)
R_K=base.R_Kvon Klitzing constant (RK)
'cd'=base.'cd'candela
k_F=base.k_FFaraday constant (kF)
R=base.Rgas constant
m=10973731.568539*rad/R_infmeter
s=299792458*m/csecond
Wb=483597.870e9/k_Jweber
S=25812.8074434/R_Ksiemen
mol=96485.3365*Wb*S/k_Fmole
K=8.3144621*(Wb*rad)^2*S/(s*mol*R)kelvin
V=Wb*rad/svolt
A=V*Sampere
C=A*scoulomb
J=V*Cjoule
Gy=(m/s)^2gray
kg=J/Svkilogram
yotta=1e24yotta (Y)
zetta=1e21zetta (Z)
exa=1e18exa (E)
peta=1e15peta (P)
tera=1e12tera (T)
giga=1e9giga (G)
mega=1e6mega (M)
kilo=1e3kilo (k)
hecto=1e2hecto (h)
deca=1e1deca (da)
deci=1e-1deci (d)
centi=1e-2centi (c)
milli=1e-3milli (m)
micro=1e-6micro (u)
nano=1e-9nano (n)
pico=1e-12pico (p)
femto=1e-15femto (f)
atto=1e-18atto (a)
zepto=1e-21zepto (z)
yocto=1e-24yocto (y)
cyc=2*pi*radcycle
Hz=cyc/shertz
sr=rad^2steradian
N=J/mnewton
Pa=N/m^2pascal
W=J/swatt
F=C/Vfarad
ohm=1/Sohm (Ω)
H=V*s/Ahenry
T=Wb/m^2tesla
lm='cd'*srlumen
lx=lm/m^2lux
Bq=Hzbecquerel
Sv=Gysievert
kat=mol/skatal
g=kg/kilogram
min=60*sminute
hr=60*minhour
day=24*hrday
degree=2*pi*rad/360degree (°)
L=(deci*m)^3liter (L or l)
G_0=2/R_Kconductance quantum (G0)
Phi_0=1/k_Jmagnetic flux quantum (Φ0)
q=G_0*Phi_0elementary charge
h=2*q*Phi_0Planck constant
alpha=pi*1e-7*c*s*G_0/(m*S)fine-structure constant (α)
Z_0=2*R_K*alphacharacteristic impedance of vacuum (Z0)
mu_0=Z_0/cmagnetic constant (μ0)
epsilon_0=1/(Z_0*c)electric constant (ε0)
k_A=mu_0/(4*pi)magnetic force constant (kA)
k_e=k_A*c^2Coulomb constant (ke)
E_h=2*R_inf*h*cHartree energy (Eh)
eV=q*Velectron volt
N_A=k_F/qAvogadro constant (NA)
k_B=R/N_ABoltzmann constant (kB)
c_1=cyc*h*c^2first radiation constant (c1)
c_2=h*c/k_Bsecond radiation constant (c2)
c_3_lambda=c_2/4.965114231744276Wien wavelength displacement law constant (c3 λ)
c_3_f=2.821439372122079*k_B/hWien frequency displacement law constant (c3 f)
sigma=2*pi*(k_B*pi)^4/(15*(h*rad)^3*c^2)Stefan-Boltzmann constant (σ)
bar=1e5*Pabar
angstrom=0.1*nano*mangstrom (Å)
atm=101325*Paatmosphere
kPa=kilo*Pakilopascal
kJ=kilo*Jkilojoule
cm=centi*mcentimeter
mm=milli*mmillimeter
mA=milli*Amilliampere
um=micro*mmicrometer
ms=milli*smillisecond
'%'=centipercent (%)
M=U.mol/U.Lmolar
cc=U.cm^3cubic centimeter

Types and constants

final constant Q.Number pi=2*acos(0) "pi (π)";
replaceable constant Bases.LH base constrainedby Bases.SImols 
  "Scalable base constants and units";
final constant Q.Angle rad=base.rad "radian";
final constant Q.Wavenumber R_inf=base.R_inf 
  "Rydberg constant (R)";
final constant Q.Velocity c=base.c 
  "speed of light in vacuum (c)";
final constant Q.MagneticFluxReciprocal k_J=base.k_J 
  "Josephson constant (kJ)";
final constant Q.ResistanceElectrical R_K=base.R_K 
  "von Klitzing constant (RK)";
final constant Q.PowerRadiant 'cd'=base.'cd' "candela";
final constant Q.Number k_F=base.k_F 
  "Faraday constant (kF)";
final constant Q.Number R=base.R "gas constant";
constant Q.Length m=10973731.568539*rad/R_inf "meter";
constant Q.Time s=299792458*m/c "second";
constant Q.MagneticFlux Wb=483597.870e9/k_J "weber";
constant Q.ConductanceElectrical S=25812.8074434/R_K "siemen";
constant Q.Amount mol=96485.3365*Wb*S/k_F "mole";
constant Q.Potential K=8.3144621*(Wb*rad)^2*S/(s*mol*R) "kelvin";
final constant Q.Potential V=Wb*rad/s "volt";
final constant Q.Current A=V*S "ampere";
final constant Q.Amount C=A*s "coulomb";
final constant Q.Energy J=V*C "joule";
final constant Q.Velocity2 Gy=(m/s)^2 "gray";
final constant Q.Mass kg=J/Sv "kilogram ";
final constant Q.Number yotta=1e24 "yotta (Y)";
final constant Q.Number zetta=1e21 "zetta (Z)";
final constant Q.Number exa=1e18 "exa (E)";
final constant Q.Number peta=1e15 "peta (P)";
final constant Q.Number tera=1e12 "tera (T)";
final constant Q.Number giga=1e9 "giga (G)";
final constant Q.Number mega=1e6 "mega (M)";
final constant Q.Number kilo=1e3 "kilo (k)";
final constant Q.Number hecto=1e2 "hecto (h)";
final constant Q.Number deca=1e1 "deca (da)";
final constant Q.Number deci=1e-1 "deci (d)";
final constant Q.Number centi=1e-2 "centi (c)";
final constant Q.Number milli=1e-3 "milli (m)";
final constant Q.Number micro=1e-6 "micro (u)";
final constant Q.Number nano=1e-9 "nano (n)";
final constant Q.Number pico=1e-12 "pico (p)";
final constant Q.Number femto=1e-15 "femto (f)";
final constant Q.Number atto=1e-18 "atto (a)";
final constant Q.Number zepto=1e-21 "zepto (z)";
final constant Q.Number yocto=1e-24 "yocto (y)";
final constant Q.Angle cyc=2*pi*rad "cycle";
final constant Q.Frequency Hz=cyc/s "hertz";
final constant Q.Angle2 sr=rad^2 "steradian";
final constant Q.Force N=J/m "newton";
final constant Q.Pressure Pa=N/m^2 "pascal";
final constant Q.Power W=J/s "watt";
final constant Q.Capacitance F=C/V "farad";
final constant Q.ResistanceElectrical ohm=1/S "ohm (Ω)";
final constant Q.Inductance H=V*s/A "henry";
final constant Q.MagneticFluxAreic T=Wb/m^2 "tesla";
final constant Q.Power lm='cd'*sr "lumen";
final constant Q.PowerAreic lx=lm/m^2 "lux";
final constant Q.Frequency Bq=Hz "becquerel";
final constant Q.Velocity2 Sv=Gy "sievert";
final constant Q.Current kat=mol/s "katal";
final constant Q.Mass g=kg/kilo "gram";
final constant Q.Time min=60*s "minute";
final constant Q.Time hr=60*min "hour";
final constant Q.Time day=24*hr "day";
final constant Q.Angle degree=2*pi*rad/360 "degree (°)";
final constant Q.Volume L=(deci*m)^3 "liter (L or l)";
final constant Q.ConductanceElectrical G_0=2/R_K 
  "conductance quantum (G0)";
final constant Q.MagneticFlux Phi_0=1/k_J 
  "magnetic flux quantum (Φ0)";
final constant Q.Amount q=G_0*Phi_0 "elementary charge";
final constant Q.MomentumRotational h=2*q*Phi_0 "Planck constant";
final constant Q.Number alpha=pi*1e-7*c*s*G_0/(m*S) 
  "fine-structure constant (α)";
final constant Q.ResistanceElectrical Z_0=2*R_K*alpha 
  "characteristic impedance of vacuum (Z0)";
final constant Q.Permeability mu_0=Z_0/c 
  "magnetic constant (μ0)";
final constant Q.Permittivity epsilon_0=1/(Z_0*c) 
  "electric constant (ε0)";
final constant Q.Permeability k_A=mu_0/(4*pi) 
  "magnetic force constant (kA)";
final constant Q.PermittivityReciprocal k_e=k_A*c^2 
  "Coulomb constant (ke)";
final constant Q.Energy E_h=2*R_inf*h*c 
  "Hartree energy (Eh)";
final constant Q.Energy eV=q*V "electron volt";
final constant Q.AmountReciprocal N_A=k_F/q 
  "Avogadro constant (NA)";
final constant Q.Amount k_B=R/N_A 
  "Boltzmann constant (kB)";
final constant Q.PowerArea c_1=cyc*h*c^2 
  "first radiation constant (c1)";
final constant Q.PotentialPerWavenumber c_2=h*c/k_B 
  "second radiation constant (c2)";
final constant Q.PotentialPerWavenumber c_3_lambda=c_2/4.965114231744276 
  "Wien wavelength displacement law constant (c3 λ)";
final constant Q.MagneticFluxReciprocal c_3_f=2.821439372122079*k_B/h 
  "Wien frequency displacement law constant (c3 f)";
final constant Q.PowerAreicPerPotential4 sigma=2*pi*(k_B*pi)^4/(15*(h*rad)^3*c^2)
  "Stefan-Boltzmann constant (σ)";
final constant Q.Pressure bar=1e5*Pa "bar";
final constant Q.Length angstrom=0.1*nano*m "angstrom (Å)";
final constant Q.Pressure atm=101325*Pa "atmosphere";
final constant Q.Pressure kPa=kilo*Pa "kilopascal";
final constant Q.Energy kJ=kilo*J "kilojoule";
final constant Q.Length cm=centi*m "centimeter";
final constant Q.Length mm=milli*m "millimeter";
final constant Q.Current mA=milli*A "milliampere";
final constant Q.Length um=micro*m "micrometer";
final constant Q.Time ms=milli*s "millisecond";
final constant Q.Number '%'=centi "percent (%)";
final constant Q.Density M=U.mol/U.L "molar";
final constant Q.Volume cc=U.cm^3 "cubic centimeter";

FCSys.Units.setup FCSys.Units.setup

Establish conversions to display quantities in units

Information

This function has no inputs or outputs. It is essentially a script. The defineDefaultDisplayUnit and defineUnitConversion functions used here are not defined in the Modelica language (as of version 3.3) but are recognized by Dymola.

For more information, please see the documentation for the Units package.

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

Modelica definition

function setup "Establish conversions to display quantities in units"
  import Modelica.Utilities.Streams.print;
  extends Modelica.Icons.Function;

algorithm 
  print("Establishing display units…");

  // ------------------------------------------------------------------------
  // Default display units
  // ------------------------------------------------------------------------
  // If units other than those in the displayUnit attribute of the
  // quantities in FCSys.Quantities should be used by default, then specify
  // them here.  Be sure that the desired unit is included in a
  // defineUnitConversion command below.
  // Generated from FCSys/Resources/quantities.xls, 2013-11-4
  defineDefaultDisplayUnit("L2.M/(A2.T3)", "cd") "Radiant power";
  defineDefaultDisplayUnit("L.T2/M", "1/kPa") "Reciprocal of pressure";
  defineDefaultDisplayUnit("1/L", "1/m") "Reciprocal of length";
  defineDefaultDisplayUnit("1/N", "1/mol") "Reciprocal of amount";
  defineDefaultDisplayUnit("L.T/M", "1/(Pa.s)") "Fluidity";
  defineDefaultDisplayUnit("A.N.T/(L2.M)", "1/Wb") 
    "Reciprocal of magnetic flux";
  defineDefaultDisplayUnit("N/s", "A") "for derivative of amount in Dymola";
  defineDefaultDisplayUnit("N/T", "A") "Current";
  defineDefaultDisplayUnit("N/(L2.T)", "A/cm2") "Areic current";
  defineDefaultDisplayUnit("N/(T.s)", "A/s") 
    "for derivative of current in Dymola";
  defineDefaultDisplayUnit("N/T2", "A/s") "Rate of current";
  defineDefaultDisplayUnit("L.N/T", "A.m") "VelocityAmount";
  defineDefaultDisplayUnit("N", "C") "Amount";
  defineDefaultDisplayUnit("N/L3", "C/cc") "Density";
  defineDefaultDisplayUnit("N/(L3.T)", "C/(cc.s)") "Rate of concentration";
  defineDefaultDisplayUnit("L3", "cc") "Volume";
  defineDefaultDisplayUnit("L3/N", "cc/C") "Specific volume";
  defineDefaultDisplayUnit("L3/(N.s)", "cc/(C.s)") 
    "for derivative of specific volume in Dymola";
  defineDefaultDisplayUnit("L3/(N.T)", "cc/(C.s)") "Rate of specific volume";
  defineDefaultDisplayUnit("L", "cm") "Length";
  defineDefaultDisplayUnit("L.T/N", "cm/A") "Resistivity";
  defineDefaultDisplayUnit("L/T", "cm/s") "Velocity";
  defineDefaultDisplayUnit("L/(T.s)", "cm/s2") 
    "for derivative of velocity in Dymola";
  defineDefaultDisplayUnit("L2", "cm2") "Area";
  defineDefaultDisplayUnit("L2/N", "cm2/C") "Specific area";
  defineDefaultDisplayUnit("L2/T", "cm2/s") "Diffusivity";
  defineDefaultDisplayUnit("N.T/M", "cm2/(V.s)") "Mobility";
  defineDefaultDisplayUnit("N2.T2/(L3.M)", "F/m") "Permittivity";
  defineDefaultDisplayUnit("M", "g") "Mass";
  defineDefaultDisplayUnit("M/L3", "g/cc") "Volumic mass";
  defineDefaultDisplayUnit("M/N", "g/mol") "Specific mass";
  defineDefaultDisplayUnit("M/s", "g/s") "for derivative of mass in Dymola";
  defineDefaultDisplayUnit("L.M/N2", "H/m") "Permeability";
  defineDefaultDisplayUnit("L2.M/T2", "J") "Energy";
  defineDefaultDisplayUnit("L2.M/(A.T)", "J.s/rad") "Rotational momentum";
  defineDefaultDisplayUnit("T/N", "K/W") "Thermal resistance";
  defineDefaultDisplayUnit("M/(L.T2)", "kPa") "Pressure";
  defineDefaultDisplayUnit("L3/T", "L/min") "Rate of volume";
  defineDefaultDisplayUnit("L3.M/(N2.T2)", "m/H") "Reciprocal of permittivity";
  defineDefaultDisplayUnit("L/N", "m/mol") "Specific length";
  defineDefaultDisplayUnit("L/T2", "m/s2") "Acceleration";
  defineDefaultDisplayUnit("N2.T2/(L2.M)", "uF") "Capacitance";
  defineDefaultDisplayUnit("L2.M/N2", "uH") "Inductance";
  defineDefaultDisplayUnit("L.M/T2", "N") "Force";
  defineDefaultDisplayUnit("L.M/(N.T)", "N/A") "Continuity";
  defineDefaultDisplayUnit("L2.M/(N2.T)", "ohm") "Electrical resistance";
  defineDefaultDisplayUnit("M/(L.T2.s)", "Pa/s") 
    "for derivative of pressure in Dymola";
  defineDefaultDisplayUnit("M/(L.T3)", "Pa/s") "Rate of pressure";
  defineDefaultDisplayUnit("M/(T2.L.s)", "Pa/s") 
    "for derivative of pressure in Dymola";
  defineDefaultDisplayUnit("A", "rad") "Angle";
  defineDefaultDisplayUnit("A/L", "rad/m") "Wavenumber";
  defineDefaultDisplayUnit("A/T", "rad/s") "Frequency";
  defineDefaultDisplayUnit("N2.T/(L2.M)", "S") "Electrical conductance";
  defineDefaultDisplayUnit("T", "s") "Time";
  defineDefaultDisplayUnit("N2.T/(L3.M)", "S/m") "Electrical conductivity";
  defineDefaultDisplayUnit("T/L", "s/m") "Lineic time";
  defineDefaultDisplayUnit("A2", "sr") "Solid angle";
  defineDefaultDisplayUnit("L2/T2", "Sv") "Squared velocity";
  defineDefaultDisplayUnit("M/(A.N.T)", "T") "Areic magnetic flux";
  defineDefaultDisplayUnit("L2.M/(N.T2)", "V") "Potential";
  defineDefaultDisplayUnit("L.M/(N.T2)", "V/m") "Specific force";
  defineDefaultDisplayUnit("L2.M/(N.T2.s)", "V/s") 
    "for derivative of potential in Dymola";
  defineDefaultDisplayUnit("L2.M/(N.T3)", "V/s") "Rate of potential";
  defineDefaultDisplayUnit("L3.M/(A.N.T2)", "V.m/rad") 
    "Potential per wavenumber";
  defineDefaultDisplayUnit("L2.M/T3", "W") "Power";
  defineDefaultDisplayUnit("M/T3", "W/m2") "Areic power";
  defineDefaultDisplayUnit("M.T5/L8", "W/(m2.K4)") 
    "Areic power per 4th power of potential";
  defineDefaultDisplayUnit("L4.M/T3", "W.m2") "Power times area";
  defineDefaultDisplayUnit("L2.M/(A.N.T)", "Wb") "Magnetic flux";

  // ------------------------------------------------------------------------
  // Conversions to display quantities in units
  // ------------------------------------------------------------------------
  defineUnitConversion(
    "1",
    "%",
    1/'%') "Number";
  defineUnitConversion(
    "L2.M/(A2.T3)",
    "cd",
    1/'cd') "Radiant power";
  defineUnitConversion(
    "1/N",
    "1/C",
    C) "Reciprocal of amount";
  defineUnitConversion(
    "L.T2/M",
    "1/kPa",
    kilo*Pa) "Reciprocal of pressure";
  defineUnitConversion(
    "1/L",
    "1/m",
    m) "Reciprocal of length";
  defineUnitConversion(
    "1/N",
    "1/mol",
    mol) "Reciprocal of amount";
  defineUnitConversion(
    "L.T2/M",
    "1/Pa",
    Pa) "Reciprocal of pressure";
  defineUnitConversion(
    "L.T/M",
    "1/(Pa.s)",
    Pa*s) "Fluidity";
  defineUnitConversion(
    "A.N.T/(L2.M)",
    "1/Wb",
    Wb) "Reciprocal of magnetic flux";
  defineUnitConversion(
    "N/s",
    "A",
    1/C) "for derivative of amount in Dymola";
  defineUnitConversion(
    "N/T",
    "A",
    1/A) "Current";
  defineUnitConversion(
    "N/(L2.T)",
    "A/cm2",
    cm^2/A) "Areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "A/cm2",
    cm^2/A) "Absolute areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "A/m2",
    m^2/A) "Areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "A/m2",
    m^2/A) "Absolute areic current";
  defineUnitConversion(
    "N/(T.s)",
    "A/s",
    1/A) "for derivative of current in Dymola";
  defineUnitConversion(
    "N/T2",
    "A/s",
    s/A) "Rate of current";
  defineUnitConversion(
    "L.N/T",
    "A.m",
    1/(A*m)) "VelocityAmount";
  defineUnitConversion(
    "M/(L.T2)",
    "atm",
    1/atm) "Pressure";
  defineUnitConversion(
    "M/(L.T2)",
    "bar",
    1/bar) "Pressure";
  defineUnitConversion(
    "N",
    "C",
    1/C) "Amount";
  defineUnitConversion(
    "N/L3",
    "C/cc",
    cc/C) "Density";
  defineUnitConversion(
    "N/(L3.T)",
    "C/(cc.s)",
    cc*s/C) "Rate of concentration";
  defineUnitConversion(
    "N/L3",
    "C/m3",
    m^3/C) "Density";
  defineUnitConversion(
    "N/(L3.T)",
    "C/(m3.s)",
    m^3*s/C) "Rate of concentration";
  defineUnitConversion(
    "N.T/M",
    "C.s/g",
    g/(C*s)) "Mobility";
  defineUnitConversion(
    "L3",
    "cc",
    1/cc) "Volume";
  defineUnitConversion(
    "L3/N",
    "cc/C",
    C/cc) "Specific volume";
  defineUnitConversion(
    "L3/(N.s)",
    "cc/(C.s)",
    C*s/cc) "for derivative of specific volume in Dymola";
  defineUnitConversion(
    "L3/(N.T)",
    "cc/(C.s)",
    C*s/cc) "Rate of specific volume";
  defineUnitConversion(
    "L3/s",
    "cc/s",
    1/(centi*m)^3) "for derivative of volume in Dymola";
  defineUnitConversion(
    "L3/T",
    "cc/s",
    s/cc) "Rate of volume";
  defineUnitConversion(
    "L",
    "cm",
    1/cm) "Length";
  defineUnitConversion(
    "L.T/N",
    "cm/A",
    A/cm) "Resistivity";
  defineUnitConversion(
    "L/T",
    "cm/s",
    s/cm) "Velocity";
  defineUnitConversion(
    "L/(T.s)",
    "cm/s2",
    s/cm) "for derivative of velocity in Dymola";
  defineUnitConversion(
    "L/T2",
    "cm/s2",
    s^2/cm) "Acceleration";
  defineUnitConversion(
    "L.T/M",
    "cm.s/g",
    g/(cm*s)) "Fluidity";
  defineUnitConversion(
    "L2",
    "cm2",
    1/cm^2) "Area";
  defineUnitConversion(
    "L2/N",
    "cm2/C",
    C/cm^2) "Specific area";
  defineUnitConversion(
    "L2/N",
    "cm2/mol",
    mol/cm^2) "Specific area";
  defineUnitConversion(
    "L2/T",
    "cm2/s",
    s/cm^2) "Diffusivity";
  defineUnitConversion(
    "N.T/M",
    "cm2/(V.s)",
    V*s/cm^2) "Mobility";
  defineUnitConversion(
    "A/L",
    "cyc/m",
    m/cyc) "Wavenumber";
  defineUnitConversion(
    "T",
    "day",
    1/day) "Time";
  defineUnitConversion(
    "A",
    "degree",
    1/degree) "Angle";
  defineUnitConversion(
    "N2.T2/(L2.M)",
    "F",
    1/F) "Capacitance";
  defineUnitConversion(
    "N2.T2/(L3.M)",
    "F/m",
    m/F) "Permittivity";
  defineUnitConversion(
    "M",
    "g",
    1/g) "Mass";
  defineUnitConversion(
    "M/L3",
    "g/cc",
    cc/g) "Volumic mass";
  defineUnitConversion(
    "M/N",
    "g/mol",
    mol/g) "Specific mass";
  defineUnitConversion(
    "M/s",
    "g/s",
    s/g) "for derivative of mass in Dymola";
  defineUnitConversion(
    "L2.M/N2",
    "H",
    1/H) "Inductance";
  defineUnitConversion(
    "L.M/N2",
    "H/m",
    m/H) "Permeability";
  defineUnitConversion(
    "T",
    "hr",
    1/hr) "Time";
  defineUnitConversion(
    "A/T",
    "Hz",
    1/Hz) "Frequency";
  defineUnitConversion(
    "L2.M/T2",
    "J",
    1/J) "Energy";
  defineUnitConversion(
    "N",
    "J/K",
    K/J) "Amount";
  defineUnitConversion(
    "N/L3",
    "J/(m3.K)",
    m^3*K/J) "Density";
  defineUnitConversion(
    "L2.M/(N.T2)",
    "J/mol",
    mol/J) "Potential";
  defineUnitConversion(
    "1",
    "J/(mol.K)",
    mol*K/J) "Absolute number";
  defineUnitConversion(
    "L2.M/(A.T)",
    "J.s/rad",
    rad/(J*s)) "Rotational momentum";
  defineUnitConversion(
    "L2.M/(N.T2)",
    "K",
    1/K) "Absolute potential";
  defineUnitConversion(
    "L2.M/(N.T2.s)",
    "K/s",
    1/K) "for derivative of potential in Dymola";
  defineUnitConversion(
    "L2.M/(N.T3)",
    "K/s",
    s/K) "Rate of potential";
  defineUnitConversion(
    "T/N",
    "K/W",
    W/K) "Thermal resistance";
  defineUnitConversion(
    "L3.M/(A.N.T2)",
    "K.m/rad",
    rad/(K*m)) "Potential per wavenumber";
  defineUnitConversion(
    "N/s",
    "kat",
    1/mol) "for derivative of amount in Dymola";
  defineUnitConversion(
    "N/T",
    "kat",
    1/kat) "Current";
  defineUnitConversion(
    "N/(L2.T)",
    "kat/cm2",
    cm^2/kat) "Areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "kat/cm2",
    cm^2/kat) "Absolute areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "kat/m2",
    m^2/kat) "Areic current";
  defineUnitConversion(
    "N/(L2.T)",
    "kat/m2",
    m^2/kat) "Absolute areic current";
  defineUnitConversion(
    "N/(T.s)",
    "kat/s",
    1/kat) "for derivative of current in Dymola";
  defineUnitConversion(
    "N/T2",
    "kat/s",
    s/kat) "Rate of current";
  defineUnitConversion(
    "M",
    "kg",
    1/kg) "Mass";
  defineUnitConversion(
    "M/L3",
    "kg/m3",
    m^3/kg) "Volumic mass";
  defineUnitConversion(
    "M/N",
    "kg/mol",
    mol/kg) "Specific mass";
  defineUnitConversion(
    "M/s",
    "kg/s",
    s/kg) "for derivative of mass in Dymola";
  defineUnitConversion(
    "M/(L.T2)",
    "kPa",
    1/(kilo*Pa)) "Pressure";
  defineUnitConversion(
    "L2.M/(N.T2)",
    "kJ/mol",
    mol/kJ) "Potential";
  defineUnitConversion(
    "L3",
    "L",
    1/L) "Volume";
  defineUnitConversion(
    "L3/s",
    "L/min",
    min/(L*s)) "for derivative of volume in Dymola";
  defineUnitConversion(
    "L3/T",
    "L/min",
    min/L) "Rate of volume";
  defineUnitConversion(
    "L2.M/T3",
    "lm",
    1/lm) "Power";
  defineUnitConversion(
    "M/T3",
    "lm/m2",
    m^2/lm) "Areic power";
  defineUnitConversion(
    "N/L3",
    "M",
    1/M) "Density";
  defineUnitConversion(
    "L",
    "m",
    1/m) "Length";
  defineUnitConversion(
    "L.T/N",
    "m/A",
    A/m) "Resistivity";
  defineUnitConversion(
    "L/N",
    "m/C",
    C/m) "Specific length";
  defineUnitConversion(
    "L3.M/(N2.T2)",
    "m/H",
    H/m) "Reciprocal of permittivity";
  defineUnitConversion(
    "L/N",
    "m/mol",
    mol/m) "Specific length";
  defineUnitConversion(
    "N/(L3.T)",
    "M/s",
    s/M) "Rate of concentration";
  defineUnitConversion(
    "L/T",
    "m/s",
    s/m) "Velocity";
  defineUnitConversion(
    "L/(T.s)",
    "m/s2",
    s/m) "for derivative of velocity in Dymola";
  defineUnitConversion(
    "L/T2",
    "m/s2",
    s^2/m) "Acceleration";
  defineUnitConversion(
    "L.T/N",
    "m.K/W",
    W/(m*K)) "Resistivity";
  defineUnitConversion(
    "L2",
    "m2",
    1/m^2) "Area";
  defineUnitConversion(
    "L2/N",
    "m2/C",
    C/m^2) "Specific area";
  defineUnitConversion(
    "L2/N",
    "m2/mol",
    mol/m^2) "Specific area";
  defineUnitConversion(
    "L2/T",
    "m2/s",
    s/m^2) "Diffusivity";
  defineUnitConversion(
    "L3",
    "m3",
    1/m^3) "Volume";
  defineUnitConversion(
    "L3/N",
    "m3/C",
    C/m^3) "Specific volume";

  defineUnitConversion(
    "L3/(N.s)",
    "m3/(C.s)",
    C*s/m^3) "for derivative of specific volume in Dymola";
  defineUnitConversion(
    "L3/(N.T)",
    "m3/(C.s)",
    C*s/m^3) "Rate of specific volume";
  defineUnitConversion(
    "L3/N",
    "m3/mol",
    mol/m^3) "Specific volume";
  defineUnitConversion(
    "L3/(N.s)",
    "m3/(mol.s)",
    mol*s/m^3) "for derivative of specific volume in Dymola";
  defineUnitConversion(
    "L3/(N.T)",
    "m3/(mol.s)",
    mol*s/m^3) "Rate of specific volume";
  defineUnitConversion(
    "L3/s",
    "m3/s",
    1/m^3) "for derivative of volume in Dymola";
  defineUnitConversion(
    "L3/T",
    "m3/s",
    s/m^3) "Rate of volume";
  defineUnitConversion(
    "N2.T2/(L2.M)",
    "uF",
    1/(micro*F)) "Capacitance";
  defineUnitConversion(
    "M/N",
    "ug/C",
    C/(micro*g)) "Specific mass";
  defineUnitConversion(
    "L2.M/N2",
    "uH",
    1/(micro*H)) "Inductance";
  defineUnitConversion(
    "T",
    "us",
    1/(micro*s)) "Time";
  defineUnitConversion(
    "T",
    "ms",
    1/(milli*s)) "Time";
  defineUnitConversion(
    "T",
    "min",
    1/min) "Time";
  defineUnitConversion(
    "L",
    "mm",
    1/mm) "Length";
  defineUnitConversion(
    "L/T",
    "mm/s",
    s/mm) "Velocity";
  defineUnitConversion(
    "L2/T",
    "mm2/s",
    s/mm^2) "Diffusivity";
  defineUnitConversion(
    "N",
    "mol",
    1/mol) "Amount";
  defineUnitConversion(
    "L.N/T",
    "mol.m/s",
    s/(mol*m)) "VelocityAmount";
  defineUnitConversion(
    "L.M/T2",
    "N",
    1/N) "Force";
  defineUnitConversion(
    "L.M/(N.T)",
    "N/A",
    A/N) "Continuity";
  defineUnitConversion(
    "L2.M/(N2.T)",
    "ohm",
    1/ohm) "Electrical resistance";
  defineUnitConversion(
    "M/(L.T2)",
    "Pa",
    1/Pa) "Pressure";
  defineUnitConversion(
    "M/(L.T2.s)",
    "Pa/s",
    1/Pa) "for derivative of pressure in Dymola";
  defineUnitConversion(
    "M/(L.T3)",
    "Pa/s",
    s/Pa) "Rate of pressure";
  defineUnitConversion(
    "M/(T2.L.s)",
    "Pa/s",
    1/Pa) "for derivative of pressure in Dymola";
  defineUnitConversion(
    "N",
    "q",
    1/q) "Amount";
  defineUnitConversion(
    "A",
    "rad",
    1/rad) "Angle";
  defineUnitConversion(
    "A/L",
    "rad/m",
    m/rad) "Wavenumber";
  defineUnitConversion(
    "A/T",
    "rad/s",
    s/rad) "Frequency";
  defineUnitConversion(
    "N2.T/(L2.M)",
    "S",
    1/S) "Electrical conductance";
  defineUnitConversion(
    "T",
    "s",
    1/s) "Time";
  defineUnitConversion(
    "N2.T/(L3.M)",
    "S/m",
    m/S) "Electrical conductivity";
  defineUnitConversion(
    "T/L",
    "s/m",
    m/s) "Lineic time";
  defineUnitConversion(
    "A2",
    "sr",
    1/sr) "Solid angle";
  defineUnitConversion(
    "L2/T2",
    "Sv",
    1/Sv) "Squared velocity";
  defineUnitConversion(
    "M/(A.N.T)",
    "T",
    1/T) "Areic magnetic flux";
  defineUnitConversion(
    "L/T",
    "um/s",
    s/um) "Velocity";
  defineUnitConversion(
    "L2.M/(N.T2)",
    "V",
    1/V) "Potential";
  defineUnitConversion(
    "L.M/(N.T2)",
    "V/m",
    m/V) "Specific force";
  defineUnitConversion(
    "L2.M/(N.T2.s)",
    "V/s",
    1/V) "for derivative of potential in Dymola";
  defineUnitConversion(
    "L2.M/(N.T3)",
    "V/s",
    s/V) "Rate of potential";
  defineUnitConversion(
    "L3.M/(A.N.T2)",
    "V.m/rad",
    rad/(V*m)) "Potential per wavenumber";
  defineUnitConversion(
    "L2.M/T3",
    "W",
    1/W) "Power";
  defineUnitConversion(
    "N/s",
    "W/K",
    K/J) "for derivative of amount in Dymola";
  defineUnitConversion(
    "N/T",
    "W/K",
    K/W) "Current";
  defineUnitConversion(
    "M/T3",
    "W/m2",
    m^2/W) "Areic power";
  defineUnitConversion(
    "M.T5/L8",
    "W/(m2.K4)",
    m^2*K^4/W) "Areic power per 4th power of potential";
  defineUnitConversion(
    "L2.M/(A2.T3)",
    "W/sr",
    sr/W) "Radiant power";
  defineUnitConversion(
    "L4.M/T3",
    "W.m2",
    1/(W*m^2)) "Power times area";
  defineUnitConversion(
    "L2.M/(A.N.T)",
    "Wb",
    1/Wb) "Magnetic flux";
  // -------- end from FCSys/Resources/quantities.xls

  // ------------------------------------------------------------------------
  // Conversions with offsets
  // ------------------------------------------------------------------------

  defineUnitConversion(
    "L2.M/(N.T2)",
    "degC",
    1/K,
    -273.15) "Temperature in degC";
  defineUnitConversion(
    "L2.M/(N.T2)",
    "degF",
    9/(5*K),
    32 - (9/5)*273.15) "Temperature in degF";
  defineUnitConversion(
    "M/(L.T2)",
    "kPag",
    1/kPa,
    -atm/kPa) "Pressure in kPag";

  print("Done.");
end setup;

FCSys.Units.from_degC FCSys.Units.from_degC

Convert from temperature in degree Celsius to temperature as a quantity

Information

Extends from Modelica.SIunits.Icons.Conversion (Base icon for conversion functions).

Inputs

TypeNameDefaultDescription
RealT_degC Temperature in degree Celsius

Outputs

TypeNameDescription
TemperatureAbsoluteTThermodynamic temperature [L2.M/(N.T2)]

Modelica definition

function from_degC 
  "Convert from temperature in degree Celsius to temperature as a quantity"
  extends Modelica.SIunits.Icons.Conversion;

  input Real T_degC "Temperature in degree Celsius";
  output Q.TemperatureAbsolute T "Thermodynamic temperature";

algorithm 
  T := (T_degC + 273.15)*K;
end from_degC;

FCSys.Units.to_degC FCSys.Units.to_degC

Convert from temperature as a quantity to temperature in degree Celsius

Information

Extends from Modelica.SIunits.Icons.Conversion (Base icon for conversion functions).

Inputs

TypeNameDefaultDescription
TemperatureAbsoluteT Thermodynamic temperature [L2.M/(N.T2)]

Outputs

TypeNameDescription
RealT_degCTemperature in degree Celsius

Modelica definition

function to_degC 
  "Convert from temperature as a quantity to temperature in degree Celsius"
  extends Modelica.SIunits.Icons.Conversion;

  input Q.TemperatureAbsolute T "Thermodynamic temperature";
  output Real T_degC "Temperature in degree Celsius";

algorithm 
  T_degC := T/K - 273.15;
end to_degC;

FCSys.Units.from_kPag FCSys.Units.from_kPag

Convert from gauge pressure in kilopascals to absolute pressure as a quantity

Information

Extends from Modelica.SIunits.Icons.Conversion (Base icon for conversion functions).

Inputs

TypeNameDefaultDescription
Realp_kPag Gauge pressure in kilopascals

Outputs

TypeNameDescription
PressureAbsolutepAbsolute pressure [M/(L.T2)]

Modelica definition

function from_kPag 
  "Convert from gauge pressure in kilopascals to absolute pressure as a quantity"
  extends Modelica.SIunits.Icons.Conversion;

  input Real p_kPag "Gauge pressure in kilopascals";
  output Q.PressureAbsolute p "Absolute pressure";

algorithm 
  p := p_kPag*kPa + atm;
end from_kPag;

FCSys.Units.to_kPag FCSys.Units.to_kPag

Convert from absolute pressure as a quantity to gauge pressure in kilopascals

Information

Extends from Modelica.SIunits.Icons.Conversion (Base icon for conversion functions).

Inputs

TypeNameDefaultDescription
PressureAbsolutep Absolute pressure [M/(L.T2)]

Outputs

TypeNameDescription
Realp_kPagGauge pressure in kilopascals

Modelica definition

function to_kPag 
  "Convert from absolute pressure as a quantity to gauge pressure in kilopascals"
  extends Modelica.SIunits.Icons.Conversion;

  input Q.PressureAbsolute p "Absolute pressure";
  output Real p_kPag "Gauge pressure in kilopascals";

algorithm 
  p_kPag := (p - atm)/kPa;
end to_kPag;