15. Constants and unit conversion — MDAnalysis.units

The base units of MDAnalysis trajectories are the Å (ångström) for length and ps (pico second) for time. By default, all positions are in Å and all times are in ps, regardless of how the MD code stored trajectory data. By default, MDAnalysis converts automatically to the MDAnalysis units when reading trajectories and converts back when writing. This makes it possible to write scripts that can be agnostic of the specifics of how a particular MD code stores trajectory data. Other base units are listed in the table on Base units in MDAnalysis.

Base units in MDAnalysis
quantity unit SI units
length Å \(10^{-10}\) m
time ps \(10^{-12}\) s
energy kJ/mol \(1.66053892103219 \times 10^{-21}\) J
charge \(e\) \(1.602176565 \times 10^{-19}\) As
force kJ/(mol·Å) \(1.66053892103219 \times 10^{-11}\) J/m
speed Å/ps \(100\) m/s

15.1. Implementation notes

All conversions with convert() are carried out in a simple fashion: the conversion factor \(f_{b,b'}\) from the base unit \(b\) to another unit \(b'\) is precomputed and stored (see Data). The numerical value of a quantity in unit \(b\) is \(X/b\) (e.g. for \(X = 1.23\,\mathrm{ps}\) the numerical value is \(X/\mathrm{ps} = 1.23\)). [1]

The new numerical value \(X'/b'\) of the quantity (in units of \(b'\)) is then

\[X'/b' = f_{b,b'} X/b\]

The function get_conversion_factor() returns the appropriate factor \(f_{b,b'}\).

Conversion between different units is always carried out via the base unit as an intermediate step:

x is in u1: from u1 to b:  x'  = x  / factor[u1]
            from b  to u2: x'' = x' * factor[u2]
so f[u1,u2] = factor[u2]/factor[u1]

15.1.1. Conversions

Examples for how to calculate some of the conversion factors that are hard-coded in units (see Data).

density:

Base unit is \(\mathrm{Å}^{-3}\):

n/x = n/A**3 * densityUnit_factor[x]

Example for how to calculate the conversion factor \(f_{\mathrm{Å}^{-3},\mathrm{nm}^{-3}}\) from \(\mathrm{Å^-3}\) to \(\mathrm{nm}^{-3}\):

\[f_{\mathrm{Å}^{-3},\mathrm{nm}^{-3}} = \frac{1\,\mathrm{nm}^{-3}}{1\,\mathrm{Å}^{-3}} = \frac{(10\,\mathrm{Å})^{-3}}{1\,\mathrm{Å}^{-3}} = 10^{-3}\]
concentration:

Example for how to convert the conversion factor to Molar (mol/l):

factor = 1 A**-3 / (N_Avogadro * (10**-9 dm)**-3)

relative to a density rho0 in g/cm^3:

M(H2O) = 18 g/mol   Molar mass of water

factor = 1/(1e-24 * N_Avogadro / M(H2O))

from rho/rho0 = n/(N_A * M**-1) / rho0 where [n] = 1/Volume, [rho] = mass/Volume

Note

In the future me might move towards using the Quantities package or scipy.constants.

15.2. Functions

MDAnalysis.units.get_conversion_factor(unit_type, u1, u2)[source]

generate the conversion factor u1 -> u2 by using the base unit as an intermediate

f[u1 -> u2] = factor[u2]/factor[u1]

Conversion of X (in u1) to X’ (in u2):

X’ = conversion_factor * X
MDAnalysis.units.convert(x, u1, u2)[source]

Convert value x in unit u1 to new value in u2.

Returns:Converted value.
Return type:float
Raises:ValueError – The units are not known or if one attempts to convert between incompatible units.

15.3. Data

MDAnalysis.units.constants = {'N_Avogadro': 6.02214129e+23, 'calorie': 4.184, 'elementary_charge': 1.602176565e-19}

Values of physical constants are taken from CODATA 2010 at NIST. The thermochemical calorie is defined in the ISO 80000-5:2007 standard and is also listed in the NIST Guide to SI: Appendix B.8: Factors for Units.

New in version 0.9.0.

MDAnalysis.units.lengthUnit_factor = {'A': 1.0, 'Angstrom': 1.0, 'angstrom': 1.0, 'femtometer': 100000.0, 'fm': 100000.0, 'nanometer': 0.1, 'nm': 0.1, 'picometer': 100.0, 'pm': 100.0, 'Å': 1.0}

The basic unit of length in MDAnalysis is the Angstrom. Conversion factors between the base unit and other lengthUnits x are stored. Conversions follow L/x = L/Angstrom * lengthUnit_factor[x]. x can be nm/nanometer or fm.

MDAnalysis.units.N_Avogadro = 6.02214129e+23

Avogadro’s constant in mol**-1.

Deprecated since version 0.9.0: Use N_Avogadro in dict constants instead.

Changed in version 0.9.0: Use CODATA 2010 value, which differs from the previously used value 6.02214179e+23 mol**-1 by -5.00000000e+16 mol**-1.

MDAnalysis.units.water = {'MolarMass': 18.016, 'SPC': 0.985, 'TIP3P': 1.002, 'TIP4P': 1.001, 'exp': 0.997}
water density values at T=298K, P=1atm [Jorgensen1998]
model g cm**-3
SPC 0.985(1)
TIP3P 1.002(1)
TIP4P 1.001(1)
exp 0.997

and molar mass 18.016 g mol**-1.

MDAnalysis.units.densityUnit_factor = {'A^{-3}': 1.0, 'Angstrom^{-3}': 1.0, 'Molar': 1660.5389210321898, 'SPC': 30.37184690488927, 'TIP3P': 29.85655608913766, 'TIP4P': 29.886382818497438, 'nanometer^{-3}': 1000.0, 'nm^{-3}': 1000.0, 'water': 30.00628806551247, 'Å^{-3}': 1.0}

The basic unit for densities is Angstroem**(-3), i.e. the volume per molecule in A**3. Especially for water it can be convenient to measure the density relative to bulk, and hence a number of values are pre-stored in water.

MDAnalysis.units.timeUnit_factor = {'AKMA': 20.45482949774598, 'femtosecond': 1000.0, 'fs': 1000.0, 'nanosecond': 0.001, 'ns': 0.001, 'picosecond': 1.0, 'ps': 1.0, 's': 1e-12, 'sec': 1e-12, 'second': 1e-12}

For time, the basic unit is ps; in particular CHARMM’s 1 AKMA time unit = 4.888821E-14 sec is supported.

MDAnalysis.units.speedUnit_factor = {'A/ps': 1.0, 'Angstrom/AKMA': 0.04888821, 'Angstrom/femtosecond': 1000.0, 'Angstrom/fs': 1000.0, 'Angstrom/picosecond': 1.0, 'Angstrom/ps': 1.0, 'angstrom/femtosecond': 1000.0, 'angstrom/fs': 1000.0, 'angstrom/picosecond': 1.0, 'm/s': 100.0, 'nanometer/picosecond': 0.1, 'nanometer/ps': 0.1, 'nm/ns': 100.0, 'nm/ps': 0.1, 'pm/ps': 100.0, 'Å/ps': 1.0}

For speed, the basic unit is Angstrom/ps.

MDAnalysis.units.forceUnit_factor = {'J/m': 1.66053892103219e-11, 'N': 1.66053892103219e-11, 'Newton': 1.66053892103219e-11, 'kJ/(mol*A)': 1.0, 'kJ/(mol*Angstrom)': 1.0, 'kJ/(mol*nm)': 10.0, 'kJ/(mol*Å)': 1.0, 'kcal/(mol*Angstrom)': 0.2390057361376673}

For force the basic unit is kJ/(mol*Angstrom).

MDAnalysis.units.chargeUnit_factor = {'Amber': 18.2223, 'As': 1.602176565e-19, 'C': 1.602176565e-19, 'e': 1.0}

Charge is measured in multiples of the electron charge e, with the value elementary_charge in constants. The conversion factor to Amber charge units is 18.2223.

Changed in version 0.9.0: Use CODATA 2010 value for elementary charge, which differs from the previously used value e = 1.602176487 x 10**(-19) C by 7.8000000e-27 C.

MDAnalysis.units.conversion_factor = {'charge': {'Amber': 18.2223, 'As': 1.602176565e-19, 'C': 1.602176565e-19, 'e': 1.0}, 'density': {'A^{-3}': 1.0, 'Angstrom^{-3}': 1.0, 'Molar': 1660.5389210321898, 'SPC': 30.37184690488927, 'TIP3P': 29.85655608913766, 'TIP4P': 29.886382818497438, 'nanometer^{-3}': 1000.0, 'nm^{-3}': 1000.0, 'water': 30.00628806551247, 'Å^{-3}': 1.0}, 'energy': {'J': 1.66053892103219e-21, 'eV': 0.01036426919046959, 'kJ/mol': 1.0, 'kcal/mol': 0.2390057361376673}, 'force': {'J/m': 1.66053892103219e-11, 'N': 1.66053892103219e-11, 'Newton': 1.66053892103219e-11, 'kJ/(mol*A)': 1.0, 'kJ/(mol*Angstrom)': 1.0, 'kJ/(mol*nm)': 10.0, 'kJ/(mol*Å)': 1.0, 'kcal/(mol*Angstrom)': 0.2390057361376673}, 'length': {'A': 1.0, 'Angstrom': 1.0, 'angstrom': 1.0, 'femtometer': 100000.0, 'fm': 100000.0, 'nanometer': 0.1, 'nm': 0.1, 'picometer': 100.0, 'pm': 100.0, 'Å': 1.0}, 'speed': {'A/ps': 1.0, 'Angstrom/AKMA': 0.04888821, 'Angstrom/femtosecond': 1000.0, 'Angstrom/fs': 1000.0, 'Angstrom/picosecond': 1.0, 'Angstrom/ps': 1.0, 'angstrom/femtosecond': 1000.0, 'angstrom/fs': 1000.0, 'angstrom/picosecond': 1.0, 'm/s': 100.0, 'nanometer/picosecond': 0.1, 'nanometer/ps': 0.1, 'nm/ns': 100.0, 'nm/ps': 0.1, 'pm/ps': 100.0, 'Å/ps': 1.0}, 'time': {'AKMA': 20.45482949774598, 'femtosecond': 1000.0, 'fs': 1000.0, 'nanosecond': 0.001, 'ns': 0.001, 'picosecond': 1.0, 'ps': 1.0, 's': 1e-12, 'sec': 1e-12, 'second': 1e-12}}

any observable with a unit (i.e. one with an entry in the unit attribute) needs an entry in conversion_factor

Type:conversion_factor is used by get_conversion_factor()
Type:Note
MDAnalysis.units.unit_types = {'A': 'length', 'A/ps': 'speed', 'AKMA': 'time', 'A^{-3}': 'density', 'Amber': 'charge', 'Angstrom': 'length', 'Angstrom/AKMA': 'speed', 'Angstrom/femtosecond': 'speed', 'Angstrom/fs': 'speed', 'Angstrom/picosecond': 'speed', 'Angstrom/ps': 'speed', 'Angstrom^{-3}': 'density', 'As': 'charge', 'C': 'charge', 'J': 'energy', 'J/m': 'force', 'Molar': 'density', 'N': 'force', 'Newton': 'force', 'SPC': 'density', 'TIP3P': 'density', 'TIP4P': 'density', 'angstrom': 'length', 'angstrom/femtosecond': 'speed', 'angstrom/fs': 'speed', 'angstrom/picosecond': 'speed', 'e': 'charge', 'eV': 'energy', 'femtometer': 'length', 'femtosecond': 'time', 'fm': 'length', 'fs': 'time', 'kJ/(mol*A)': 'force', 'kJ/(mol*Angstrom)': 'force', 'kJ/(mol*nm)': 'force', 'kJ/(mol*Å)': 'force', 'kJ/mol': 'energy', 'kcal/(mol*Angstrom)': 'force', 'kcal/mol': 'energy', 'm/s': 'speed', 'nanometer': 'length', 'nanometer/picosecond': 'speed', 'nanometer/ps': 'speed', 'nanometer^{-3}': 'density', 'nanosecond': 'time', 'nm': 'length', 'nm/ns': 'speed', 'nm/ps': 'speed', 'nm^{-3}': 'density', 'ns': 'time', 'picometer': 'length', 'picosecond': 'time', 'pm': 'length', 'pm/ps': 'speed', 'ps': 'time', 's': 'time', 'sec': 'time', 'second': 'time', 'water': 'density', 'Å': 'length', 'Å/ps': 'speed', 'Å^{-3}': 'density'}

returns the type of unit for a known input unit. Note: Any unit must be unique because this dict is used to guess the unit type.

Type:Generated lookup table (dict)

15.4. References and footnotes

[Jorgensen1998]
  1. Jorgensen, C. Jenson, J Comp Chem 19 (1998), 1179-1186

Footnotes

[1]

One can also consider the conversion factor to carry units \(b'/b\), in which case the conversion formula would become

\[X' = f_{b,b'} X\]