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