6be1c9d2cff94bdfdaf23c22f3e92a25a41339d9
[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/domdec/localtopologychecker.h"
48 #include "gromacs/fileio/checkpoint.h"
49 #include "gromacs/fileio/xtcio.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.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/mdtypes/observablesreducer.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/futil.h"
65 #include "gromacs/utility/smalloc.h"
66
67 typedef struct gmx_global_stat
68 {
69     t_bin* rb;
70     int*   itc0;
71     int*   itc1;
72 } t_gmx_global_stat;
73
74 gmx_global_stat_t global_stat_init(const t_inputrec* ir)
75 {
76     gmx_global_stat_t gs;
77
78     snew(gs, 1);
79
80     gs->rb = mk_bin();
81     snew(gs->itc0, ir->opts.ngtc);
82     snew(gs->itc1, ir->opts.ngtc);
83
84     return gs;
85 }
86
87 void global_stat_destroy(gmx_global_stat_t gs)
88 {
89     destroy_bin(gs->rb);
90     sfree(gs->itc0);
91     sfree(gs->itc1);
92     sfree(gs);
93 }
94
95 static int filter_enerdterm(const real* afrom, gmx_bool bToBuffer, real* ato, 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             case F_ETOT:
129             case F_ECONSERVED:
130                 // Don't reduce total and conserved energy
131                 // because they are computed later (see #4301)
132                 break;
133             default:
134                 if (bEner)
135                 {
136                     ato[to++] = afrom[from++];
137                 }
138                 break;
139         }
140     }
141
142     return to;
143 }
144
145 void global_stat(const gmx_global_stat&   gs,
146                  const t_commrec*         cr,
147                  gmx_enerdata_t*          enerd,
148                  tensor                   fvir,
149                  tensor                   svir,
150                  const t_inputrec&        inputrec,
151                  gmx_ekindata_t*          ekind,
152                  gmx::ArrayRef<real>      constraintsRmsdData,
153                  t_vcm*                   vcm,
154                  gmx::ArrayRef<real>      sig,
155                  bool                     bSumEkinhOld,
156                  int                      flags,
157                  int64_t                  step,
158                  gmx::ObservablesReducer* observablesReducer)
159 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
160 {
161     int ie = 0, ifv = 0, isv = 0, irmsd = 0;
162     int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
163     int isig = -1;
164     int icj = -1, ici = -1, icx = -1;
165
166     bool checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
167     bool bVV                             = EI_VV(inputrec.eI);
168     bool bTemp                           = ((flags & CGLO_TEMPERATURE) != 0);
169     bool bEner                           = ((flags & CGLO_ENERGY) != 0);
170     bool bPres                           = ((flags & CGLO_PRESSURE) != 0);
171     bool bConstrVir                      = ((flags & CGLO_CONSTRAINT) != 0);
172     bool bEkinAveVel                     = (inputrec.eI == IntegrationAlgorithm::VV
173                         || (inputrec.eI == IntegrationAlgorithm::VVAK && bPres));
174     bool bReadEkin                       = ((flags & CGLO_READEKIN) != 0);
175
176     // This structure implements something akin to a vector. As
177     // modules add their data into it with add_bin[rd], they save the
178     // index it returns, which allows them to look up their data later
179     // after the reduction with extract_bin[rd]. The various index
180     // variables are mostly named following the pattern
181     // "i<abbreviation for module>".
182     t_bin* rb   = gs.rb;
183     int*   itc0 = gs.itc0;
184     int*   itc1 = gs.itc1;
185
186
187     reset_bin(rb);
188     /* This routine copies all the data to be summed to one big buffer
189      * using the t_bin struct.
190      */
191
192     /* First, we neeed to identify which enerd->term should be
193        communicated.  Temperature and pressure terms should only be
194        communicated and summed when they need to be, to avoid repeating
195        the sums and overcounting. */
196
197     std::array<real, F_NRE> copyenerd;
198     int nener = filter_enerdterm(enerd->term.data(), TRUE, copyenerd.data(), bTemp, bPres, bEner);
199
200     /* First, the data that needs to be communicated with velocity verlet every time
201        This is just the constraint virial.*/
202     if (bConstrVir)
203     {
204         isv = add_binr(rb, DIM * DIM, svir[0]);
205     }
206
207     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
208     if (bTemp || !bVV)
209     {
210         if (ekind)
211         {
212             for (int j = 0; (j < inputrec.opts.ngtc); j++)
213             {
214                 if (bSumEkinhOld)
215                 {
216                     itc0[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
217                 }
218                 if (bEkinAveVel && !bReadEkin)
219                 {
220                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinf[0]);
221                 }
222                 else if (!bReadEkin)
223                 {
224                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh[0]);
225                 }
226             }
227             /* these probably need to be put into one of these categories */
228             idedl = add_binr(rb, 1, &(ekind->dekindl));
229             if (bSumEkinhOld)
230             {
231                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
232             }
233             if (ekind->cosacc.cos_accel != 0)
234             {
235                 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
236             }
237         }
238     }
239
240     if (bPres)
241     {
242         ifv = add_binr(rb, DIM * DIM, fvir[0]);
243     }
244
245     gmx::EnumerationArray<NonBondedEnergyTerms, int> inn;
246     if (bEner)
247     {
248         ie = add_binr(rb, nener, copyenerd.data());
249         if (!constraintsRmsdData.empty())
250         {
251             irmsd = add_binr(rb, 2, constraintsRmsdData.data());
252         }
253         for (auto key : gmx::keysOf(inn))
254         {
255             inn[key] = add_binr(rb, enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
256         }
257         if (inputrec.efep != FreeEnergyPerturbationType::No)
258         {
259             idvdll  = add_bind(rb, enerd->dvdl_lin);
260             idvdlnl = add_bind(rb, enerd->dvdl_nonlin);
261             if (enerd->foreignLambdaTerms.numLambdas() > 0)
262             {
263                 iepl = add_bind(rb,
264                                 enerd->foreignLambdaTerms.energies().size(),
265                                 enerd->foreignLambdaTerms.energies().data());
266             }
267         }
268     }
269
270     if (vcm)
271     {
272         icm   = add_binr(rb, DIM * vcm->nr, vcm->group_p[0]);
273         imass = add_binr(rb, vcm->nr, vcm->group_mass.data());
274         if (vcm->mode == ComRemovalAlgorithm::Angular)
275         {
276             icj = add_binr(rb, DIM * vcm->nr, vcm->group_j[0]);
277             icx = add_binr(rb, DIM * vcm->nr, vcm->group_x[0]);
278             ici = add_binr(rb, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
279         }
280     }
281
282     double nb;
283     if (checkNumberOfBondedInteractions)
284     {
285         GMX_RELEASE_ASSERT(DOMAINDECOMP(cr),
286                            "No need to check number of bonded interactions when not using domain "
287                            "decomposition");
288         nb  = cr->dd->localTopologyChecker->numBondedInteractions();
289         inb = add_bind(rb, 1, &nb);
290     }
291     if (!sig.empty())
292     {
293         isig = add_binr(rb, sig);
294     }
295
296     gmx::ArrayRef<double> observablesReducerBuffer = observablesReducer->communicationBuffer();
297     int                   tbinIndexForObservablesReducer = 0;
298     if (!observablesReducerBuffer.empty())
299     {
300         tbinIndexForObservablesReducer =
301                 add_bind(rb, observablesReducerBuffer.ssize(), observablesReducerBuffer.data());
302     }
303
304     sum_bin(rb, cr);
305
306     /* Extract all the data locally */
307
308     if (bConstrVir)
309     {
310         extract_binr(rb, isv, DIM * DIM, svir[0]);
311     }
312
313     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
314     if (bTemp || !bVV)
315     {
316         if (ekind)
317         {
318             for (int j = 0; (j < inputrec.opts.ngtc); j++)
319             {
320                 if (bSumEkinhOld)
321                 {
322                     extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
323                 }
324                 if (bEkinAveVel && !bReadEkin)
325                 {
326                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
327                 }
328                 else if (!bReadEkin)
329                 {
330                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
331                 }
332             }
333             extract_binr(rb, idedl, 1, &(ekind->dekindl));
334             if (bSumEkinhOld)
335             {
336                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
337             }
338             if (ekind->cosacc.cos_accel != 0)
339             {
340                 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
341             }
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.data());
352         if (!constraintsRmsdData.empty())
353         {
354             extract_binr(rb, irmsd, constraintsRmsdData);
355         }
356         for (auto key : gmx::keysOf(inn))
357         {
358             extract_binr(rb, inn[key], enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
359         }
360         if (inputrec.efep != FreeEnergyPerturbationType::No)
361         {
362             extract_bind(rb, idvdll, enerd->dvdl_lin);
363             extract_bind(rb, idvdlnl, enerd->dvdl_nonlin);
364             if (enerd->foreignLambdaTerms.numLambdas() > 0)
365             {
366                 extract_bind(rb,
367                              iepl,
368                              enerd->foreignLambdaTerms.energies().size(),
369                              enerd->foreignLambdaTerms.energies().data());
370             }
371         }
372
373         filter_enerdterm(copyenerd.data(), FALSE, enerd->term.data(), bTemp, bPres, bEner);
374     }
375
376     if (vcm)
377     {
378         extract_binr(rb, icm, DIM * vcm->nr, vcm->group_p[0]);
379         extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
380         if (vcm->mode == ComRemovalAlgorithm::Angular)
381         {
382             extract_binr(rb, icj, DIM * vcm->nr, vcm->group_j[0]);
383             extract_binr(rb, icx, DIM * vcm->nr, vcm->group_x[0]);
384             extract_binr(rb, ici, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
385         }
386     }
387
388     if (checkNumberOfBondedInteractions)
389     {
390         extract_bind(rb, inb, 1, &nb);
391         GMX_RELEASE_ASSERT(DOMAINDECOMP(cr),
392                            "No need to check number of bonded interactions when not using domain "
393                            "decomposition");
394         cr->dd->localTopologyChecker->setNumberOfBondedInteractionsOverAllDomains(gmx::roundToInt(nb));
395     }
396
397     if (!sig.empty())
398     {
399         extract_binr(rb, isig, sig);
400     }
401
402     if (!observablesReducerBuffer.empty())
403     {
404         extract_bind(rb,
405                      tbinIndexForObservablesReducer,
406                      observablesReducerBuffer.ssize(),
407                      observablesReducerBuffer.data());
408         observablesReducer->reductionComplete(step);
409     }
410 }