Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / gen_maxwell_velocities.cpp
blob350e1ad95176a82c4f356cf4e4cd1e7e346be393
1 /*
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, 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.
37 #include "gmxpre.h"
39 #include "gen_maxwell_velocities.h"
41 #include <cmath>
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,
54 gmx_mtop_t *mtop, rvec v[], gmx::ThreeFry2x64<> * rng)
56 int i, m, nrdf;
57 real boltz, sd;
58 real ekin, temp, mass, scal;
59 gmx_mtop_atomloop_all_t aloop;
60 const t_atom *atom;
61 gmx::TabulatedNormalDistribution<real> normalDist;
63 boltz = BOLTZ*tempi;
64 ekin = 0.0;
65 nrdf = 0;
66 aloop = gmx_mtop_atomloop_all_init(mtop);
67 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
69 mass = atom->m;
70 if (mass > 0)
72 rng->restart(i, 0);
73 sd = std::sqrt(boltz/mass);
74 for (m = 0; (m < DIM); m++)
76 v[i][m] = sd*normalDist(*rng);
77 ekin += 0.5*mass*v[i][m]*v[i][m];
79 nrdf += DIM;
82 temp = (2.0*ekin)/(nrdf*BOLTZ);
83 if (temp > 0)
85 scal = std::sqrt(tempi/temp);
86 for (i = 0; (i < mtop->natoms); i++)
88 for (m = 0; (m < DIM); m++)
90 v[i][m] *= scal;
94 fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n",
95 tempi);
96 if (debug)
98 fprintf(debug,
99 "Velocities were taken from a Maxwell distribution\n"
100 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
101 temp, tempi);
105 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t *mtop, rvec v[])
108 if (seed == 0)
110 seed = static_cast<int>(gmx::makeRandomSeed());
111 fprintf(stderr, "Using random seed %u for generating velocities\n", seed);
113 gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
115 low_mspeed(tempi, mtop, v, &rng);
118 static real calc_cm(int natoms, real mass[], rvec x[], rvec v[],
119 rvec xcm, rvec vcm, rvec acm, matrix L)
121 rvec dx, a0;
122 real tm, m0;
123 int i, m;
125 clear_rvec(xcm);
126 clear_rvec(vcm);
127 clear_rvec(acm);
128 tm = 0.0;
129 for (i = 0; (i < natoms); i++)
131 m0 = mass[i];
132 tm += m0;
133 cprod(x[i], v[i], a0);
134 for (m = 0; (m < DIM); m++)
136 xcm[m] += m0*x[i][m]; /* c.o.m. position */
137 vcm[m] += m0*v[i][m]; /* c.o.m. velocity */
138 acm[m] += m0*a0[m]; /* rotational velocity around c.o.m. */
141 cprod(xcm, vcm, a0);
142 for (m = 0; (m < DIM); m++)
144 xcm[m] /= tm;
145 vcm[m] /= tm;
146 acm[m] -= a0[m]/tm;
149 #define PVEC(str, v) fprintf(log, \
150 "%s[X]: %10.5e %s[Y]: %10.5e %s[Z]: %10.5e\n", \
151 str, v[0], str, v[1], str, v[2])
152 #ifdef DEBUG
153 PVEC("xcm", xcm);
154 PVEC("acm", acm);
155 PVEC("vcm", vcm);
156 #endif
158 clear_mat(L);
159 for (i = 0; (i < natoms); i++)
161 m0 = mass[i];
162 for (m = 0; (m < DIM); m++)
164 dx[m] = x[i][m]-xcm[m];
166 L[XX][XX] += dx[XX]*dx[XX]*m0;
167 L[XX][YY] += dx[XX]*dx[YY]*m0;
168 L[XX][ZZ] += dx[XX]*dx[ZZ]*m0;
169 L[YY][YY] += dx[YY]*dx[YY]*m0;
170 L[YY][ZZ] += dx[YY]*dx[ZZ]*m0;
171 L[ZZ][ZZ] += dx[ZZ]*dx[ZZ]*m0;
173 #ifdef DEBUG
174 PVEC("L-x", L[XX]);
175 PVEC("L-y", L[YY]);
176 PVEC("L-z", L[ZZ]);
177 #endif
179 return tm;
182 void stop_cm(FILE gmx_unused *log, int natoms, real mass[], rvec x[], rvec v[])
184 rvec xcm, vcm, acm;
185 tensor L;
186 int i, m;
188 #ifdef DEBUG
189 fprintf(log, "stopping center of mass motion...\n");
190 #endif
191 (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
193 /* Subtract center of mass velocity */
194 for (i = 0; (i < natoms); i++)
196 for (m = 0; (m < DIM); m++)
198 v[i][m] -= vcm[m];
202 #ifdef DEBUG
203 (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);
204 #endif