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/mdtypes/observablesreducer.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
66 typedef struct gmx_global_stat
73 gmx_global_stat_t global_stat_init(const t_inputrec* ir)
80 snew(gs->itc0, ir->opts.ngtc);
81 snew(gs->itc1, ir->opts.ngtc);
86 void global_stat_destroy(gmx_global_stat_t gs)
94 static int filter_enerdterm(const real* afrom, gmx_bool bToBuffer, real* ato, gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
100 for (i = 0; i < F_NRE; i++)
117 ato[to++] = afrom[from++];
124 ato[to++] = afrom[from++];
129 // Don't reduce total and conserved energy
130 // because they are computed later (see #4301)
135 ato[to++] = afrom[from++];
144 void global_stat(const gmx_global_stat& gs,
146 gmx_enerdata_t* enerd,
149 const t_inputrec& inputrec,
150 gmx_ekindata_t* ekind,
152 gmx::ArrayRef<real> sig,
156 gmx::ObservablesReducer* observablesReducer)
157 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
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;
162 int icj = -1, ici = -1, icx = -1;
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);
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>".
185 /* This routine copies all the data to be summed to one big buffer
186 * using the t_bin struct.
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. */
194 std::array<real, F_NRE> copyenerd;
195 int nener = filter_enerdterm(enerd->term.data(), TRUE, copyenerd.data(), bTemp, bPres, bEner);
197 /* First, the data that needs to be communicated with velocity verlet every time
198 This is just the constraint virial.*/
201 isv = add_binr(rb, DIM * DIM, svir[0]);
204 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
209 for (int j = 0; (j < inputrec.opts.ngtc); j++)
213 itc0[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
215 if (bEkinAveVel && !bReadEkin)
217 itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinf[0]);
221 itc1[j] = add_binr(rb, DIM * DIM, ekind->tcstat[j].ekinh[0]);
224 /* these probably need to be put into one of these categories */
225 idedl = add_binr(rb, 1, &(ekind->dekindl));
228 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
230 if (ekind->cosacc.cos_accel != 0)
232 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
239 ifv = add_binr(rb, DIM * DIM, fvir[0]);
242 gmx::EnumerationArray<NonBondedEnergyTerms, int> inn;
245 ie = add_binr(rb, nener, copyenerd.data());
246 for (auto key : gmx::keysOf(inn))
248 inn[key] = add_binr(rb, enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
250 if (inputrec.efep != FreeEnergyPerturbationType::No)
252 idvdll = add_bind(rb, enerd->dvdl_lin);
253 idvdlnl = add_bind(rb, enerd->dvdl_nonlin);
254 if (enerd->foreignLambdaTerms.numLambdas() > 0)
257 enerd->foreignLambdaTerms.energies().size(),
258 enerd->foreignLambdaTerms.energies().data());
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)
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]);
277 isig = add_binr(rb, sig);
280 gmx::ArrayRef<double> observablesReducerBuffer = observablesReducer->communicationBuffer();
281 int tbinIndexForObservablesReducer = 0;
282 if (!observablesReducerBuffer.empty())
284 tbinIndexForObservablesReducer =
285 add_bind(rb, observablesReducerBuffer.ssize(), observablesReducerBuffer.data());
290 /* Extract all the data locally */
294 extract_binr(rb, isv, DIM * DIM, svir[0]);
297 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
302 for (int j = 0; (j < inputrec.opts.ngtc); j++)
306 extract_binr(rb, itc0[j], DIM * DIM, ekind->tcstat[j].ekinh_old[0]);
308 if (bEkinAveVel && !bReadEkin)
310 extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinf[0]);
314 extract_binr(rb, itc1[j], DIM * DIM, ekind->tcstat[j].ekinh[0]);
317 extract_binr(rb, idedl, 1, &(ekind->dekindl));
320 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
322 if (ekind->cosacc.cos_accel != 0)
324 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
330 extract_binr(rb, ifv, DIM * DIM, fvir[0]);
335 extract_binr(rb, ie, nener, copyenerd.data());
336 for (auto key : gmx::keysOf(inn))
338 extract_binr(rb, inn[key], enerd->grpp.nener, enerd->grpp.energyGroupPairTerms[key].data());
340 if (inputrec.efep != FreeEnergyPerturbationType::No)
342 extract_bind(rb, idvdll, enerd->dvdl_lin);
343 extract_bind(rb, idvdlnl, enerd->dvdl_nonlin);
344 if (enerd->foreignLambdaTerms.numLambdas() > 0)
348 enerd->foreignLambdaTerms.energies().size(),
349 enerd->foreignLambdaTerms.energies().data());
353 filter_enerdterm(copyenerd.data(), FALSE, enerd->term.data(), bTemp, bPres, bEner);
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)
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]);
370 extract_binr(rb, isig, sig);
373 if (!observablesReducerBuffer.empty())
376 tbinIndexForObservablesReducer,
377 observablesReducerBuffer.ssize(),
378 observablesReducerBuffer.data());
379 observablesReducer->reductionComplete(step);