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