Rework -Weverything
[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,2021, 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                                   ekin, temp;
59     gmx::TabulatedNormalDistribution<real> normalDist;
60
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(gmx::c_boltz * tempi / 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 * gmx::c_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     GMX_LOG(logger.info)
93             .asParagraph()
94             .appendTextFormatted("Velocities were taken from a Maxwell distribution at %g K", tempi);
95     if (debug)
96     {
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,
101                 tempi);
102     }
103 }
104
105 void maxwell_speed(real tempi, unsigned int seed, gmx_mtop_t* mtop, rvec v[], const gmx::MDLogger& logger)
106 {
107
108     if (seed == 0)
109     {
110         seed = static_cast<int>(gmx::makeRandomSeed());
111         GMX_LOG(logger.info)
112                 .asParagraph()
113                 .appendTextFormatted("Using random seed %u for generating velocities", seed);
114     }
115     gmx::ThreeFry2x64<> rng(seed, gmx::RandomDomain::MaxwellVelocities);
116
117     low_mspeed(tempi, mtop, v, &rng, logger);
118 }
119
120 static real calc_cm(int natoms, const real mass[], rvec x[], rvec v[], rvec xcm, rvec vcm, rvec acm, matrix L)
121 {
122     rvec dx, a0;
123     real tm, m0;
124     int  i, m;
125
126     clear_rvec(xcm);
127     clear_rvec(vcm);
128     clear_rvec(acm);
129     tm = 0.0;
130     for (i = 0; (i < natoms); i++)
131     {
132         m0 = mass[i];
133         tm += m0;
134         cprod(x[i], v[i], a0);
135         for (m = 0; (m < DIM); m++)
136         {
137             xcm[m] += m0 * x[i][m]; /* c.o.m. position */
138             vcm[m] += m0 * v[i][m]; /* c.o.m. velocity */
139             acm[m] += m0 * a0[m];   /* rotational velocity around c.o.m. */
140         }
141     }
142     cprod(xcm, vcm, a0);
143     for (m = 0; (m < DIM); m++)
144     {
145         xcm[m] /= tm;
146         vcm[m] /= tm;
147         acm[m] -= a0[m] / tm;
148     }
149
150     clear_mat(L);
151     for (i = 0; (i < natoms); i++)
152     {
153         m0 = mass[i];
154         for (m = 0; (m < DIM); m++)
155         {
156             dx[m] = x[i][m] - xcm[m];
157         }
158         L[XX][XX] += dx[XX] * dx[XX] * m0;
159         L[XX][YY] += dx[XX] * dx[YY] * m0;
160         L[XX][ZZ] += dx[XX] * dx[ZZ] * m0;
161         L[YY][YY] += dx[YY] * dx[YY] * m0;
162         L[YY][ZZ] += dx[YY] * dx[ZZ] * m0;
163         L[ZZ][ZZ] += dx[ZZ] * dx[ZZ] * m0;
164     }
165
166     return tm;
167 }
168
169 void stop_cm(const gmx::MDLogger gmx_unused& logger, int natoms, real mass[], rvec x[], rvec v[])
170 {
171     rvec   xcm, vcm, acm;
172     tensor L;
173     int    i, m;
174
175 #ifdef DEBUG
176     GMX_LOG(logger.info).asParagraph().appendTextFormatted("stopping center of mass motion...");
177 #endif
178     (void)calc_cm(natoms, mass, x, v, xcm, vcm, acm, L);
179
180     /* Subtract center of mass velocity */
181     for (i = 0; (i < natoms); i++)
182     {
183         for (m = 0; (m < DIM); m++)
184         {
185             v[i][m] -= vcm[m];
186         }
187     }
188 }