Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gen_maxwell_velocities.cpp
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,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.
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 "gen_maxwell_velocities.h"
40
41 #include <cmath>
42
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"
52
53 static void low_mspeed(real tempi, gmx_mtop_t* mtop, rvec v[], gmx::ThreeFry2x64<>* rng)
54 {
55     int                                    nrdf;
56     real                                   boltz;
57     real                                   ekin, temp;
58     gmx::TabulatedNormalDistribution<real> normalDist;
59
60     boltz = BOLTZ * tempi;
61     ekin  = 0.0;
62     nrdf  = 0;
63     for (const AtomProxy atomP : AtomRange(*mtop))
64     {
65         const t_atom& local = atomP.atom();
66         int           i     = atomP.globalAtomNumber();
67         real          mass  = local.m;
68         if (mass > 0)
69         {
70             rng->restart(i, 0);
71             real sd = std::sqrt(boltz / mass);
72             for (int m = 0; (m < DIM); m++)
73             {
74                 v[i][m] = sd * normalDist(*rng);
75                 ekin += 0.5 * mass * v[i][m] * v[i][m];
76             }
77             nrdf += DIM;
78         }
79     }
80     temp = (2.0 * ekin) / (nrdf * BOLTZ);
81     if (temp > 0)
82     {
83         real scal = std::sqrt(tempi / temp);
84         for (int i = 0; (i < mtop->natoms); i++)
85         {
86             for (int m = 0; (m < DIM); m++)
87             {
88                 v[i][m] *= scal;
89             }
90         }
91     }
92     fprintf(stderr, "Velocities were taken from a Maxwell distribution at %g K\n", tempi);
93     if (debug)
94     {
95         fprintf(debug,
96                 "Velocities were taken from a Maxwell distribution\n"
97                 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
98                 temp, tempi);
99     }
100 }
101
102 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t* mtop, rvec v[])
103 {
104
105     if (seed == 0)
106     {
107         seed = static_cast<int>(gmx::makeRandomSeed());
108         fprintf(stderr, "Using random seed %u for generating velocities\n", seed);
109     }
110     gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
111
112     low_mspeed(tempi, mtop, v, &rng);
113 }
114
115 static real calc_cm(int natoms, const real mass[], rvec x[], rvec v[], rvec xcm, rvec vcm, rvec acm, matrix L)
116 {
117     rvec dx, a0;
118     real tm, m0;
119     int  i, m;
120
121     clear_rvec(xcm);
122     clear_rvec(vcm);
123     clear_rvec(acm);
124     tm = 0.0;
125     for (i = 0; (i < natoms); i++)
126     {
127         m0 = mass[i];
128         tm += m0;
129         cprod(x[i], v[i], a0);
130         for (m = 0; (m < DIM); m++)
131         {
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. */
135         }
136     }
137     cprod(xcm, vcm, a0);
138     for (m = 0; (m < DIM); m++)
139     {
140         xcm[m] /= tm;
141         vcm[m] /= tm;
142         acm[m] -= a0[m] / tm;
143     }
144
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])
147 #ifdef DEBUG
148     PVEC("xcm", xcm);
149     PVEC("acm", acm);
150     PVEC("vcm", vcm);
151 #endif
152
153     clear_mat(L);
154     for (i = 0; (i < natoms); i++)
155     {
156         m0 = mass[i];
157         for (m = 0; (m < DIM); m++)
158         {
159             dx[m] = x[i][m] - xcm[m];
160         }
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;
167     }
168 #ifdef DEBUG
169     PVEC("L-x", L[XX]);
170     PVEC("L-y", L[YY]);
171     PVEC("L-z", L[ZZ]);
172 #endif
173
174     return tm;
175 }
176
177 void stop_cm(FILE gmx_unused* log, int natoms, real mass[], rvec x[], rvec v[])
178 {
179     rvec   xcm, vcm, acm;
180     tensor L;
181     int    i, m;
182
183 #ifdef DEBUG
184     fprintf(log, "stopping center of mass motion...\n");
185 #endif
186     (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
187
188     /* Subtract center of mass velocity */
189     for (i = 0; (i < natoms); i++)
190     {
191         for (m = 0; (m < DIM); m++)
192         {
193             v[i][m] -= vcm[m];
194         }
195     }
196
197 #ifdef DEBUG
198     (void)calc_cm(log, natoms, mass, x, v, xcm, vcm, acm, L);
199 #endif
200 }