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