c1222041dfe63680f2f98a6c4f5fe5e39de4e9db
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "stat.h"
41
42 #include <cstdio>
43 #include <cstring>
44
45 #include "gromacs/domdec/domdec.h"
46 #include "gromacs/domdec/domdec_struct.h"
47 #include "gromacs/fileio/checkpoint.h"
48 #include "gromacs/fileio/xtcio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/constr.h"
53 #include "gromacs/mdlib/md_support.h"
54 #include "gromacs/mdlib/rbin.h"
55 #include "gromacs/mdlib/tgroup.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/enerdata.h"
59 #include "gromacs/mdtypes/group.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.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(const 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(const real* afrom, gmx_bool bToBuffer, real* ato, 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,
141                  gmx_enerdata_t*         enerd,
142                  tensor                  fvir,
143                  tensor                  svir,
144                  rvec                    mu_tot,
145                  const t_inputrec*       inputrec,
146                  gmx_ekindata_t*         ekind,
147                  const gmx::Constraints* constr,
148                  t_vcm*                  vcm,
149                  int                     nsig,
150                  real*                   sig,
151                  int*                    totalNumberOfBondedInteractions,
152                  gmx_bool                bSumEkinhOld,
153                  int                     flags)
154 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
155 {
156     t_bin* rb;
157     int *  itc0, *itc1;
158     int    ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
159     int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
160     int      isig = -1;
161     int      icj = -1, ici = -1, icx = -1;
162     int      inn[egNR];
163     real     copyenerd[F_NRE];
164     int      nener, j;
165     double   nb;
166     gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
167     bool checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
168
169     bVV         = EI_VV(inputrec->eI);
170     bTemp       = ((flags & CGLO_TEMPERATURE) != 0);
171     bEner       = ((flags & CGLO_ENERGY) != 0);
172     bPres       = ((flags & CGLO_PRESSURE) != 0);
173     bConstrVir  = ((flags & CGLO_CONSTRAINT) != 0);
174     bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
175     bReadEkin   = ((flags & CGLO_READEKIN) != 0);
176
177     rb   = gs->rb;
178     itc0 = gs->itc0;
179     itc1 = gs->itc1;
180
181
182     reset_bin(rb);
183     /* This routine copies all the data to be summed to one big buffer
184      * using the t_bin struct.
185      */
186
187     /* First, we neeed to identify which enerd->term should be
188        communicated.  Temperature and pressure terms should only be
189        communicated and summed when they need to be, to avoid repeating
190        the sums and overcounting. */
191
192     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
193
194     /* First, the data that needs to be communicated with velocity verlet every time
195        This is just the constraint virial.*/
196     if (bConstrVir)
197     {
198         isv = add_binr(rb, DIM * DIM, svir[0]);
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             idedl = add_binr(rb, 1, &(ekind->dekindl));
223             if (bSumEkinhOld)
224             {
225                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
226             }
227             if (ekind->cosacc.cos_accel != 0)
228             {
229                 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
230             }
231         }
232     }
233
234     if (bPres)
235     {
236         ifv = add_binr(rb, DIM * DIM, fvir[0]);
237     }
238
239     gmx::ArrayRef<real> rmsdData;
240     if (bEner)
241     {
242         ie = add_binr(rb, nener, copyenerd);
243         if (constr)
244         {
245             rmsdData = constr->rmsdData();
246             if (!rmsdData.empty())
247             {
248                 irmsd = add_binr(rb, 2, rmsdData.data());
249             }
250         }
251         if (!inputrecNeedMutot(inputrec))
252         {
253             imu = add_binr(rb, DIM, mu_tot);
254         }
255
256         for (j = 0; (j < egNR); j++)
257         {
258             inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j].data());
259         }
260         if (inputrec->efep != efepNO)
261         {
262             idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
263             idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
264             if (!enerd->enerpart_lambda.empty())
265             {
266                 iepl = add_bind(rb, enerd->enerpart_lambda.size(), enerd->enerpart_lambda.data());
267             }
268         }
269     }
270
271     if (vcm)
272     {
273         icm   = add_binr(rb, DIM * vcm->nr, vcm->group_p[0]);
274         imass = add_binr(rb, vcm->nr, vcm->group_mass.data());
275         if (vcm->mode == ecmANGULAR)
276         {
277             icj = add_binr(rb, DIM * vcm->nr, vcm->group_j[0]);
278             icx = add_binr(rb, DIM * vcm->nr, vcm->group_x[0]);
279             ici = add_binr(rb, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
280         }
281     }
282
283     if (checkNumberOfBondedInteractions)
284     {
285         nb  = cr->dd->nbonded_local;
286         inb = add_bind(rb, 1, &nb);
287     }
288     if (nsig > 0)
289     {
290         isig = add_binr(rb, nsig, sig);
291     }
292
293     /* Global sum it all */
294     if (debug)
295     {
296         fprintf(debug, "Summing %d energies\n", rb->maxreal);
297     }
298     sum_bin(rb, cr);
299
300     /* Extract all the data locally */
301
302     if (bConstrVir)
303     {
304         extract_binr(rb, isv, DIM * DIM, svir[0]);
305     }
306
307     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
308     if (bTemp || !bVV)
309     {
310         if (ekind)
311         {
312             for (j = 0; (j < inputrec->opts.ngtc); j++)
313             {
314                 if (bSumEkinhOld)
315                 {
316                     extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
317                 }
318                 if (bEkinAveVel && !bReadEkin)
319                 {
320                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
321                 }
322                 else if (!bReadEkin)
323                 {
324                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
325                 }
326             }
327             extract_binr(rb, idedl, 1, &(ekind->dekindl));
328             if (bSumEkinhOld)
329             {
330                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
331             }
332             if (ekind->cosacc.cos_accel != 0)
333             {
334                 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
335             }
336         }
337     }
338     if (bPres)
339     {
340         extract_binr(rb, ifv, DIM * DIM, fvir[0]);
341     }
342
343     if (bEner)
344     {
345         extract_binr(rb, ie, nener, copyenerd);
346         if (!rmsdData.empty())
347         {
348             extract_binr(rb, irmsd, rmsdData);
349         }
350         if (!inputrecNeedMutot(inputrec))
351         {
352             extract_binr(rb, imu, DIM, mu_tot);
353         }
354
355         for (j = 0; (j < egNR); j++)
356         {
357             extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j].data());
358         }
359         if (inputrec->efep != efepNO)
360         {
361             extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
362             extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
363             if (!enerd->enerpart_lambda.empty())
364             {
365                 extract_bind(rb, iepl, enerd->enerpart_lambda.size(), enerd->enerpart_lambda.data());
366             }
367         }
368
369         filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
370     }
371
372     if (vcm)
373     {
374         extract_binr(rb, icm, DIM * vcm->nr, vcm->group_p[0]);
375         extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
376         if (vcm->mode == ecmANGULAR)
377         {
378             extract_binr(rb, icj, DIM * vcm->nr, vcm->group_j[0]);
379             extract_binr(rb, icx, DIM * vcm->nr, vcm->group_x[0]);
380             extract_binr(rb, ici, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
381         }
382     }
383
384     if (checkNumberOfBondedInteractions)
385     {
386         extract_bind(rb, inb, 1, &nb);
387         *totalNumberOfBondedInteractions = gmx::roundToInt(nb);
388     }
389
390     if (nsig > 0)
391     {
392         extract_binr(rb, isig, nsig, sig);
393     }
394 }