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