OCL: Make variables const
[gromacs.git] / src / gromacs / mdlib / stat.cpp
blobba9a135ecf37b075814de343b817d81d77b8e7d1
1 /*
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,2018, 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.
37 #include "gmxpre.h"
39 #include <cstdio>
40 #include <cstring>
42 #include "gromacs/domdec/domdec.h"
43 #include "gromacs/domdec/domdec_struct.h"
44 #include "gromacs/fileio/checkpoint.h"
45 #include "gromacs/fileio/xtcio.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/constr.h"
50 #include "gromacs/mdlib/md_support.h"
51 #include "gromacs/mdlib/mdrun.h"
52 #include "gromacs/mdlib/rbin.h"
53 #include "gromacs/mdlib/sim_util.h"
54 #include "gromacs/mdlib/tgroup.h"
55 #include "gromacs/mdlib/vcm.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/group.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.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
66 t_bin *rb;
67 int *itc0;
68 int *itc1;
69 } t_gmx_global_stat;
71 gmx_global_stat_t global_stat_init(const t_inputrec *ir)
73 gmx_global_stat_t gs;
75 snew(gs, 1);
77 gs->rb = mk_bin();
78 snew(gs->itc0, ir->opts.ngtc);
79 snew(gs->itc1, ir->opts.ngtc);
81 return gs;
84 void global_stat_destroy(gmx_global_stat_t gs)
86 destroy_bin(gs->rb);
87 sfree(gs->itc0);
88 sfree(gs->itc1);
89 sfree(gs);
92 static int filter_enerdterm(const real *afrom, gmx_bool bToBuffer, real *ato,
93 gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
95 int i, to, from;
97 from = 0;
98 to = 0;
99 for (i = 0; i < F_NRE; i++)
101 if (bToBuffer)
103 from = i;
105 else
107 to = i;
109 switch (i)
111 case F_EKIN:
112 case F_TEMP:
113 case F_DKDL:
114 if (bTemp)
116 ato[to++] = afrom[from++];
118 break;
119 case F_PRES:
120 case F_PDISPCORR:
121 if (bPres)
123 ato[to++] = afrom[from++];
125 break;
126 default:
127 if (bEner)
129 ato[to++] = afrom[from++];
131 break;
135 return to;
138 void global_stat(const gmx_global_stat *gs,
139 const t_commrec *cr, gmx_enerdata_t *enerd,
140 tensor fvir, tensor svir, rvec mu_tot,
141 const t_inputrec *inputrec,
142 gmx_ekindata_t *ekind, const gmx::Constraints *constr,
143 t_vcm *vcm,
144 int nsig, real *sig,
145 int *totalNumberOfBondedInteractions,
146 gmx_bool bSumEkinhOld, int flags)
147 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
149 t_bin *rb;
150 int *itc0, *itc1;
151 int ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
152 int idedl = 0, idedlo = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
153 int isig = -1;
154 int icj = -1, ici = -1, icx = -1;
155 int inn[egNR];
156 real copyenerd[F_NRE];
157 int nener, j;
158 double nb;
159 gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
160 bool checkNumberOfBondedInteractions = (flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0;
162 bVV = EI_VV(inputrec->eI);
163 bTemp = ((flags & CGLO_TEMPERATURE) != 0);
164 bEner = ((flags & CGLO_ENERGY) != 0);
165 bPres = ((flags & CGLO_PRESSURE) != 0);
166 bConstrVir = ((flags & CGLO_CONSTRAINT) != 0);
167 bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
168 bReadEkin = ((flags & CGLO_READEKIN) != 0);
170 rb = gs->rb;
171 itc0 = gs->itc0;
172 itc1 = gs->itc1;
175 reset_bin(rb);
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.*/
189 if (bConstrVir)
191 isv = add_binr(rb, DIM*DIM, svir[0]);
194 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
195 if (bTemp || !bVV)
197 if (ekind)
199 for (j = 0; (j < inputrec->opts.ngtc); j++)
201 if (bSumEkinhOld)
203 itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
205 if (bEkinAveVel && !bReadEkin)
207 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
209 else if (!bReadEkin)
211 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
214 /* these probably need to be put into one of these categories */
215 idedl = add_binr(rb, 1, &(ekind->dekindl));
216 if (bSumEkinhOld)
218 idedlo = add_binr(rb, 1, &(ekind->dekindl_old));
220 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
224 if (bPres)
226 ifv = add_binr(rb, DIM*DIM, fvir[0]);
229 gmx::ArrayRef<real> rmsdData;
230 if (bEner)
232 ie = add_binr(rb, nener, copyenerd);
233 if (constr)
235 rmsdData = constr->rmsdData();
236 if (!rmsdData.empty())
238 irmsd = add_binr(rb, 2, rmsdData.data());
241 if (!inputrecNeedMutot(inputrec))
243 imu = add_binr(rb, DIM, mu_tot);
246 for (j = 0; (j < egNR); j++)
248 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
250 if (inputrec->efep != efepNO)
252 idvdll = add_bind(rb, efptNR, enerd->dvdl_lin);
253 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
254 if (enerd->n_lambda > 0)
256 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
261 if (vcm)
263 icm = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
264 imass = add_binr(rb, vcm->nr, vcm->group_mass);
265 if (vcm->mode == ecmANGULAR)
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]);
273 if (checkNumberOfBondedInteractions)
275 nb = cr->dd->nbonded_local;
276 inb = add_bind(rb, 1, &nb);
278 if (nsig > 0)
280 isig = add_binr(rb, nsig, sig);
283 /* Global sum it all */
284 if (debug)
286 fprintf(debug, "Summing %d energies\n", rb->maxreal);
288 sum_bin(rb, cr);
290 /* Extract all the data locally */
292 if (bConstrVir)
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 */
298 if (bTemp || !bVV)
300 if (ekind)
302 for (j = 0; (j < inputrec->opts.ngtc); j++)
304 if (bSumEkinhOld)
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]);
312 else if (!bReadEkin)
314 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
317 extract_binr(rb, idedl, 1, &(ekind->dekindl));
318 if (bSumEkinhOld)
320 extract_binr(rb, idedlo, 1, &(ekind->dekindl_old));
322 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
325 if (bPres)
327 extract_binr(rb, ifv, DIM*DIM, fvir[0]);
330 if (bEner)
332 extract_binr(rb, ie, nener, copyenerd);
333 if (!rmsdData.empty())
335 extract_binr(rb, irmsd, rmsdData);
337 if (!inputrecNeedMutot(inputrec))
339 extract_binr(rb, imu, DIM, mu_tot);
342 for (j = 0; (j < egNR); j++)
344 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
346 if (inputrec->efep != efepNO)
348 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
349 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
350 if (enerd->n_lambda > 0)
352 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
356 filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
359 if (vcm)
361 extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
362 extract_binr(rb, imass, vcm->nr, vcm->group_mass);
363 if (vcm->mode == ecmANGULAR)
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]);
371 if (checkNumberOfBondedInteractions)
373 extract_bind(rb, inb, 1, &nb);
374 *totalNumberOfBondedInteractions = gmx::roundToInt(nb);
377 if (nsig > 0)
379 extract_binr(rb, isig, nsig, sig);
383 bool do_per_step(int64_t step, int64_t nstep)
385 if (nstep != 0)
387 return (step % nstep) == 0;
389 else
391 return false;