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