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