Use ObservablesReducer for LINCS RMSD computation
[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/mdtypes/observablesreducer.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             case F_ETOT:
128             case F_ECONSERVED:
129                 // Don't reduce total and conserved energy
130                 // because they are computed later (see #4301)
131                 break;
132             default:
133                 if (bEner)
134                 {
135                     ato[to++] = afrom[from++];
136                 }
137                 break;
138         }
139     }
140
141     return to;
142 }
143
144 void global_stat(const gmx_global_stat&   gs,
145                  const t_commrec*         cr,
146                  gmx_enerdata_t*          enerd,
147                  tensor                   fvir,
148                  tensor                   svir,
149                  const t_inputrec&        inputrec,
150                  gmx_ekindata_t*          ekind,
151                  t_vcm*                   vcm,
152                  gmx::ArrayRef<real>      sig,
153                  bool                     bSumEkinhOld,
154                  int                      flags,
155                  int64_t                  step,
156                  gmx::ObservablesReducer* observablesReducer)
157 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
158 {
159     int ie = 0, ifv = 0, isv = 0;
160     int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0;
161     int isig = -1;
162     int icj = -1, ici = -1, icx = -1;
163
164     bool bVV         = EI_VV(inputrec.eI);
165     bool bTemp       = ((flags & CGLO_TEMPERATURE) != 0);
166     bool bEner       = ((flags & CGLO_ENERGY) != 0);
167     bool bPres       = ((flags & CGLO_PRESSURE) != 0);
168     bool bConstrVir  = ((flags & CGLO_CONSTRAINT) != 0);
169     bool bEkinAveVel = (inputrec.eI == IntegrationAlgorithm::VV
170                         || (inputrec.eI == IntegrationAlgorithm::VVAK && bPres));
171     bool bReadEkin   = ((flags & CGLO_READEKIN) != 0);
172
173     // This structure implements something akin to a vector. As
174     // modules add their data into it with add_bin[rd], they save the
175     // index it returns, which allows them to look up their data later
176     // after the reduction with extract_bin[rd]. The various index
177     // variables are mostly named following the pattern
178     // "i<abbreviation for module>".
179     t_bin* rb   = gs.rb;
180     int*   itc0 = gs.itc0;
181     int*   itc1 = gs.itc1;
182
183
184     reset_bin(rb);
185     /* This routine copies all the data to be summed to one big buffer
186      * using the t_bin struct.
187      */
188
189     /* First, we neeed to identify which enerd->term should be
190        communicated.  Temperature and pressure terms should only be
191        communicated and summed when they need to be, to avoid repeating
192        the sums and overcounting. */
193
194     std::array<real, F_NRE> copyenerd;
195     int nener = filter_enerdterm(enerd->term.data(), TRUE, copyenerd.data(), bTemp, bPres, bEner);
196
197     /* First, the data that needs to be communicated with velocity verlet every time
198        This is just the constraint virial.*/
199     if (bConstrVir)
200     {
201         isv = add_binr(rb, DIM * DIM, svir[0]);
202     }
203
204     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
205     if (bTemp || !bVV)
206     {
207         if (ekind)
208         {
209             for (int j = 0; (j < inputrec.opts.ngtc); j++)
210             {
211                 if (bSumEkinhOld)
212                 {
213                     itc0[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
214                 }
215                 if (bEkinAveVel && !bReadEkin)
216                 {
217                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinf[0]);
218                 }
219                 else if (!bReadEkin)
220                 {
221                     itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh[0]);
222                 }
223             }
224             /* these probably need to be put into one of these categories */
225             idedl = add_binr(rb, 1, &(ekind->dekindl));
226             if (bSumEkinhOld)
227             {
228                 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
229             }
230             if (ekind->cosacc.cos_accel != 0)
231             {
232                 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
233             }
234         }
235     }
236
237     if (bPres)
238     {
239         ifv = add_binr(rb, DIM * DIM, fvir[0]);
240     }
241
242     gmx::EnumerationArray<NonBondedEnergyTerms, int> inn;
243     if (bEner)
244     {
245         ie = add_binr(rb, nener, copyenerd.data());
246         for (auto key : gmx::keysOf(inn))
247         {
248             inn[key] = add_binr(rb, enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].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 (!sig.empty())
276     {
277         isig = add_binr(rb, sig);
278     }
279
280     gmx::ArrayRef<double> observablesReducerBuffer = observablesReducer->communicationBuffer();
281     int                   tbinIndexForObservablesReducer = 0;
282     if (!observablesReducerBuffer.empty())
283     {
284         tbinIndexForObservablesReducer =
285                 add_bind(rb, observablesReducerBuffer.ssize(), observablesReducerBuffer.data());
286     }
287
288     sum_bin(rb, cr);
289
290     /* Extract all the data locally */
291
292     if (bConstrVir)
293     {
294         extract_binr(rb, isv, DIM * DIM, svir[0]);
295     }
296
297     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
298     if (bTemp || !bVV)
299     {
300         if (ekind)
301         {
302             for (int j = 0; (j < inputrec.opts.ngtc); j++)
303             {
304                 if (bSumEkinhOld)
305                 {
306                     extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
307                 }
308                 if (bEkinAveVel && !bReadEkin)
309                 {
310                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
311                 }
312                 else if (!bReadEkin)
313                 {
314                     extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
315                 }
316             }
317             extract_binr(rb, idedl, 1, &(ekind->dekindl));
318             if (bSumEkinhOld)
319             {
320                 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
321             }
322             if (ekind->cosacc.cos_accel != 0)
323             {
324                 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
325             }
326         }
327     }
328     if (bPres)
329     {
330         extract_binr(rb, ifv, DIM * DIM, fvir[0]);
331     }
332
333     if (bEner)
334     {
335         extract_binr(rb, ie, nener, copyenerd.data());
336         for (auto key : gmx::keysOf(inn))
337         {
338             extract_binr(rb, inn[key], enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
339         }
340         if (inputrec.efep != FreeEnergyPerturbationType::No)
341         {
342             extract_bind(rb, idvdll, enerd->dvdl_lin);
343             extract_bind(rb, idvdlnl, enerd->dvdl_nonlin);
344             if (enerd->foreignLambdaTerms.numLambdas() > 0)
345             {
346                 extract_bind(rb,
347                              iepl,
348                              enerd->foreignLambdaTerms.energies().size(),
349                              enerd->foreignLambdaTerms.energies().data());
350             }
351         }
352
353         filter_enerdterm(copyenerd.data(), FALSE, enerd->term.data(), bTemp, bPres, bEner);
354     }
355
356     if (vcm)
357     {
358         extract_binr(rb, icm, DIM * vcm->nr, vcm->group_p[0]);
359         extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
360         if (vcm->mode == ComRemovalAlgorithm::Angular)
361         {
362             extract_binr(rb, icj, DIM * vcm->nr, vcm->group_j[0]);
363             extract_binr(rb, icx, DIM * vcm->nr, vcm->group_x[0]);
364             extract_binr(rb, ici, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
365         }
366     }
367
368     if (!sig.empty())
369     {
370         extract_binr(rb, isig, sig);
371     }
372
373     if (!observablesReducerBuffer.empty())
374     {
375         extract_bind(rb,
376                      tbinIndexForObservablesReducer,
377                      observablesReducerBuffer.ssize(),
378                      observablesReducerBuffer.data());
379         observablesReducer->reductionComplete(step);
380     }
381 }