865fd31175994d1a426977af9d2afbca082ce88f
[alexxy/gromacs.git] / src / gromacs / math / units.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifndef GMX_MATH_UNITS_H
38 #define GMX_MATH_UNITS_H
39
40 /*
41  * Physical constants to be used in Gromacs.
42  * No constants (apart from 0, 1 or 2) should
43  * be anywhere else in the code.
44  */
45
46 #include "gromacs/math/utilities.h"
47
48 #define ANGSTROM (1e-10)           /* Old...    */
49 #define KILO (1e3)                 /* Thousand  */
50 #define NANO (1e-9)                /* A Number  */
51 #define PICO (1e-12)               /* A Number  */
52 #define A2NM (ANGSTROM / NANO)     /* NANO              */
53 #define NM2A (NANO / ANGSTROM)     /* 10.0              */
54 #define RAD2DEG (180.0 / M_PI)     /* Conversion        */
55 #define DEG2RAD (M_PI / 180.0)     /* id                */
56 #define CAL2JOULE (4.184)          /* Exact definition of the calorie */
57 #define E_CHARGE (1.602176634e-19) /* Exact definition, Coulomb NIST 2018 CODATA */
58
59 #define AMU (1.66053906660e-27)                     /* kg, NIST 2018 CODATA  */
60 #define BOLTZMANN (1.380649e-23)                    /* (J/K, Exact definition, NIST 2018 CODATA */
61 #define AVOGADRO (6.02214076e23)                    /* 1/mol, Exact definition, NIST 2018 CODATA */
62 #define RGAS (BOLTZMANN * AVOGADRO)                 /* (J/(mol K))  */
63 #define BOLTZ (RGAS / KILO)                         /* (kJ/(mol K)) */
64 #define FARADAY (E_CHARGE * AVOGADRO)               /* (C/mol)      */
65 #define ELECTRONVOLT (E_CHARGE * AVOGADRO / KILO)   /* (kJ/mol)   */
66 #define PLANCK1 (6.62607015e-34)                    /* J/Hz, Exact definition, NIST 2018 CODATA */
67 #define PLANCK (PLANCK1 * AVOGADRO / (PICO * KILO)) /* (kJ/mol) ps */
68
69 #define EPSILON0_SI (8.8541878128e-12) /* F/m,  NIST 2018 CODATA */
70 /* Epsilon in our MD units: (e^2 / Na (kJ nm)) == (e^2 mol/(kJ nm)) */
71 #define EPSILON0 ((EPSILON0_SI * NANO * KILO) / (E_CHARGE * E_CHARGE * AVOGADRO))
72
73 #define SPEED_OF_LIGHT (2.99792458e05)  /* units of nm/ps, Exact definition, NIST 2018 CODATA */
74 #define ATOMICMASS_keV (931494.10242)   /* Atomic mass in keV, NIST 2018 CODATA   */
75 #define ELECTRONMASS_keV (510.99895000) /* Electron mas in keV, NIST 2018 CODATA  */
76
77 #define RYDBERG (1.0973731568160e-02) /* nm^-1, NIST 2018 CODATA */
78
79 #define ONE_4PI_EPS0 (1.0 / (4.0 * M_PI * EPSILON0))
80 #define FACEL (10.0 * ONE_4PI_EPS0)
81
82 /* Pressure in MD units is:
83  * 1 bar = 1e5 Pa = 1e5 kg m^-1 s^-2 = 1e-28 kg nm^-1 ps^-2 = 1e-28 / AMU amu nm^1 ps ^2
84  */
85 #define BAR_MDUNITS (1e5 * NANO * PICO * PICO / AMU)
86 #define PRESFAC (1.0 / BAR_MDUNITS)
87
88 /* DEBYE2ENM should be (1e-21*PICO)/(SPEED_OF_LIGHT*E_CHARGE*NANO*NANO),
89  * but we need to factor out some of the exponents to avoid single-precision overflows.
90  */
91 #define DEBYE2ENM (1e-15 / (SPEED_OF_LIGHT * E_CHARGE))
92 #define ENM2DEBYE (1.0 / DEBYE2ENM)
93
94 /* to convert from a acceleration in (e V)/(amu nm) */
95 /* FIELDFAC is also Faraday's constant and E_CHARGE/(1e6 AMU) */
96 #define FIELDFAC (FARADAY / KILO)
97
98 /* to convert AU to MD units: */
99 #define HARTREE2KJ ((2.0 * RYDBERG * PLANCK * SPEED_OF_LIGHT) / AVOGADRO)
100 #define BOHR2NM (0.0529177210903) /* nm^-1, NIST 2018 CODATA */
101 #define HARTREE_BOHR2MD (HARTREE2KJ * AVOGADRO / BOHR2NM)
102
103
104 /* The four basic units */
105 #define unit_length "nm"
106 #define unit_time "ps"
107 #define unit_mass "u"
108 #define unit_energy "kJ/mol"
109
110 /* Temperature unit, T in this unit times BOLTZ give energy in unit_energy */
111 #define unit_temp_K "K"
112
113 /* Charge unit, electron charge, involves ONE_4PI_EPS0 */
114 #define unit_charge_e "e"
115
116 /* Pressure unit, pressure in basic units times PRESFAC gives this unit */
117 #define unit_pres_bar "bar"
118
119 /* Dipole unit, debye, conversion from the unit_charge_e involves ENM2DEBYE */
120 #define unit_dipole_D "D"
121
122 /* Derived units from basic units only */
123 #define unit_vel unit_length "/" unit_time
124 #define unit_volume unit_length "^3"
125 #define unit_invtime "1/" unit_time
126
127 /* Other derived units */
128 #define unit_surft_bar unit_pres_bar " " unit_length
129
130 /* SI units, conversion from basic units involves NANO, PICO and AMU */
131 #define unit_length_SI "m"
132 #define unit_time_SI "s"
133 #define unit_mass_SI "kg"
134
135 #define unit_density_SI unit_mass_SI "/" unit_length_SI "^3"
136 #define unit_invvisc_SI unit_length_SI " " unit_time_SI "/" unit_mass_SI
137
138 /* The routines below can be used for converting units from or to GROMACS
139    internal units. */
140 enum
141 {
142     eg2cAngstrom,
143     eg2cNm,
144     eg2cBohr,
145     eg2cKcal_Mole,
146     eg2cHartree,
147     eg2cHartree_e,
148     eg2cAngstrom3,
149     eg2cCoulomb,
150     eg2cDebye,
151     eg2cElectron,
152     eg2cBuckingham,
153     eg2cNR
154 };
155
156 /* Convert value x to GROMACS units. Energy -> Energy, Length -> Length etc.
157    The type of x is deduced from unit,
158    which should be taken from the enum above. */
159 extern double convert2gmx(double x, int unit);
160
161 /* Convert value x from GROMACS units to the desired one.
162    The type of return value is deduced from unit, see above */
163 extern double gmx2convert(double x, int unit);
164
165 /* Convert the string to one of the units supported. Returns -1 if not found. */
166 extern int string2unit(char* string);
167
168 /* Convert the unit to a string. Return NULL when unit is out of range. */
169 extern const char* unit2string(int unit);
170
171 #endif