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