Revert "Eliminate mdlib/mdrun.h"
[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,2017,2018,2019, 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 <cstdio>
40 #include <cstring>
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/md_support.h"
51 #include "gromacs/mdlib/mdrun.h"
52 #include "gromacs/mdlib/rbin.h"
53 #include "gromacs/mdlib/sim_util.h"
54 #include "gromacs/mdlib/tgroup.h"
55 #include "gromacs/mdlib/vcm.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/enerdata.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(const 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(const 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(const gmx_global_stat *gs,
140                  const t_commrec *cr, gmx_enerdata_t *enerd,
141                  tensor fvir, tensor svir, rvec mu_tot,
142                  const t_inputrec *inputrec,
143                  gmx_ekindata_t *ekind, const gmx::Constraints *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     double     nb;
160     gmx_bool   bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
161     bool       checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
162
163     bVV           = EI_VV(inputrec->eI);
164     bTemp         = ((flags & CGLO_TEMPERATURE) != 0);
165     bEner         = ((flags & CGLO_ENERGY) != 0);
166     bPres         = ((flags & CGLO_PRESSURE) != 0);
167     bConstrVir    = ((flags & CGLO_CONSTRAINT) != 0);
168     bEkinAveVel   = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
169     bReadEkin     = ((flags & CGLO_READEKIN) != 0);
170
171     rb   = gs->rb;
172     itc0 = gs->itc0;
173     itc1 = gs->itc1;
174
175
176     reset_bin(rb);
177     /* This routine copies all the data to be summed to one big buffer
178      * using the t_bin struct.
179      */
180
181     /* First, we neeed to identify which enerd->term should be
182        communicated.  Temperature and pressure terms should only be
183        communicated and summed when they need to be, to avoid repeating
184        the sums and overcounting. */
185
186     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
187
188     /* First, the data that needs to be communicated with velocity verlet every time
189        This is just the constraint virial.*/
190     if (bConstrVir)
191     {
192         isv = add_binr(rb, DIM*DIM, svir[0]);
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             idedl = add_binr(rb, 1, &(ekind->dekindl));
217             if (bSumEkinhOld)
218             {
219                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
220             }
221             ica   = add_binr(rb, 1, &(ekind->cosacc.mvcos));
222         }
223     }
224
225     if (bPres)
226     {
227         ifv = add_binr(rb, DIM*DIM, fvir[0]);
228     }
229
230     gmx::ArrayRef<real> rmsdData;
231     if (bEner)
232     {
233         ie  = add_binr(rb, nener, copyenerd);
234         if (constr)
235         {
236             rmsdData = constr->rmsdData();
237             if (!rmsdData.empty())
238             {
239                 irmsd = add_binr(rb, 2, rmsdData.data());
240             }
241         }
242         if (!inputrecNeedMutot(inputrec))
243         {
244             imu = add_binr(rb, DIM, mu_tot);
245         }
246
247         for (j = 0; (j < egNR); j++)
248         {
249             inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
250         }
251         if (inputrec->efep != efepNO)
252         {
253             idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
254             idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
255             if (enerd->n_lambda > 0)
256             {
257                 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
258             }
259         }
260     }
261
262     if (vcm)
263     {
264         icm   = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
265         imass = add_binr(rb, vcm->nr, vcm->group_mass.data());
266         if (vcm->mode == ecmANGULAR)
267         {
268             icj   = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
269             icx   = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
270             ici   = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
271         }
272     }
273
274     if (checkNumberOfBondedInteractions)
275     {
276         nb  = cr->dd->nbonded_local;
277         inb = add_bind(rb, 1, &nb);
278     }
279     if (nsig > 0)
280     {
281         isig = add_binr(rb, nsig, sig);
282     }
283
284     /* Global sum it all */
285     if (debug)
286     {
287         fprintf(debug, "Summing %d energies\n", rb->maxreal);
288     }
289     sum_bin(rb, cr);
290
291     /* Extract all the data locally */
292
293     if (bConstrVir)
294     {
295         extract_binr(rb, isv, DIM*DIM, svir[0]);
296     }
297
298     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
299     if (bTemp || !bVV)
300     {
301         if (ekind)
302         {
303             for (j = 0; (j < inputrec->opts.ngtc); j++)
304             {
305                 if (bSumEkinhOld)
306                 {
307                     extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
308                 }
309                 if (bEkinAveVel && !bReadEkin)
310                 {
311                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
312                 }
313                 else if (!bReadEkin)
314                 {
315                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
316                 }
317             }
318             extract_binr(rb, idedl, 1, &(ekind->dekindl));
319             if (bSumEkinhOld)
320             {
321                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
322             }
323             extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
324         }
325     }
326     if (bPres)
327     {
328         extract_binr(rb, ifv, DIM*DIM, fvir[0]);
329     }
330
331     if (bEner)
332     {
333         extract_binr(rb, ie, nener, copyenerd);
334         if (!rmsdData.empty())
335         {
336             extract_binr(rb, irmsd, rmsdData);
337         }
338         if (!inputrecNeedMutot(inputrec))
339         {
340             extract_binr(rb, imu, DIM, mu_tot);
341         }
342
343         for (j = 0; (j < egNR); j++)
344         {
345             extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
346         }
347         if (inputrec->efep != efepNO)
348         {
349             extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
350             extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
351             if (enerd->n_lambda > 0)
352             {
353                 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
354             }
355         }
356
357         filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
358     }
359
360     if (vcm)
361     {
362         extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
363         extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
364         if (vcm->mode == ecmANGULAR)
365         {
366             extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
367             extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
368             extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
369         }
370     }
371
372     if (checkNumberOfBondedInteractions)
373     {
374         extract_bind(rb, inb, 1, &nb);
375         *totalNumberOfBondedInteractions = gmx::roundToInt(nb);
376     }
377
378     if (nsig > 0)
379     {
380         extract_binr(rb, isig, nsig, sig);
381     }
382 }
383
384 bool do_per_step(int64_t step, int64_t nstep)
385 {
386     if (nstep != 0)
387     {
388         return (step % nstep) == 0;
389     }
390     else
391     {
392         return false;
393     }
394 }