Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / stat.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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <string.h>
42 #include <stdio.h>
43 #include "typedefs.h"
44 #include "types/commrec.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "network.h"
47 #include "txtdump.h"
48 #include "names.h"
49 #include "physics.h"
50 #include "vec.h"
51 #include "gromacs/math/utilities.h"
52 #include "force.h"
53 #include "vcm.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/futil.h"
56 #include "network.h"
57 #include "rbin.h"
58 #include "tgroup.h"
59 #include "gromacs/fileio/xtcio.h"
60 #include "gromacs/fileio/gmxfio.h"
61 #include "gromacs/fileio/trnio.h"
62 #include "domdec.h"
63 #include "constr.h"
64 #include "checkpoint.h"
65 #include "md_support.h"
66 #include "mdrun.h"
67 #include "sim_util.h"
68
69 typedef struct gmx_global_stat
70 {
71     t_bin *rb;
72     int   *itc0;
73     int   *itc1;
74 } t_gmx_global_stat;
75
76 gmx_global_stat_t global_stat_init(t_inputrec *ir)
77 {
78     gmx_global_stat_t gs;
79
80     snew(gs, 1);
81
82     gs->rb = mk_bin();
83     snew(gs->itc0, ir->opts.ngtc);
84     snew(gs->itc1, ir->opts.ngtc);
85
86     return gs;
87 }
88
89 void global_stat_destroy(gmx_global_stat_t gs)
90 {
91     destroy_bin(gs->rb);
92     sfree(gs->itc0);
93     sfree(gs->itc1);
94     sfree(gs);
95 }
96
97 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
98                             gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
99 {
100     int i, to, from;
101
102     from = 0;
103     to   = 0;
104     for (i = 0; i < F_NRE; i++)
105     {
106         if (bToBuffer)
107         {
108             from = i;
109         }
110         else
111         {
112             to = i;
113         }
114         switch (i)
115         {
116             case F_EKIN:
117             case F_TEMP:
118             case F_DKDL:
119                 if (bTemp)
120                 {
121                     ato[to++] = afrom[from++];
122                 }
123                 break;
124             case F_PRES:
125             case F_PDISPCORR:
126                 if (bPres)
127                 {
128                     ato[to++] = afrom[from++];
129                 }
130                 break;
131             default:
132                 if (bEner)
133                 {
134                     ato[to++] = afrom[from++];
135                 }
136                 break;
137         }
138     }
139
140     return to;
141 }
142
143 void global_stat(FILE *fplog, gmx_global_stat_t gs,
144                  t_commrec *cr, gmx_enerdata_t *enerd,
145                  tensor fvir, tensor svir, rvec mu_tot,
146                  t_inputrec *inputrec,
147                  gmx_ekindata_t *ekind, gmx_constr_t constr,
148                  t_vcm *vcm,
149                  int nsig, real *sig,
150                  gmx_mtop_t *top_global, t_state *state_local,
151                  gmx_bool bSumEkinhOld, int flags)
152 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
153 {
154     t_bin     *rb;
155     int       *itc0, *itc1;
156     int        ie    = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
157     int        idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
158     int        isig  = -1;
159     int        icj   = -1, ici = -1, icx = -1;
160     int        inn[egNR];
161     real       copyenerd[F_NRE];
162     int        nener, j;
163     real      *rmsd_data = NULL;
164     double     nb;
165     gmx_bool   bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bFirstIterate, bReadEkin;
166
167     bVV           = EI_VV(inputrec->eI);
168     bTemp         = flags & CGLO_TEMPERATURE;
169     bEner         = flags & CGLO_ENERGY;
170     bPres         = (flags & CGLO_PRESSURE);
171     bConstrVir    = (flags & CGLO_CONSTRAINT);
172     bFirstIterate = (flags & CGLO_FIRSTITERATE);
173     bEkinAveVel   = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
174     bReadEkin     = (flags & CGLO_READEKIN);
175
176     rb   = gs->rb;
177     itc0 = gs->itc0;
178     itc1 = gs->itc1;
179
180
181     reset_bin(rb);
182     /* This routine copies all the data to be summed to one big buffer
183      * using the t_bin struct.
184      */
185
186     /* First, we neeed to identify which enerd->term should be
187        communicated.  Temperature and pressure terms should only be
188        communicated and summed when they need to be, to avoid repeating
189        the sums and overcounting. */
190
191     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
192
193     /* First, the data that needs to be communicated with velocity verlet every time
194        This is just the constraint virial.*/
195     if (bConstrVir)
196     {
197         isv = add_binr(rb, DIM*DIM, svir[0]);
198         where();
199     }
200
201 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
202     if (bTemp || !bVV)
203     {
204         if (ekind)
205         {
206             for (j = 0; (j < inputrec->opts.ngtc); j++)
207             {
208                 if (bSumEkinhOld)
209                 {
210                     itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
211                 }
212                 if (bEkinAveVel && !bReadEkin)
213                 {
214                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
215                 }
216                 else if (!bReadEkin)
217                 {
218                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
219                 }
220             }
221             /* these probably need to be put into one of these categories */
222             where();
223             idedl = add_binr(rb, 1, &(ekind->dekindl));
224             where();
225             ica   = add_binr(rb, 1, &(ekind->cosacc.mvcos));
226             where();
227         }
228     }
229     where();
230
231     if ((bPres || !bVV) && bFirstIterate)
232     {
233         ifv = add_binr(rb, DIM*DIM, fvir[0]);
234     }
235
236
237     if (bEner)
238     {
239         where();
240         if (bFirstIterate)
241         {
242             ie  = add_binr(rb, nener, copyenerd);
243         }
244         where();
245         if (constr)
246         {
247             rmsd_data = constr_rmsd_data(constr);
248             if (rmsd_data)
249             {
250                 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
251             }
252         }
253         if (!NEED_MUTOT(*inputrec))
254         {
255             imu = add_binr(rb, DIM, mu_tot);
256             where();
257         }
258
259         if (bFirstIterate)
260         {
261             for (j = 0; (j < egNR); j++)
262             {
263                 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
264             }
265             where();
266             if (inputrec->efep != efepNO)
267             {
268                 idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
269                 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
270                 if (enerd->n_lambda > 0)
271                 {
272                     iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
273                 }
274             }
275         }
276     }
277
278     if (vcm)
279     {
280         icm   = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
281         where();
282         imass = add_binr(rb, vcm->nr, vcm->group_mass);
283         where();
284         if (vcm->mode == ecmANGULAR)
285         {
286             icj   = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
287             where();
288             icx   = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
289             where();
290             ici   = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
291             where();
292         }
293     }
294
295     if (DOMAINDECOMP(cr))
296     {
297         nb  = cr->dd->nbonded_local;
298         inb = add_bind(rb, 1, &nb);
299     }
300     where();
301     if (nsig > 0)
302     {
303         isig = add_binr(rb, nsig, sig);
304     }
305
306     /* Global sum it all */
307     if (debug)
308     {
309         fprintf(debug, "Summing %d energies\n", rb->maxreal);
310     }
311     sum_bin(rb, cr);
312     where();
313
314     /* Extract all the data locally */
315
316     if (bConstrVir)
317     {
318         extract_binr(rb, isv, DIM*DIM, svir[0]);
319     }
320
321     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
322     if (bTemp || !bVV)
323     {
324         if (ekind)
325         {
326             for (j = 0; (j < inputrec->opts.ngtc); j++)
327             {
328                 if (bSumEkinhOld)
329                 {
330                     extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
331                 }
332                 if (bEkinAveVel && !bReadEkin)
333                 {
334                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
335                 }
336                 else if (!bReadEkin)
337                 {
338                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
339                 }
340             }
341             extract_binr(rb, idedl, 1, &(ekind->dekindl));
342             extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
343             where();
344         }
345     }
346     if ((bPres || !bVV) && bFirstIterate)
347     {
348         extract_binr(rb, ifv, DIM*DIM, fvir[0]);
349     }
350
351     if (bEner)
352     {
353         if (bFirstIterate)
354         {
355             extract_binr(rb, ie, nener, copyenerd);
356             if (rmsd_data)
357             {
358                 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
359             }
360             if (!NEED_MUTOT(*inputrec))
361             {
362                 extract_binr(rb, imu, DIM, mu_tot);
363             }
364
365             for (j = 0; (j < egNR); j++)
366             {
367                 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
368             }
369             if (inputrec->efep != efepNO)
370             {
371                 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
372                 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
373                 if (enerd->n_lambda > 0)
374                 {
375                     extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
376                 }
377             }
378             if (DOMAINDECOMP(cr))
379             {
380                 extract_bind(rb, inb, 1, &nb);
381                 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
382                 {
383                     dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
384                 }
385             }
386             where();
387
388             filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
389         }
390     }
391
392     if (vcm)
393     {
394         extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
395         where();
396         extract_binr(rb, imass, vcm->nr, vcm->group_mass);
397         where();
398         if (vcm->mode == ecmANGULAR)
399         {
400             extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
401             where();
402             extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
403             where();
404             extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
405             where();
406         }
407     }
408
409     if (nsig > 0)
410     {
411         extract_binr(rb, isig, nsig, sig);
412     }
413     where();
414 }
415
416 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
417 {
418     if (nstep != 0)
419     {
420         return ((step % nstep) == 0);
421     }
422     else
423     {
424         return 0;
425     }
426 }