2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
42 #include "gromacs/random/random.h"
44 #include "gromacs/utility/smalloc.h"
48 #include "gen_maxwell_velocities.h"
49 #include "mtop_util.h"
51 static void low_mspeed(real tempi,
52 gmx_mtop_t *mtop, rvec v[], gmx_rng_t rng)
56 real ekin, temp, mass, scal;
57 gmx_mtop_atomloop_all_t aloop;
63 aloop = gmx_mtop_atomloop_all_init(mtop);
64 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
69 sd = sqrt(boltz/mass);
70 for (m = 0; (m < DIM); m++)
72 v[i][m] = sd*gmx_rng_gaussian_real(rng);
73 ekin += 0.5*mass*v[i][m]*v[i][m];
78 temp = (2.0*ekin)/(nrdf*BOLTZ);
81 scal = sqrt(tempi/temp);
82 for (i = 0; (i < mtop->natoms); i++)
84 for (m = 0; (m < DIM); m++)
90 fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n",
95 "Velocities were taken from a Maxwell distribution\n"
96 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
101 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t *mtop, rvec v[])
109 seed = gmx_rng_make_seed();
110 fprintf(stderr, "Using random seed %d for generating velocities\n", seed);
113 rng = gmx_rng_init(seed);
115 low_mspeed(tempi, mtop, v, rng);
117 gmx_rng_destroy(rng);
120 static real calc_cm(int natoms, real mass[], rvec x[], rvec v[],
121 rvec xcm, rvec vcm, rvec acm, matrix L)
131 for (i = 0; (i < natoms); i++)
135 cprod(x[i], v[i], a0);
136 for (m = 0; (m < DIM); m++)
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. */
144 for (m = 0; (m < DIM); m++)
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])
161 for (i = 0; (i < natoms); i++)
164 for (m = 0; (m < DIM); m++)
166 dx[m] = x[i][m]-xcm[m];
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;
184 void stop_cm(FILE gmx_unused *log, int natoms, real mass[], rvec x[], rvec v[])
191 fprintf(log, "stopping center of mass motion...\n");
193 (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
195 /* Subtract center of mass velocity */
196 for (i = 0; (i < natoms); i++)
198 for (m = 0; (m < DIM); m++)
205 (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);