Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / gmxpreprocess / gen_maxwell_velocities.cpp
blobdee74291c0827e0214df50dbc41b4341ff3eb1a1
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,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.
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 nrdf;
57 real boltz;
58 real ekin, temp;
59 gmx::TabulatedNormalDistribution<real> normalDist;
61 boltz = BOLTZ*tempi;
62 ekin = 0.0;
63 nrdf = 0;
64 for (const AtomProxy atomP : AtomRange(*mtop))
66 const t_atom &local = atomP.atom();
67 int i = atomP.globalAtomNumber();
68 real mass = local.m;
69 if (mass > 0)
71 rng->restart(i, 0);
72 real sd = std::sqrt(boltz/mass);
73 for (int m = 0; (m < DIM); m++)
75 v[i][m] = sd*normalDist(*rng);
76 ekin += 0.5*mass*v[i][m]*v[i][m];
78 nrdf += DIM;
81 temp = (2.0*ekin)/(nrdf*BOLTZ);
82 if (temp > 0)
84 real scal = std::sqrt(tempi/temp);
85 for (int i = 0; (i < mtop->natoms); i++)
87 for (int m = 0; (m < DIM); m++)
89 v[i][m] *= scal;
93 fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n",
94 tempi);
95 if (debug)
97 fprintf(debug,
98 "Velocities were taken from a Maxwell distribution\n"
99 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
100 temp, tempi);
104 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t *mtop, rvec v[])
107 if (seed == 0)
109 seed = static_cast<int>(gmx::makeRandomSeed());
110 fprintf(stderr, "Using random seed %u for generating velocities\n", seed);
112 gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
114 low_mspeed(tempi, mtop, v, &rng);
117 static real calc_cm(int natoms, const real mass[], rvec x[], rvec v[],
118 rvec xcm, rvec vcm, rvec acm, matrix L)
120 rvec dx, a0;
121 real tm, m0;
122 int i, m;
124 clear_rvec(xcm);
125 clear_rvec(vcm);
126 clear_rvec(acm);
127 tm = 0.0;
128 for (i = 0; (i < natoms); i++)
130 m0 = mass[i];
131 tm += m0;
132 cprod(x[i], v[i], a0);
133 for (m = 0; (m < DIM); m++)
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. */
140 cprod(xcm, vcm, a0);
141 for (m = 0; (m < DIM); m++)
143 xcm[m] /= tm;
144 vcm[m] /= tm;
145 acm[m] -= a0[m]/tm;
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
157 clear_mat(L);
158 for (i = 0; (i < natoms); i++)
160 m0 = mass[i];
161 for (m = 0; (m < DIM); m++)
163 dx[m] = x[i][m]-xcm[m];
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;
172 #ifdef DEBUG
173 PVEC("L-x", L[XX]);
174 PVEC("L-y", L[YY]);
175 PVEC("L-z", L[ZZ]);
176 #endif
178 return tm;
181 void stop_cm(FILE gmx_unused *log, int natoms, real mass[], rvec x[], rvec v[])
183 rvec xcm, vcm, acm;
184 tensor L;
185 int i, m;
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);
192 /* Subtract center of mass velocity */
193 for (i = 0; (i < natoms); i++)
195 for (m = 0; (m < DIM); m++)
197 v[i][m] -= vcm[m];
201 #ifdef DEBUG
202 (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);
203 #endif