Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / 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  * Copyright (c) 2013, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42
43 #include <math.h>
44 #include "macros.h"
45 #include "main.h"
46 #include "smalloc.h"
47 #include "gromacs/fileio/futil.h"
48 #include "tgroup.h"
49 #include "vec.h"
50 #include "network.h"
51 #include "smalloc.h"
52 #include "update.h"
53 #include "rbin.h"
54 #include "mtop_util.h"
55 #include "gmx_omp_nthreads.h"
56
57 static void init_grptcstat(int ngtc, t_grp_tcstat tcstat[])
58 {
59     int i, j;
60
61     for (i = 0; (i < ngtc); i++)
62     {
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(gmx_mtop_t *mtop, int ngacc, t_grp_acc gstat[])
71 {
72     gmx_groups_t           *groups;
73     gmx_mtop_atomloop_all_t aloop;
74     int                     i, grp;
75     t_atom                 *atom;
76
77     if (ngacc > 0)
78     {
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 gmx_unused *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     snew(ekind->dekindl_work, nthread);
132 #pragma omp parallel for num_threads(nthread) schedule(static)
133     for (thread = 0; thread < nthread; thread++)
134     {
135 #define EKIN_WORK_BUFFER_SIZE 2
136         /* Allocate 2 extra elements on both sides, so in single
137          * precision we have
138          * EKIN_WORK_BUFFER_SIZE*DIM*DIM*sizeof(real) = 72/144 bytes
139          * buffer on both sides to avoid cache pollution.
140          */
141         snew(ekind->ekin_work_alloc[thread], ekind->ngtc+2*EKIN_WORK_BUFFER_SIZE);
142         ekind->ekin_work[thread] = ekind->ekin_work_alloc[thread] + EKIN_WORK_BUFFER_SIZE;
143         /* Nasty hack so we can have the per-thread accumulation
144          * variable for dekindl in the same thread-local cache lines
145          * as the per-thread accumulation tensors for ekin[fh],
146          * because they are accumulated in the same loop. */
147         ekind->dekindl_work[thread] = &(ekind->ekin_work[thread][ekind->ngtc][0][0]);
148 #undef EKIN_WORK_BUFFER_SIZE
149     }
150
151     ekind->ngacc = opts->ngacc;
152     snew(ekind->grpstat, opts->ngacc);
153     init_grpstat(mtop, opts->ngacc, ekind->grpstat);
154 }
155
156 void accumulate_u(t_commrec *cr, t_grpopts *opts, gmx_ekindata_t *ekind)
157 {
158     /* This routine will only be called when it's necessary */
159     t_bin *rb;
160     int    g;
161
162     rb = mk_bin();
163
164     for (g = 0; (g < opts->ngacc); g++)
165     {
166         add_binr(rb, DIM, ekind->grpstat[g].u);
167     }
168     sum_bin(rb, cr);
169
170     for (g = 0; (g < opts->ngacc); g++)
171     {
172         extract_binr(rb, DIM*g, DIM, ekind->grpstat[g].u);
173     }
174     destroy_bin(rb);
175 }
176
177 /* I don't think accumulate_ekin is used anymore? */
178
179 #if 0
180 static void accumulate_ekin(t_commrec *cr, t_grpopts *opts,
181                             gmx_ekindata_t *ekind)
182 {
183     int g;
184
185     if (PAR(cr))
186     {
187         for (g = 0; (g < opts->ngtc); g++)
188         {
189             gmx_sum(DIM*DIM, ekind->tcstat[g].ekinf[0], cr);
190         }
191     }
192 }
193 #endif
194
195 void update_ekindata(int start, int homenr, gmx_ekindata_t *ekind,
196                      t_grpopts *opts, rvec v[], t_mdatoms *md, real lambda)
197 {
198     int  d, g, n;
199     real mv;
200
201     /* calculate mean velocities at whole timestep */
202     for (g = 0; (g < opts->ngtc); g++)
203     {
204         ekind->tcstat[g].T = 0;
205     }
206
207     if (ekind->bNEMD)
208     {
209         for (g = 0; (g < opts->ngacc); g++)
210         {
211             clear_rvec(ekind->grpstat[g].u);
212         }
213
214         g = 0;
215         for (n = start; (n < start+homenr); n++)
216         {
217             if (md->cACC)
218             {
219                 g = md->cACC[n];
220             }
221             for (d = 0; (d < DIM); d++)
222             {
223                 mv                      = md->massT[n]*v[n][d];
224                 ekind->grpstat[g].u[d] += mv;
225             }
226         }
227
228         for (g = 0; (g < opts->ngacc); g++)
229         {
230             for (d = 0; (d < DIM); d++)
231             {
232                 ekind->grpstat[g].u[d] /=
233                     (1-lambda)*ekind->grpstat[g].mA + lambda*ekind->grpstat[g].mB;
234             }
235         }
236     }
237 }
238
239 real sum_ekin(t_grpopts *opts, gmx_ekindata_t *ekind, real *dekindlambda,
240               gmx_bool bEkinAveVel, gmx_bool bScaleEkin)
241 {
242     int           i, j, m, ngtc;
243     real          T, ek;
244     t_grp_tcstat *tcstat;
245     real          nrdf, nd, *ndf;
246
247     ngtc = opts->ngtc;
248     ndf  = opts->nrdf;
249
250     T    = 0;
251     nrdf = 0;
252
253     clear_mat(ekind->ekin);
254
255     for (i = 0; (i < ngtc); i++)
256     {
257
258         nd     = ndf[i];
259         tcstat = &ekind->tcstat[i];
260         /* Sometimes a group does not have degrees of freedom, e.g.
261          * when it consists of shells and virtual sites, then we just
262          * set the temperatue to 0 and also neglect the kinetic
263          * energy, which should be  zero anyway.
264          */
265
266         if (nd > 0)
267         {
268             if (bEkinAveVel)
269             {
270                 if (!bScaleEkin)
271                 {
272                     /* in this case, kinetic energy is from the current velocities already */
273                     msmul(tcstat->ekinf, tcstat->ekinscalef_nhc, tcstat->ekinf);
274                 }
275             }
276             else
277             {
278                 /* Calculate the full step Ekin as the average of the half steps */
279                 for (j = 0; (j < DIM); j++)
280                 {
281                     for (m = 0; (m < DIM); m++)
282                     {
283                         tcstat->ekinf[j][m] =
284                             0.5*(tcstat->ekinh[j][m]*tcstat->ekinscaleh_nhc + tcstat->ekinh_old[j][m]);
285                     }
286                 }
287             }
288             m_add(tcstat->ekinf, ekind->ekin, ekind->ekin);
289
290             tcstat->Th = calc_temp(trace(tcstat->ekinh), nd);
291             tcstat->T  = calc_temp(trace(tcstat->ekinf), nd);
292
293             /* after the scaling factors have been multiplied in, we can remove them */
294             if (bEkinAveVel)
295             {
296                 tcstat->ekinscalef_nhc = 1.0;
297             }
298             else
299             {
300                 tcstat->ekinscaleh_nhc = 1.0;
301             }
302         }
303         else
304         {
305             tcstat->T  = 0;
306             tcstat->Th = 0;
307         }
308         T    += nd*tcstat->T;
309         nrdf += nd;
310     }
311     if (nrdf > 0)
312     {
313         T /= nrdf;
314     }
315     if (dekindlambda)
316     {
317         if (bEkinAveVel)
318         {
319             *dekindlambda = ekind->dekindl;
320         }
321         else
322         {
323             *dekindlambda = 0.5*(ekind->dekindl + ekind->dekindl_old);
324         }
325     }
326     return T;
327 }