Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / gmxlib / random.c
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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "physics.h"
46 #include "typedefs.h"
47 #include "vec.h"
48 #include "gmx_random.h"
49 #include "random.h"
50 #include "mtop_util.h"
51
52 static void low_mspeed(real tempi,
53                        gmx_mtop_t *mtop,rvec v[], gmx_rng_t rng)
54 {
55   int  i,m,nrdf;
56   real boltz,sd;
57   real ekin,temp,mass,scal;
58   gmx_mtop_atomloop_all_t aloop;
59   t_atom *atom;
60
61   boltz=BOLTZ*tempi;
62   ekin=0.0;
63   nrdf=0;
64   aloop = gmx_mtop_atomloop_all_init(mtop);
65   while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
66     mass = atom->m;
67     if (mass > 0) {
68       sd=sqrt(boltz/mass);
69       for (m=0; (m<DIM); m++) {
70         v[i][m]=sd*gmx_rng_gaussian_real(rng);
71         ekin += 0.5*mass*v[i][m]*v[i][m];
72       }
73       nrdf += DIM;
74     }
75   }
76   temp=(2.0*ekin)/(nrdf*BOLTZ);
77   if (temp > 0) {
78     scal=sqrt(tempi/temp);
79     for(i=0; (i<mtop->natoms); i++)
80       for(m=0; (m<DIM); m++)
81         v[i][m]*=scal;
82   }
83   fprintf(stderr,"Velocities were taken from a Maxwell distribution at %g K\n",
84           tempi);
85   if (debug) {
86     fprintf(debug,
87             "Velocities were taken from a Maxwell distribution\n"
88             "Initial generated temperature: %12.5e (scaled to: %12.5e)\n",
89             temp,tempi);
90   }
91 }
92
93 void maxwell_speed(real tempi,int seed,gmx_mtop_t *mtop, rvec v[])
94 {
95   atom_id *dummy;
96   int     i;
97   gmx_rng_t rng;
98   
99   if (seed == -1) {
100     seed = make_seed();
101     fprintf(stderr,"Using random seed %d for generating velocities\n",seed);
102   }
103   
104   rng = gmx_rng_init(seed);
105
106   low_mspeed(tempi,mtop,v,rng);
107
108   gmx_rng_destroy(rng);
109 }
110
111 real calc_cm(FILE *log,int natoms,real mass[],rvec x[],rvec v[],
112              rvec xcm,rvec vcm,rvec acm,matrix L)
113 {
114   rvec dx,a0;
115   real tm,m0;
116   int  i,m;
117
118   clear_rvec(xcm);
119   clear_rvec(vcm);
120   clear_rvec(acm);
121   tm=0.0;
122   for(i=0; (i<natoms); i++) {
123     m0=mass[i];
124     tm+=m0;
125     cprod(x[i],v[i],a0);
126     for(m=0; (m<DIM); m++) {
127       xcm[m]+=m0*x[i][m]; /* c.o.m. position */
128       vcm[m]+=m0*v[i][m]; /* c.o.m. velocity */
129       acm[m]+=m0*a0[m];   /* rotational velocity around c.o.m. */
130     }
131   }
132   cprod(xcm,vcm,a0);
133   for(m=0; (m<DIM); m++) {
134     xcm[m]/=tm;
135     vcm[m]/=tm;
136     acm[m]-=a0[m]/tm;
137   }
138
139 #define PVEC(str,v) fprintf(log,\
140                             "%s[X]: %10.5e  %s[Y]: %10.5e  %s[Z]: %10.5e\n", \
141                            str,v[0],str,v[1],str,v[2])
142 #ifdef DEBUG
143   PVEC("xcm",xcm);
144   PVEC("acm",acm);
145   PVEC("vcm",vcm);
146 #endif
147   
148   clear_mat(L);
149   for(i=0; (i<natoms); i++) {
150     m0=mass[i];
151     for(m=0; (m<DIM); m++)
152       dx[m]=x[i][m]-xcm[m];
153     L[XX][XX]+=dx[XX]*dx[XX]*m0;
154     L[XX][YY]+=dx[XX]*dx[YY]*m0;
155     L[XX][ZZ]+=dx[XX]*dx[ZZ]*m0;
156     L[YY][YY]+=dx[YY]*dx[YY]*m0;
157     L[YY][ZZ]+=dx[YY]*dx[ZZ]*m0;
158     L[ZZ][ZZ]+=dx[ZZ]*dx[ZZ]*m0;
159   }
160 #ifdef DEBUG
161   PVEC("L-x",L[XX]);
162   PVEC("L-y",L[YY]);
163   PVEC("L-z",L[ZZ]);
164 #endif
165
166   return tm;
167 }
168
169 void stop_cm(FILE *log,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   fprintf(log,"stopping center of mass motion...\n");
177 #endif
178   (void)calc_cm(log,natoms,mass,x,v,xcm,vcm,acm,L);
179   
180   /* Subtract center of mass velocity */
181   for(i=0; (i<natoms); i++)
182     for(m=0; (m<DIM); m++)
183       v[i][m]-=vcm[m];
184
185 #ifdef DEBUG
186   (void)calc_cm(log,natoms,mass,x,v,xcm,vcm,acm,L);
187 #endif
188 }