Moved random number generator code to separate dir.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gen_maxwell_velocities.c
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) 2013,2014, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "gromacs/random/random.h"
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "physics.h"
46 #include "typedefs.h"
47 #include "vec.h"
48 #include "gen_maxwell_velocities.h"
49 #include "mtop_util.h"
50
51 static void low_mspeed(real tempi,
52                        gmx_mtop_t *mtop, rvec v[], gmx_rng_t rng)
53 {
54     int                     i, m, nrdf;
55     real                    boltz, sd;
56     real                    ekin, temp, mass, scal;
57     gmx_mtop_atomloop_all_t aloop;
58     t_atom                 *atom;
59
60     boltz = BOLTZ*tempi;
61     ekin  = 0.0;
62     nrdf  = 0;
63     aloop = gmx_mtop_atomloop_all_init(mtop);
64     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
65     {
66         mass = atom->m;
67         if (mass > 0)
68         {
69             sd = sqrt(boltz/mass);
70             for (m = 0; (m < DIM); m++)
71             {
72                 v[i][m] = sd*gmx_rng_gaussian_real(rng);
73                 ekin   += 0.5*mass*v[i][m]*v[i][m];
74             }
75             nrdf += DIM;
76         }
77     }
78     temp = (2.0*ekin)/(nrdf*BOLTZ);
79     if (temp > 0)
80     {
81         scal = sqrt(tempi/temp);
82         for (i = 0; (i < mtop->natoms); i++)
83         {
84             for (m = 0; (m < DIM); m++)
85             {
86                 v[i][m] *= scal;
87             }
88         }
89     }
90     fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n",
91             tempi);
92     if (debug)
93     {
94         fprintf(debug,
95                 "Velocities were taken from a Maxwell distribution\n"
96                 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
97                 temp, tempi);
98     }
99 }
100
101 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t *mtop, rvec v[])
102 {
103     atom_id  *dummy;
104     int       i;
105     gmx_rng_t rng;
106
107     if (seed == 0)
108     {
109         seed = gmx_rng_make_seed();
110         fprintf(stderr, "Using random seed %d for generating velocities\n", seed);
111     }
112
113     rng = gmx_rng_init(seed);
114
115     low_mspeed(tempi, mtop, v, rng);
116
117     gmx_rng_destroy(rng);
118 }
119
120 static real calc_cm(int natoms, real mass[], rvec x[], rvec v[],
121                     rvec xcm, rvec vcm, rvec acm, matrix L)
122 {
123     rvec dx, a0;
124     real tm, m0;
125     int  i, m;
126
127     clear_rvec(xcm);
128     clear_rvec(vcm);
129     clear_rvec(acm);
130     tm = 0.0;
131     for (i = 0; (i < natoms); i++)
132     {
133         m0  = mass[i];
134         tm += m0;
135         cprod(x[i], v[i], a0);
136         for (m = 0; (m < DIM); m++)
137         {
138             xcm[m] += m0*x[i][m]; /* c.o.m. position */
139             vcm[m] += m0*v[i][m]; /* c.o.m. velocity */
140             acm[m] += m0*a0[m];   /* rotational velocity around c.o.m. */
141         }
142     }
143     cprod(xcm, vcm, a0);
144     for (m = 0; (m < DIM); m++)
145     {
146         xcm[m] /= tm;
147         vcm[m] /= tm;
148         acm[m] -= a0[m]/tm;
149     }
150
151 #define PVEC(str, v) fprintf(log, \
152                              "%s[X]: %10.5e  %s[Y]: %10.5e  %s[Z]: %10.5e\n", \
153                              str, v[0], str, v[1], str, v[2])
154 #ifdef DEBUG
155     PVEC("xcm", xcm);
156     PVEC("acm", acm);
157     PVEC("vcm", vcm);
158 #endif
159
160     clear_mat(L);
161     for (i = 0; (i < natoms); i++)
162     {
163         m0 = mass[i];
164         for (m = 0; (m < DIM); m++)
165         {
166             dx[m] = x[i][m]-xcm[m];
167         }
168         L[XX][XX] += dx[XX]*dx[XX]*m0;
169         L[XX][YY] += dx[XX]*dx[YY]*m0;
170         L[XX][ZZ] += dx[XX]*dx[ZZ]*m0;
171         L[YY][YY] += dx[YY]*dx[YY]*m0;
172         L[YY][ZZ] += dx[YY]*dx[ZZ]*m0;
173         L[ZZ][ZZ] += dx[ZZ]*dx[ZZ]*m0;
174     }
175 #ifdef DEBUG
176     PVEC("L-x", L[XX]);
177     PVEC("L-y", L[YY]);
178     PVEC("L-z", L[ZZ]);
179 #endif
180
181     return tm;
182 }
183
184 void stop_cm(FILE gmx_unused *log, int natoms, real mass[], rvec x[], rvec v[])
185 {
186     rvec   xcm, vcm, acm;
187     tensor L;
188     int    i, m;
189
190 #ifdef DEBUG
191     fprintf(log, "stopping center of mass motion...\n");
192 #endif
193     (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
194
195     /* Subtract center of mass velocity */
196     for (i = 0; (i < natoms); i++)
197     {
198         for (m = 0; (m < DIM); m++)
199         {
200             v[i][m] -= vcm[m];
201         }
202     }
203
204 #ifdef DEBUG
205     (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);
206 #endif
207 }