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