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,2015,2016,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.
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.
39 #include "gen_maxwell_velocities.h"
43 #include "gromacs/math/units.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/random/tabulatednormaldistribution.h"
47 #include "gromacs/random/threefry.h"
48 #include "gromacs/topology/mtop_util.h"
49 #include "gromacs/topology/topology.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 static void low_mspeed(real tempi, gmx_mtop_t* mtop, rvec v[], gmx::ThreeFry2x64<>* rng)
58 gmx::TabulatedNormalDistribution<real> normalDist;
60 boltz = BOLTZ * tempi;
63 for (const AtomProxy atomP : AtomRange(*mtop))
65 const t_atom& local = atomP.atom();
66 int i = atomP.globalAtomNumber();
71 real sd = std::sqrt(boltz / mass);
72 for (int m = 0; (m < DIM); m++)
74 v[i][m] = sd * normalDist(*rng);
75 ekin += 0.5 * mass * v[i][m] * v[i][m];
80 temp = (2.0 * ekin) / (nrdf * BOLTZ);
83 real scal = std::sqrt(tempi / temp);
84 for (int i = 0; (i < mtop->natoms); i++)
86 for (int m = 0; (m < DIM); m++)
92 fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n", tempi);
96 "Velocities were taken from a Maxwell distribution\n"
97 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
102 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t* mtop, rvec v[])
107 seed = static_cast<int>(gmx::makeRandomSeed());
108 fprintf(stderr, "Using random seed %u for generating velocities\n", seed);
110 gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
112 low_mspeed(tempi, mtop, v, &rng);
115 static real calc_cm(int natoms, const real mass[], rvec x[], rvec v[], rvec xcm, rvec vcm, rvec acm, matrix L)
125 for (i = 0; (i < natoms); i++)
129 cprod(x[i], v[i], a0);
130 for (m = 0; (m < DIM); m++)
132 xcm[m] += m0 * x[i][m]; /* c.o.m. position */
133 vcm[m] += m0 * v[i][m]; /* c.o.m. velocity */
134 acm[m] += m0 * a0[m]; /* rotational velocity around c.o.m. */
138 for (m = 0; (m < DIM); m++)
142 acm[m] -= a0[m] / tm;
145 #define PVEC(str, v) \
146 fprintf(log, "%s[X]: %10.5e %s[Y]: %10.5e %s[Z]: %10.5e\n", str, (v)[0], str, (v)[1], str, (v)[2])
154 for (i = 0; (i < natoms); i++)
157 for (m = 0; (m < DIM); m++)
159 dx[m] = x[i][m] - xcm[m];
161 L[XX][XX] += dx[XX] * dx[XX] * m0;
162 L[XX][YY] += dx[XX] * dx[YY] * m0;
163 L[XX][ZZ] += dx[XX] * dx[ZZ] * m0;
164 L[YY][YY] += dx[YY] * dx[YY] * m0;
165 L[YY][ZZ] += dx[YY] * dx[ZZ] * m0;
166 L[ZZ][ZZ] += dx[ZZ] * dx[ZZ] * m0;
177 void stop_cm(FILE gmx_unused* log, int natoms, real mass[], rvec x[], rvec v[])
184 fprintf(log, "stopping center of mass motion...\n");
186 (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
188 /* Subtract center of mass velocity */
189 for (i = 0; (i < natoms); i++)
191 for (m = 0; (m < DIM); m++)
198 (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);