Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / mdlib / tgroup.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, 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 /* This file is completely threadsafe - keep it that way! */
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43
44 #include <math.h>
45 #include "macros.h"
46 #include "main.h"
47 #include "smalloc.h"
48 #include "futil.h"
49 #include "tgroup.h"
50 #include "vec.h"
51 #include "network.h"
52 #include "smalloc.h"
53 #include "update.h"
54 #include "rbin.h"
55 #include "mtop_util.h"
56 #include "gmx_omp_nthreads.h"
57
58 static void init_grptcstat(int ngtc,t_grp_tcstat tcstat[])
59
60     int i,j;
61     
62     for(i=0; (i<ngtc); i++) {
63         tcstat[i].T = 0;
64         clear_mat(tcstat[i].ekinh);
65         clear_mat(tcstat[i].ekinh_old);
66         clear_mat(tcstat[i].ekinf);
67     }
68 }
69
70 static void init_grpstat(FILE *log,
71                          gmx_mtop_t *mtop,int ngacc,t_grp_acc gstat[])
72 {
73     gmx_groups_t *groups;
74     gmx_mtop_atomloop_all_t aloop;
75     int    i,grp;
76     t_atom *atom;
77     
78     if (ngacc > 0) {
79         groups = &mtop->groups;
80         aloop = gmx_mtop_atomloop_all_init(mtop);
81         while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) 
82         {
83             grp = ggrpnr(groups,egcACC,i);
84             if ((grp < 0) && (grp >= ngacc))
85             {
86                 gmx_incons("Input for acceleration groups wrong");
87             }
88             gstat[grp].nat++;
89             /* This will not work for integrator BD */
90             gstat[grp].mA += atom->m;
91             gstat[grp].mB += atom->mB;
92         }
93     }
94 }
95
96 void init_ekindata(FILE *log,gmx_mtop_t *mtop,t_grpopts *opts,
97                    gmx_ekindata_t *ekind)
98 {
99   int i;
100   int nthread,thread;
101 #ifdef DEBUG
102   fprintf(log,"ngtc: %d, ngacc: %d, ngener: %d\n",opts->ngtc,opts->ngacc,
103           opts->ngener);
104 #endif
105
106   /* bNEMD tells if we should remove remove the COM velocity
107    * from the velocities during velocity scaling in T-coupling.
108    * Turn this on when we have multiple acceleration groups
109    * or one accelerated group.
110    */
111   ekind->bNEMD = (opts->ngacc > 1 || norm(opts->acc[0]) > 0);
112
113   ekind->ngtc = opts->ngtc;
114   snew(ekind->tcstat,opts->ngtc);
115   init_grptcstat(opts->ngtc,ekind->tcstat);
116   /* Set Berendsen tcoupl lambda's to 1, 
117    * so runs without Berendsen coupling are not affected.
118    */
119   for(i=0; i<opts->ngtc; i++) 
120   {
121       ekind->tcstat[i].lambda = 1.0;
122       ekind->tcstat[i].vscale_nhc = 1.0;
123       ekind->tcstat[i].ekinscaleh_nhc = 1.0;
124       ekind->tcstat[i].ekinscalef_nhc = 1.0;
125   }
126   
127     nthread = gmx_omp_nthreads_get(emntUpdate);
128
129     snew(ekind->ekin_work_alloc,nthread);
130     snew(ekind->ekin_work,nthread);
131 #pragma omp parallel for num_threads(nthread) schedule(static)
132     for(thread=0; thread<nthread; thread++)
133     {
134         /* Allocate 2 elements extra on both sides,
135          * so in single precision we have 2*3*3*4=72 bytes buffer
136          * on both sides to avoid cache pollution.
137          */
138         snew(ekind->ekin_work_alloc[thread],ekind->ngtc+4);
139         ekind->ekin_work[thread] = ekind->ekin_work_alloc[thread] + 2;
140     }
141
142   ekind->ngacc = opts->ngacc;
143   snew(ekind->grpstat,opts->ngacc);
144   init_grpstat(log,mtop,opts->ngacc,ekind->grpstat);
145 }
146
147 void accumulate_u(t_commrec *cr,t_grpopts *opts,gmx_ekindata_t *ekind)
148 {
149     /* This routine will only be called when it's necessary */
150     t_bin *rb;
151     int   g;
152     
153     rb = mk_bin();
154     
155     for(g=0; (g<opts->ngacc); g++) 
156     {
157         add_binr(rb,DIM,ekind->grpstat[g].u);
158     }
159     sum_bin(rb,cr);
160     
161     for(g=0; (g<opts->ngacc); g++) 
162     {
163         extract_binr(rb,DIM*g,DIM,ekind->grpstat[g].u);
164     }
165     destroy_bin(rb);
166 }       
167
168 /* I don't think accumulate_ekin is used anymore? */
169
170 #if 0
171 static void accumulate_ekin(t_commrec *cr,t_grpopts *opts,
172                             gmx_ekindata_t *ekind)
173 {
174     int g;
175
176     if(PAR(cr))
177     {
178         for(g=0; (g<opts->ngtc); g++) 
179         {
180             gmx_sum(DIM*DIM,ekind->tcstat[g].ekinf[0],cr);
181         }
182     }
183 }       
184 #endif 
185
186 void update_ekindata(int start,int homenr,gmx_ekindata_t *ekind,
187                      t_grpopts *opts,rvec v[],t_mdatoms *md,real lambda)
188 {
189   int  d,g,n;
190   real mv;
191
192   /* calculate mean velocities at whole timestep */ 
193   for(g=0; (g<opts->ngtc); g++) {
194     ekind->tcstat[g].T = 0;
195   }
196
197   if (ekind->bNEMD) {
198     for (g=0; (g<opts->ngacc); g++)
199       clear_rvec(ekind->grpstat[g].u);
200     
201     g = 0;
202     for(n=start; (n<start+homenr); n++) {
203       if (md->cACC)
204         g = md->cACC[n];
205       for(d=0; (d<DIM);d++) {
206           mv = md->massT[n]*v[n][d];
207           ekind->grpstat[g].u[d] += mv;
208       }
209     }
210
211     for (g=0; (g < opts->ngacc); g++) {
212       for(d=0; (d<DIM);d++) {
213           ekind->grpstat[g].u[d] /=
214               (1-lambda)*ekind->grpstat[g].mA + lambda*ekind->grpstat[g].mB;
215       }
216     }
217   }
218 }
219
220 real sum_ekin(t_grpopts *opts,gmx_ekindata_t *ekind,real *dekindlambda,
221               gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld, gmx_bool bScaleEkin)
222 {
223     int          i,j,m,ngtc;
224     real         T,ek;
225     t_grp_tcstat *tcstat;
226     real         nrdf,nd,*ndf;
227     
228     ngtc = opts->ngtc;
229     ndf  = opts->nrdf;
230     
231     T = 0; 
232     nrdf = 0;
233
234     clear_mat(ekind->ekin);
235
236     for(i=0; (i<ngtc); i++) 
237     {
238
239         nd = ndf[i];
240         tcstat = &ekind->tcstat[i];
241         /* Sometimes a group does not have degrees of freedom, e.g.
242          * when it consists of shells and virtual sites, then we just
243          * set the temperatue to 0 and also neglect the kinetic
244          * energy, which should be  zero anyway.
245          */
246         
247         if (nd > 0) {
248             if (bEkinAveVel)
249             {
250                 if (!bScaleEkin)
251                 {
252                     /* in this case, kinetic energy is from the current velocities already */
253                     msmul(tcstat->ekinf,tcstat->ekinscalef_nhc,tcstat->ekinf);
254                 }
255             } 
256             else
257
258             {
259                 /* Calculate the full step Ekin as the average of the half steps */
260                 for(j=0; (j<DIM); j++)
261                 {
262                     for(m=0; (m<DIM); m++)
263                     {
264                         tcstat->ekinf[j][m] =
265                             0.5*(tcstat->ekinh[j][m]*tcstat->ekinscaleh_nhc + tcstat->ekinh_old[j][m]);
266                     }
267                 }
268             }
269             m_add(tcstat->ekinf,ekind->ekin,ekind->ekin);
270             
271             tcstat->Th = calc_temp(trace(tcstat->ekinh),nd);
272             tcstat->T  = calc_temp(trace(tcstat->ekinf),nd);
273
274             /* after the scaling factors have been multiplied in, we can remove them */
275             if (bEkinAveVel) 
276             {
277                 tcstat->ekinscalef_nhc = 1.0;
278             } 
279             else 
280             {
281                 tcstat->ekinscaleh_nhc = 1.0;
282             }
283         }
284         else 
285         {
286             tcstat->T  = 0;
287             tcstat->Th = 0;
288         }
289         T    += nd*tcstat->T;
290         nrdf += nd;
291     }
292     if (nrdf > 0)
293     {
294         T/=nrdf;
295     }
296     if (dekindlambda) 
297     {
298         *dekindlambda = 0.5*(ekind->dekindl + ekind->dekindl_old);
299     }
300     return T;
301 }