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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/domdec/domdec.h"
43 #include "gromacs/fileio/xtcio.h"
44 #include "gromacs/legacyheaders/checkpoint.h"
45 #include "gromacs/legacyheaders/constr.h"
46 #include "gromacs/legacyheaders/force.h"
47 #include "gromacs/legacyheaders/md_support.h"
48 #include "gromacs/legacyheaders/mdrun.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/network.h"
51 #include "gromacs/legacyheaders/rbin.h"
52 #include "gromacs/legacyheaders/sim_util.h"
53 #include "gromacs/legacyheaders/tgroup.h"
54 #include "gromacs/legacyheaders/txtdump.h"
55 #include "gromacs/legacyheaders/typedefs.h"
56 #include "gromacs/legacyheaders/vcm.h"
57 #include "gromacs/legacyheaders/types/commrec.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 typedef struct gmx_global_stat
71 gmx_global_stat_t global_stat_init(t_inputrec *ir)
78 snew(gs->itc0, ir->opts.ngtc);
79 snew(gs->itc1, ir->opts.ngtc);
84 void global_stat_destroy(gmx_global_stat_t gs)
92 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
93 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++];
129 ato[to++] = afrom[from++];
138 void global_stat(FILE *fplog, gmx_global_stat_t gs,
139 t_commrec *cr, gmx_enerdata_t *enerd,
140 tensor fvir, tensor svir, rvec mu_tot,
141 t_inputrec *inputrec,
142 gmx_ekindata_t *ekind, gmx_constr_t constr,
145 gmx_mtop_t *top_global, t_state *state_local,
146 gmx_bool bSumEkinhOld, int flags)
147 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
151 int ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
152 int idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
154 int icj = -1, ici = -1, icx = -1;
156 real copyenerd[F_NRE];
158 real *rmsd_data = NULL;
160 gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
162 bVV = EI_VV(inputrec->eI);
163 bTemp = flags & CGLO_TEMPERATURE;
164 bEner = flags & CGLO_ENERGY;
165 bPres = (flags & CGLO_PRESSURE);
166 bConstrVir = (flags & CGLO_CONSTRAINT);
167 bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
168 bReadEkin = (flags & CGLO_READEKIN);
176 /* This routine copies all the data to be summed to one big buffer
177 * using the t_bin struct.
180 /* First, we neeed to identify which enerd->term should be
181 communicated. Temperature and pressure terms should only be
182 communicated and summed when they need to be, to avoid repeating
183 the sums and overcounting. */
185 nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
187 /* First, the data that needs to be communicated with velocity verlet every time
188 This is just the constraint virial.*/
191 isv = add_binr(rb, DIM*DIM, svir[0]);
195 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
200 for (j = 0; (j < inputrec->opts.ngtc); j++)
204 itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
206 if (bEkinAveVel && !bReadEkin)
208 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
212 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
215 /* these probably need to be put into one of these categories */
217 idedl = add_binr(rb, 1, &(ekind->dekindl));
219 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
227 ifv = add_binr(rb, DIM*DIM, fvir[0]);
234 ie = add_binr(rb, nener, copyenerd);
238 rmsd_data = constr_rmsd_data(constr);
241 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
244 if (!NEED_MUTOT(*inputrec))
246 imu = add_binr(rb, DIM, mu_tot);
250 for (j = 0; (j < egNR); j++)
252 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
255 if (inputrec->efep != efepNO)
257 idvdll = add_bind(rb, efptNR, enerd->dvdl_lin);
258 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
259 if (enerd->n_lambda > 0)
261 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
268 icm = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
270 imass = add_binr(rb, vcm->nr, vcm->group_mass);
272 if (vcm->mode == ecmANGULAR)
274 icj = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
276 icx = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
278 ici = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
283 if (DOMAINDECOMP(cr))
285 nb = cr->dd->nbonded_local;
286 inb = add_bind(rb, 1, &nb);
291 isig = add_binr(rb, nsig, sig);
294 /* Global sum it all */
297 fprintf(debug, "Summing %d energies\n", rb->maxreal);
302 /* Extract all the data locally */
306 extract_binr(rb, isv, DIM*DIM, svir[0]);
309 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
314 for (j = 0; (j < inputrec->opts.ngtc); j++)
318 extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
320 if (bEkinAveVel && !bReadEkin)
322 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
326 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
329 extract_binr(rb, idedl, 1, &(ekind->dekindl));
330 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
336 extract_binr(rb, ifv, DIM*DIM, fvir[0]);
341 extract_binr(rb, ie, nener, copyenerd);
344 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
346 if (!NEED_MUTOT(*inputrec))
348 extract_binr(rb, imu, DIM, mu_tot);
351 for (j = 0; (j < egNR); j++)
353 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
355 if (inputrec->efep != efepNO)
357 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
358 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
359 if (enerd->n_lambda > 0)
361 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
364 if (DOMAINDECOMP(cr))
366 extract_bind(rb, inb, 1, &nb);
367 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
369 dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
374 filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
379 extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
381 extract_binr(rb, imass, vcm->nr, vcm->group_mass);
383 if (vcm->mode == ecmANGULAR)
385 extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
387 extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
389 extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
396 extract_binr(rb, isig, nsig, sig);
401 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
405 return ((step % nstep) == 0);