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