fa522f3537f461d01dddf430dbaccd98fca72de7
[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,2021, 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/md_support.h"
53 #include "gromacs/mdlib/rbin.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, 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(const gmx_global_stat* gs,
139                  const t_commrec*       cr,
140                  gmx_enerdata_t*        enerd,
141                  tensor                 fvir,
142                  tensor                 svir,
143                  const t_inputrec*      inputrec,
144                  gmx_ekindata_t*        ekind,
145                  gmx::ArrayRef<real>    constraintsRmsdData,
146                  t_vcm*                 vcm,
147                  int                    nsig,
148                  real*                  sig,
149                  int*                   totalNumberOfBondedInteractions,
150                  bool                   bSumEkinhOld,
151                  int                    flags)
152 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
153 {
154     t_bin* rb;
155     int *  itc0, *itc1;
156     int    ie = 0, ifv = 0, isv = 0, irmsd = 0;
157     int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
158     int      isig = -1;
159     int      icj = -1, ici = -1, icx = -1;
160     int      inn[egNR];
161     real     copyenerd[F_NRE];
162     int      nener, j;
163     double   nb;
164     gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
165     bool checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
166
167     bVV         = EI_VV(inputrec->eI);
168     bTemp       = ((flags & CGLO_TEMPERATURE) != 0);
169     bEner       = ((flags & CGLO_ENERGY) != 0);
170     bPres       = ((flags & CGLO_PRESSURE) != 0);
171     bConstrVir  = ((flags & CGLO_CONSTRAINT) != 0);
172     bEkinAveVel = (inputrec->eI == IntegrationAlgorithm::VV
173                    || (inputrec->eI == IntegrationAlgorithm::VVAK && bPres));
174     bReadEkin   = ((flags & CGLO_READEKIN) != 0);
175
176     rb   = gs->rb;
177     itc0 = gs->itc0;
178     itc1 = gs->itc1;
179
180
181     reset_bin(rb);
182     /* This routine copies all the data to be summed to one big buffer
183      * using the t_bin struct.
184      */
185
186     /* First, we neeed to identify which enerd->term should be
187        communicated.  Temperature and pressure terms should only be
188        communicated and summed when they need to be, to avoid repeating
189        the sums and overcounting. */
190
191     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
192
193     /* First, the data that needs to be communicated with velocity verlet every time
194        This is just the constraint virial.*/
195     if (bConstrVir)
196     {
197         isv = add_binr(rb, DIM * DIM, svir[0]);
198     }
199
200     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
201     if (bTemp || !bVV)
202     {
203         if (ekind)
204         {
205             for (j = 0; (j < inputrec->opts.ngtc); j++)
206             {
207                 if (bSumEkinhOld)
208                 {
209                     itc0[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
210                 }
211                 if (bEkinAveVel && !bReadEkin)
212                 {
213                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinf[0]);
214                 }
215                 else if (!bReadEkin)
216                 {
217                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh[0]);
218                 }
219             }
220             /* these probably need to be put into one of these categories */
221             idedl = add_binr(rb, 1, &(ekind->dekindl));
222             if (bSumEkinhOld)
223             {
224                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
225             }
226             if (ekind->cosacc.cos_accel != 0)
227             {
228                 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
229             }
230         }
231     }
232
233     if (bPres)
234     {
235         ifv = add_binr(rb, DIM * DIM, fvir[0]);
236     }
237
238     if (bEner)
239     {
240         ie = add_binr(rb, nener, copyenerd);
241         if (!constraintsRmsdData.empty())
242         {
243             irmsd = add_binr(rb, 2, constraintsRmsdData.data());
244         }
245
246         for (j = 0; (j < egNR); j++)
247         {
248             inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j].data());
249         }
250         if (inputrec->efep != FreeEnergyPerturbationType::No)
251         {
252             idvdll  = add_bind(rb, enerd->dvdl_lin);
253             idvdlnl = add_bind(rb, enerd->dvdl_nonlin);
254             if (enerd->foreignLambdaTerms.numLambdas() > 0)
255             {
256                 iepl = add_bind(rb,
257                                 enerd->foreignLambdaTerms.energies().size(),
258                                 enerd->foreignLambdaTerms.energies().data());
259             }
260         }
261     }
262
263     if (vcm)
264     {
265         icm   = add_binr(rb, DIM * vcm->nr, vcm->group_p[0]);
266         imass = add_binr(rb, vcm->nr, vcm->group_mass.data());
267         if (vcm->mode == ComRemovalAlgorithm::Angular)
268         {
269             icj = add_binr(rb, DIM * vcm->nr, vcm->group_j[0]);
270             icx = add_binr(rb, DIM * vcm->nr, vcm->group_x[0]);
271             ici = add_binr(rb, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
272         }
273     }
274
275     if (checkNumberOfBondedInteractions)
276     {
277         nb  = cr->dd->nbonded_local;
278         inb = add_bind(rb, 1, &nb);
279     }
280     if (nsig > 0)
281     {
282         isig = add_binr(rb, nsig, sig);
283     }
284
285     /* Global sum it all */
286     if (debug)
287     {
288         fprintf(debug, "Summing %d energies\n", rb->maxreal);
289     }
290     sum_bin(rb, cr);
291
292     /* Extract all the data locally */
293
294     if (bConstrVir)
295     {
296         extract_binr(rb, isv, DIM * DIM, svir[0]);
297     }
298
299     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
300     if (bTemp || !bVV)
301     {
302         if (ekind)
303         {
304             for (j = 0; (j < inputrec->opts.ngtc); j++)
305             {
306                 if (bSumEkinhOld)
307                 {
308                     extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
309                 }
310                 if (bEkinAveVel && !bReadEkin)
311                 {
312                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
313                 }
314                 else if (!bReadEkin)
315                 {
316                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
317                 }
318             }
319             extract_binr(rb, idedl, 1, &(ekind->dekindl));
320             if (bSumEkinhOld)
321             {
322                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
323             }
324             if (ekind->cosacc.cos_accel != 0)
325             {
326                 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
327             }
328         }
329     }
330     if (bPres)
331     {
332         extract_binr(rb, ifv, DIM * DIM, fvir[0]);
333     }
334
335     if (bEner)
336     {
337         extract_binr(rb, ie, nener, copyenerd);
338         if (!constraintsRmsdData.empty())
339         {
340             extract_binr(rb, irmsd, constraintsRmsdData);
341         }
342
343         for (j = 0; (j < egNR); j++)
344         {
345             extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j].data());
346         }
347         if (inputrec->efep != FreeEnergyPerturbationType::No)
348         {
349             extract_bind(rb, idvdll, enerd->dvdl_lin);
350             extract_bind(rb, idvdlnl, enerd->dvdl_nonlin);
351             if (enerd->foreignLambdaTerms.numLambdas() > 0)
352             {
353                 extract_bind(rb,
354                              iepl,
355                              enerd->foreignLambdaTerms.energies().size(),
356                              enerd->foreignLambdaTerms.energies().data());
357             }
358         }
359
360         filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
361     }
362
363     if (vcm)
364     {
365         extract_binr(rb, icm, DIM * vcm->nr, vcm->group_p[0]);
366         extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
367         if (vcm->mode == ComRemovalAlgorithm::Angular)
368         {
369             extract_binr(rb, icj, DIM * vcm->nr, vcm->group_j[0]);
370             extract_binr(rb, icx, DIM * vcm->nr, vcm->group_x[0]);
371             extract_binr(rb, ici, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
372         }
373     }
374
375     if (checkNumberOfBondedInteractions)
376     {
377         extract_bind(rb, inb, 1, &nb);
378         *totalNumberOfBondedInteractions = gmx::roundToInt(nb);
379     }
380
381     if (nsig > 0)
382     {
383         extract_binr(rb, isig, nsig, sig);
384     }
385 }