Structural Bioinformatics Library
Template C++ / Python API for developping structural bioinformatics applications.
User Manual

UnitSystemTraits

Authors: F. Cazals and T. Dreyfus

The UnitSystemTraits package provides two things, namely a concept to represent a units system, together with accompanying models. It is in particular meant to ensure the compile time coherence of the dimensions of the manipulated quantities.

Introduction

The International System of Units is a coherent system of units based on seven base units–represented on the package logo. This kind of system ensures the coherence of calculations relying on different data types.

For example, in biomolecular modeling, computing the potential energy of a protein uses input files describing atomic positions expressed in Angstroms, while force fields may express distances in nanometers.

To foster the coherence of such calculations, this package provides models to represent unit systems with easy inter-conversions based on the Boost Units Library

Using existing models in existing Applications : I/O

The main unit systems used in bio-physics are summarized in the following table (see also Parmed documentation for an overview of these unit systems) :

Unit system Length unit Mass unit Time unit Charge unit Temperature unit Amount unit Energy unit
SI meters kilograms second Ampere kelvin mole joule
CGS centimeters gram second Ampere kelvin mole 1e-7 joule
MD nanometers gram second q electron kelvin mole kilojoule
Planck Planck length Planck mass Planck time Planck charge Planck temperature item Planck energy
AKMA angstroms daltons AKMA time q electron kelvin mole kilocalorie

The SBL provides the two unit systems classically used in molecular dynamics : AKMA and MD.

When the data used do not comply with the units required for the calculations undertaken, a conversion is in order. For example, if the distance between two atoms expressed in MD is printed, the expression 1 nm is printed, even if the atomic coordinates are specified in Angstrom.

Using existing models to develop novel applications

We first introduce the Boost Units Library and discuss its usage in this package. Then we describe the different models provided in this package.

Pre-requisites : Boost units

The Boost Units Library offers a full mechanism to represent existing and custom unit systems in C++. There are five important notions :

  • a dimension is the representation of a physical dimension (e.g length, mass, time, etc.) and is the base for the dimensional analysis of the library;
  • a base unit is a combination of one or more dimensions;
  • a unit system is a free ensemble of base units; by free, we mean that any combination of base units cannot be used to represent another base unit;
  • a unit is a combination of one or more dimensions in a given unit system;
  • a quantity is a value in a specified ensemble (integers, real values, etc.) associated to a specified unit.

This framework allows defining the traits class for the unit system, which contains :

  • the types for the units used in the SBL (e.g Length_unit, Mass_unit, Time_unit, etc.),
  • the types for the quantities used in the SBL (e.g Length_quantity, Mass_quantity, Time_quantity, etc.).
  • the unit constants used in the SBL (e.g 1 angstrom, 1 nanometer, 1 second, etc.).
The previous enumeration does not mention base units. Those types are internal and defined in the file SBL/Models/Unit_system_traits_base_units.hpp


The type of a constant depends on the unit system : if a constant refers to a base unit of the system, then it is a unit type; otherwise, t is a quantity type. Constants of both types can be used indifferently since they can be freely converted / cast to the corresponding quantity type.


Main models

The primary goal of this package is to ensure compile time coherence of calculations involving units. To this end, the following two classes are offered:

  • SBL::Models::T_Unit_system_traits_MD for MD.

The simplest way to use these models is to replace number types by the quantities these models define. For example, a classical dimensionless calculation of the Euclidean squared distance between two 2D points is done as follows :

double x1 = 0;
double y1 = 0;
double x2 = 1;
double y2 = 0;
double d = (x2 - x1)*(x2 - x1) + (y2 - y1)*(y2 - y1);

If we know that the coordinates are in angstroms, we can replace the double type by the Length_quantity type :

UnitSystemTraits::Length_quantity x1 = 0*UnitSystemTraits::angstrom;
UnitSystemTraits::Length_quantity y1 = 0*UnitSystemTraits::angstrom;
UnitSystemTraits::Length_quantity x2 = 1*UnitSystemTraits::angstrom;
UnitSystemTraits::Length_quantity y2 = 0*UnitSystemTraits::angstrom;
UnitSystemTraits::Length_quantity d = (x2 - x1)*(x2 - x1) + (y2 - y1)*(y2 - y1);

Note that the type UnitSystemTraits can be replaced by any model of the concept, even if the Angstrom is not a base unit of the unit system : the unit constant angstrom is defined in all models of UnitSystemTraits

Note also that it is possible to use nameless units when values are independent from the unit system to be used. For exemple, a molar energy value that has to be initialized to 1 whatever the underlying system cannot explicitly use AKMA or MD units, and will use these system-independent units :

UnitSystemTraits::Molar_energy_quantity E = 1*UnitSystemTraits::energy_unit/UnitSystemTraits::amount_unit;

The value of E will be always 1, whatever the underlaying unit system.

Dimensionless models

There are two limitations of this approach when using a third party library :

  • if the library is not generic w.r.t the number type, quantities cannot be used as a replacement of the number type;
  • even if the library is generic w.r.t the number type, only basic arithmetic operations are possible with quantities (+, -, /, *) : hence, if the library is using special functions (sqrt, cos, sin, etc.), it will not work.

To overcome these difficulties, one can use the static method value() of the concept UnitSystemTraits that returns the dimensionless value of the quantity. This package also provides dimensionless versions of AKMA and MD :

The main difference is that units and quantities are represented by a simple number type (double for example). This means also that there is no compile time check for the coherence of the dimensions. For example, one can check the coherence of dimensions of examples within Molecular_potential_energy using SBL::Models::T_Unit_system_traits_AKMA . However, the energy minimization provided in the package Real_value_function_minimizer may use libraries which are not generic w.r.t the number type .

Other models

Finally, the class SBL::Models::T_Unit_system_traits_for_potential_energy< UnitSystemTraitsBase > specializes the parameter UnitSystemTraitsBase and defines supplemental units and quantities used in the package Molecular_potential_energy

Developing new SBL models for the UnitSystemTraits concept

A model of UnitSystemTraits has to define the following types :

Length_unit
Mass_unit
Time_unit
Temperature_unit
Angle_unit
Amount_unit
Charge_unit
Energy_unit
Molar_energy_unit

Length_quantity
Mass_quantity
Time_quantity
Temperature_quantity
Angle_quantity
Amount_quantity
Charge_quantity
Energy_quantity
Molar_energy_quantity

It also has to define the following static constants :

length_unit;
angstrom;
angstroms;
nanometer;
nanometers;

mass_unit;
dalton;
daltons;

time_unit;
T;
picosecond;
picoseconds;

temperature_unit;
kelvin;
kelvins;

angle_unit;
radian;
radians;
degree;
degrees;

amount_unit;
mole;
moles;

charge_unit;
q_electron;
coulomb;
coulombs;

energy_unit;
kilocalorie;
kilocalories;
kilojoule;
kilojoules;

Note that these constants can be either of a unit type either of a quantity type indifferenty.

Finally, a model of UnitSystemTraits has to define the static method value(const Quantity& q) that returns the dimensionless value of the input quantity.