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,2016, 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/domdec/domdec_struct.h"
44 #include "gromacs/fileio/checkpoint.h"
45 #include "gromacs/fileio/xtcio.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/constr.h"
50 #include "gromacs/mdlib/force.h"
51 #include "gromacs/mdlib/md_support.h"
52 #include "gromacs/mdlib/mdrun.h"
53 #include "gromacs/mdlib/rbin.h"
54 #include "gromacs/mdlib/sim_util.h"
55 #include "gromacs/mdlib/tgroup.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/group.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
64
65 typedef struct gmx_global_stat
66 {
67     t_bin *rb;
68     int   *itc0;
69     int   *itc1;
70 } t_gmx_global_stat;
71
72 gmx_global_stat_t global_stat_init(t_inputrec *ir)
73 {
74     gmx_global_stat_t gs;
75
76     snew(gs, 1);
77
78     gs->rb = mk_bin();
79     snew(gs->itc0, ir->opts.ngtc);
80     snew(gs->itc1, ir->opts.ngtc);
81
82     return gs;
83 }
84
85 void global_stat_destroy(gmx_global_stat_t gs)
86 {
87     destroy_bin(gs->rb);
88     sfree(gs->itc0);
89     sfree(gs->itc1);
90     sfree(gs);
91 }
92
93 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
94                             gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
95 {
96     int i, to, from;
97
98     from = 0;
99     to   = 0;
100     for (i = 0; i < F_NRE; i++)
101     {
102         if (bToBuffer)
103         {
104             from = i;
105         }
106         else
107         {
108             to = i;
109         }
110         switch (i)
111         {
112             case F_EKIN:
113             case F_TEMP:
114             case F_DKDL:
115                 if (bTemp)
116                 {
117                     ato[to++] = afrom[from++];
118                 }
119                 break;
120             case F_PRES:
121             case F_PDISPCORR:
122                 if (bPres)
123                 {
124                     ato[to++] = afrom[from++];
125                 }
126                 break;
127             default:
128                 if (bEner)
129                 {
130                     ato[to++] = afrom[from++];
131                 }
132                 break;
133         }
134     }
135
136     return to;
137 }
138
139 void global_stat(gmx_global_stat_t gs,
140                  t_commrec *cr, gmx_enerdata_t *enerd,
141                  tensor fvir, tensor svir, rvec mu_tot,
142                  t_inputrec *inputrec,
143                  gmx_ekindata_t *ekind, gmx_constr_t constr,
144                  t_vcm *vcm,
145                  int nsig, real *sig,
146                  int *totalNumberOfBondedInteractions,
147                  gmx_bool bSumEkinhOld, int flags)
148 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
149 {
150     t_bin     *rb;
151     int       *itc0, *itc1;
152     int        ie    = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
153     int        idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
154     int        isig  = -1;
155     int        icj   = -1, ici = -1, icx = -1;
156     int        inn[egNR];
157     real       copyenerd[F_NRE];
158     int        nener, j;
159     real      *rmsd_data = NULL;
160     double     nb;
161     gmx_bool   bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
162     bool       checkNumberOfBondedInteractions = flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
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             if (bSumEkinhOld)
221             {
222                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
223             }
224             where();
225             ica   = add_binr(rb, 1, &(ekind->cosacc.mvcos));
226             where();
227         }
228     }
229     where();
230
231     if (bPres)
232     {
233         ifv = add_binr(rb, DIM*DIM, fvir[0]);
234     }
235
236
237     if (bEner)
238     {
239         where();
240         ie  = add_binr(rb, nener, copyenerd);
241         where();
242         if (constr)
243         {
244             rmsd_data = constr_rmsd_data(constr);
245             if (rmsd_data)
246             {
247                 irmsd = add_binr(rb, 2, rmsd_data);
248             }
249         }
250         if (!inputrecNeedMutot(inputrec))
251         {
252             imu = add_binr(rb, DIM, mu_tot);
253             where();
254         }
255
256         for (j = 0; (j < egNR); j++)
257         {
258             inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
259         }
260         where();
261         if (inputrec->efep != efepNO)
262         {
263             idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
264             idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
265             if (enerd->n_lambda > 0)
266             {
267                 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
268             }
269         }
270     }
271
272     if (vcm)
273     {
274         icm   = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
275         where();
276         imass = add_binr(rb, vcm->nr, vcm->group_mass);
277         where();
278         if (vcm->mode == ecmANGULAR)
279         {
280             icj   = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
281             where();
282             icx   = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
283             where();
284             ici   = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
285             where();
286         }
287     }
288
289     if (checkNumberOfBondedInteractions)
290     {
291         nb  = cr->dd->nbonded_local;
292         inb = add_bind(rb, 1, &nb);
293     }
294     where();
295     if (nsig > 0)
296     {
297         isig = add_binr(rb, nsig, sig);
298     }
299
300     /* Global sum it all */
301     if (debug)
302     {
303         fprintf(debug, "Summing %d energies\n", rb->maxreal);
304     }
305     sum_bin(rb, cr);
306     where();
307
308     /* Extract all the data locally */
309
310     if (bConstrVir)
311     {
312         extract_binr(rb, isv, DIM*DIM, svir[0]);
313     }
314
315     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
316     if (bTemp || !bVV)
317     {
318         if (ekind)
319         {
320             for (j = 0; (j < inputrec->opts.ngtc); j++)
321             {
322                 if (bSumEkinhOld)
323                 {
324                     extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
325                 }
326                 if (bEkinAveVel && !bReadEkin)
327                 {
328                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
329                 }
330                 else if (!bReadEkin)
331                 {
332                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
333                 }
334             }
335             extract_binr(rb, idedl, 1, &(ekind->dekindl));
336             if (bSumEkinhOld)
337             {
338                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
339             }
340             extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
341             where();
342         }
343     }
344     if (bPres)
345     {
346         extract_binr(rb, ifv, DIM*DIM, fvir[0]);
347     }
348
349     if (bEner)
350     {
351         extract_binr(rb, ie, nener, copyenerd);
352         if (rmsd_data)
353         {
354             extract_binr(rb, irmsd, 2, rmsd_data);
355         }
356         if (!inputrecNeedMutot(inputrec))
357         {
358             extract_binr(rb, imu, DIM, mu_tot);
359         }
360
361         for (j = 0; (j < egNR); j++)
362         {
363             extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
364         }
365         if (inputrec->efep != efepNO)
366         {
367             extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
368             extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
369             if (enerd->n_lambda > 0)
370             {
371                 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
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 (checkNumberOfBondedInteractions)
397     {
398         extract_bind(rb, inb, 1, &nb);
399         *totalNumberOfBondedInteractions = static_cast<int>(nb+0.5);
400     }
401
402     if (nsig > 0)
403     {
404         extract_binr(rb, isig, nsig, sig);
405     }
406     where();
407 }
408
409 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
410 {
411     if (nstep != 0)
412     {
413         return ((step % nstep) == 0);
414     }
415     else
416     {
417         return 0;
418     }
419 }