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