2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
65 typedef struct gmx_global_stat
72 gmx_global_stat_t global_stat_init(const t_inputrec* ir)
79 snew(gs->itc0, ir->opts.ngtc);
80 snew(gs->itc1, ir->opts.ngtc);
85 void global_stat_destroy(gmx_global_stat_t gs)
93 static int filter_enerdterm(const real* afrom, gmx_bool bToBuffer, real* ato, gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
99 for (i = 0; i < F_NRE; i++)
116 ato[to++] = afrom[from++];
123 ato[to++] = afrom[from++];
128 // Don't reduce total and conserved energy
129 // because they are computed later (see #4301)
134 ato[to++] = afrom[from++];
143 void global_stat(const gmx_global_stat& gs,
145 gmx_enerdata_t* enerd,
148 const t_inputrec& inputrec,
149 gmx_ekindata_t* ekind,
150 gmx::ArrayRef<real> constraintsRmsdData,
152 gmx::ArrayRef<real> sig,
155 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
157 int ie = 0, ifv = 0, isv = 0, irmsd = 0;
158 int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
160 int icj = -1, ici = -1, icx = -1;
162 bool checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
163 bool bVV = EI_VV(inputrec.eI);
164 bool bTemp = ((flags & CGLO_TEMPERATURE) != 0);
165 bool bEner = ((flags & CGLO_ENERGY) != 0);
166 bool bPres = ((flags & CGLO_PRESSURE) != 0);
167 bool bConstrVir = ((flags & CGLO_CONSTRAINT) != 0);
168 bool bEkinAveVel = (inputrec.eI == IntegrationAlgorithm::VV
169 || (inputrec.eI == IntegrationAlgorithm::VVAK && bPres));
170 bool bReadEkin = ((flags & CGLO_READEKIN) != 0);
178 /* This routine copies all the data to be summed to one big buffer
179 * using the t_bin struct.
182 /* First, we neeed to identify which enerd->term should be
183 communicated. Temperature and pressure terms should only be
184 communicated and summed when they need to be, to avoid repeating
185 the sums and overcounting. */
187 std::array<real, F_NRE> copyenerd;
188 int nener = filter_enerdterm(enerd->term.data(), TRUE, copyenerd.data(), bTemp, bPres, bEner);
190 /* First, the data that needs to be communicated with velocity verlet every time
191 This is just the constraint virial.*/
194 isv = add_binr(rb, DIM * DIM, svir[0]);
197 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
202 for (int j = 0; (j < inputrec.opts.ngtc); j++)
206 itc0[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
208 if (bEkinAveVel && !bReadEkin)
210 itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinf[0]);
214 itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh[0]);
217 /* these probably need to be put into one of these categories */
218 idedl = add_binr(rb, 1, &(ekind->dekindl));
221 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
223 if (ekind->cosacc.cos_accel != 0)
225 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
232 ifv = add_binr(rb, DIM * DIM, fvir[0]);
235 gmx::EnumerationArray<NonBondedEnergyTerms, int> inn;
238 ie = add_binr(rb, nener, copyenerd.data());
239 if (!constraintsRmsdData.empty())
241 irmsd = add_binr(rb, 2, constraintsRmsdData.data());
243 for (auto key : gmx::keysOf(inn))
245 inn[key] = add_binr(rb, enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
247 if (inputrec.efep != FreeEnergyPerturbationType::No)
249 idvdll = add_bind(rb, enerd->dvdl_lin);
250 idvdlnl = add_bind(rb, enerd->dvdl_nonlin);
251 if (enerd->foreignLambdaTerms.numLambdas() > 0)
254 enerd->foreignLambdaTerms.energies().size(),
255 enerd->foreignLambdaTerms.energies().data());
262 icm = add_binr(rb, DIM * vcm->nr, vcm->group_p[0]);
263 imass = add_binr(rb, vcm->nr, vcm->group_mass.data());
264 if (vcm->mode == ComRemovalAlgorithm::Angular)
266 icj = add_binr(rb, DIM * vcm->nr, vcm->group_j[0]);
267 icx = add_binr(rb, DIM * vcm->nr, vcm->group_x[0]);
268 ici = add_binr(rb, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
273 if (checkNumberOfBondedInteractions)
275 GMX_RELEASE_ASSERT(DOMAINDECOMP(cr),
276 "No need to check number of bonded interactions when not using domain "
278 nb = numBondedInteractions(*cr->dd);
279 inb = add_bind(rb, 1, &nb);
283 isig = add_binr(rb, sig);
288 /* Extract all the data locally */
292 extract_binr(rb, isv, DIM * DIM, svir[0]);
295 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
300 for (int j = 0; (j < inputrec.opts.ngtc); j++)
304 extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
306 if (bEkinAveVel && !bReadEkin)
308 extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
312 extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
315 extract_binr(rb, idedl, 1, &(ekind->dekindl));
318 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
320 if (ekind->cosacc.cos_accel != 0)
322 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
328 extract_binr(rb, ifv, DIM * DIM, fvir[0]);
333 extract_binr(rb, ie, nener, copyenerd.data());
334 if (!constraintsRmsdData.empty())
336 extract_binr(rb, irmsd, constraintsRmsdData);
338 for (auto key : gmx::keysOf(inn))
340 extract_binr(rb, inn[key], enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
342 if (inputrec.efep != FreeEnergyPerturbationType::No)
344 extract_bind(rb, idvdll, enerd->dvdl_lin);
345 extract_bind(rb, idvdlnl, enerd->dvdl_nonlin);
346 if (enerd->foreignLambdaTerms.numLambdas() > 0)
350 enerd->foreignLambdaTerms.energies().size(),
351 enerd->foreignLambdaTerms.energies().data());
355 filter_enerdterm(copyenerd.data(), FALSE, enerd->term.data(), bTemp, bPres, bEner);
360 extract_binr(rb, icm, DIM * vcm->nr, vcm->group_p[0]);
361 extract_binr(rb, imass, vcm->nr, vcm->group_mass.data());
362 if (vcm->mode == ComRemovalAlgorithm::Angular)
364 extract_binr(rb, icj, DIM * vcm->nr, vcm->group_j[0]);
365 extract_binr(rb, icx, DIM * vcm->nr, vcm->group_x[0]);
366 extract_binr(rb, ici, DIM * DIM * vcm->nr, vcm->group_i[0][0]);
370 if (checkNumberOfBondedInteractions)
372 extract_bind(rb, inb, 1, &nb);
373 GMX_RELEASE_ASSERT(DOMAINDECOMP(cr),
374 "No need to check number of bonded interactions when not using domain "
376 setNumberOfBondedInteractionsOverAllDomains(*cr->dd, gmx::roundToInt(nb));
381 extract_binr(rb, isig, sig);