Apply re-formatting to C++ in src/ 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 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "gen_maxwell_velocities.h"
41
42 #include <cmath>
43
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/random/tabulatednormaldistribution.h"
48 #include "gromacs/random/threefry.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/logger.h"
53 #include "gromacs/utility/smalloc.h"
54
55 static void low_mspeed(real tempi, gmx_mtop_t* mtop, rvec v[], gmx::ThreeFry2x64<>* rng, const gmx::MDLogger& logger)
56 {
57     int                                    nrdf;
58     real                                   boltz;
59     real                                   ekin, temp;
60     gmx::TabulatedNormalDistribution<real> normalDist;
61
62     boltz = BOLTZ * tempi;
63     ekin  = 0.0;
64     nrdf  = 0;
65     for (const AtomProxy atomP : AtomRange(*mtop))
66     {
67         const t_atom& local = atomP.atom();
68         int           i     = atomP.globalAtomNumber();
69         real          mass  = local.m;
70         if (mass > 0)
71         {
72             rng->restart(i, 0);
73             real sd = std::sqrt(boltz / mass);
74             for (int m = 0; (m < DIM); m++)
75             {
76                 v[i][m] = sd * normalDist(*rng);
77                 ekin += 0.5 * mass * v[i][m] * v[i][m];
78             }
79             nrdf += DIM;
80         }
81     }
82     temp = (2.0 * ekin) / (nrdf * BOLTZ);
83     if (temp > 0)
84     {
85         real scal = std::sqrt(tempi / temp);
86         for (int i = 0; (i < mtop->natoms); i++)
87         {
88             for (int m = 0; (m < DIM); m++)
89             {
90                 v[i][m] *= scal;
91             }
92         }
93     }
94     GMX_LOG(logger.info)
95             .asParagraph()
96             .appendTextFormatted("Velocities were taken from a Maxwell distribution at %g K", tempi);
97     if (debug)
98     {
99         fprintf(debug,
100                 "Velocities were taken from a Maxwell distribution\n"
101                 "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
102                 temp,
103                 tempi);
104     }
105 }
106
107 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t* mtop, rvec v[], const gmx::MDLogger& logger)
108 {
109
110     if (seed == 0)
111     {
112         seed = static_cast<int>(gmx::makeRandomSeed());
113         GMX_LOG(logger.info)
114                 .asParagraph()
115                 .appendTextFormatted("Using random seed %u for generating velocities", seed);
116     }
117     gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
118
119     low_mspeed(tempi, mtop, v, &rng, logger);
120 }
121
122 static real calc_cm(int natoms, const real mass[], rvec x[], rvec v[], rvec xcm, rvec vcm, rvec acm, matrix L)
123 {
124     rvec dx, a0;
125     real tm, m0;
126     int  i, m;
127
128     clear_rvec(xcm);
129     clear_rvec(vcm);
130     clear_rvec(acm);
131     tm = 0.0;
132     for (i = 0; (i < natoms); i++)
133     {
134         m0 = mass[i];
135         tm += m0;
136         cprod(x[i], v[i], a0);
137         for (m = 0; (m < DIM); m++)
138         {
139             xcm[m] += m0 * x[i][m]; /* c.o.m. position */
140             vcm[m] += m0 * v[i][m]; /* c.o.m. velocity */
141             acm[m] += m0 * a0[m];   /* rotational velocity around c.o.m. */
142         }
143     }
144     cprod(xcm, vcm, a0);
145     for (m = 0; (m < DIM); m++)
146     {
147         xcm[m] /= tm;
148         vcm[m] /= tm;
149         acm[m] -= a0[m] / tm;
150     }
151
152 #define PVEC(str, v) \
153     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 #ifdef DEBUG
155     PVEC("xcm", xcm);
156     PVEC("acm", acm);
157     PVEC("vcm", vcm);
158 #endif
159
160     clear_mat(L);
161     for (i = 0; (i < natoms); i++)
162     {
163         m0 = mass[i];
164         for (m = 0; (m < DIM); m++)
165         {
166             dx[m] = x[i][m] - xcm[m];
167         }
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;
174     }
175 #ifdef DEBUG
176     PVEC("L-x", L[XX]);
177     PVEC("L-y", L[YY]);
178     PVEC("L-z", L[ZZ]);
179 #endif
180
181     return tm;
182 }
183
184 void stop_cm(const gmx::MDLogger gmx_unused& logger, int natoms, real mass[], rvec x[], rvec v[])
185 {
186     rvec   xcm, vcm, acm;
187     tensor L;
188     int    i, m;
189
190 #ifdef DEBUG
191     GMX_LOG(logger.info).asParagraph().appendTextFormatted("stopping center of mass motion...");
192 #endif
193     (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
194
195     /* Subtract center of mass velocity */
196     for (i = 0; (i < natoms); i++)
197     {
198         for (m = 0; (m < DIM); m++)
199         {
200             v[i][m] -= vcm[m];
201         }
202     }
203 }