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