Removed four includes from legacyheaders/typedefs.h
[gromacs.git] / src / gromacs / mdlib / stat.cpp
blobd60ceb6f763f0ae396f5bb25fa5a92f920f497c0
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, 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 <stdio.h>
40 #include <string.h>
42 #include "gromacs/domdec/domdec.h"
43 #include "gromacs/fileio/xtcio.h"
44 #include "gromacs/legacyheaders/checkpoint.h"
45 #include "gromacs/legacyheaders/force.h"
46 #include "gromacs/legacyheaders/md_support.h"
47 #include "gromacs/legacyheaders/mdrun.h"
48 #include "gromacs/legacyheaders/names.h"
49 #include "gromacs/legacyheaders/network.h"
50 #include "gromacs/legacyheaders/rbin.h"
51 #include "gromacs/legacyheaders/sim_util.h"
52 #include "gromacs/legacyheaders/tgroup.h"
53 #include "gromacs/legacyheaders/txtdump.h"
54 #include "gromacs/legacyheaders/typedefs.h"
55 #include "gromacs/legacyheaders/vcm.h"
56 #include "gromacs/legacyheaders/types/commrec.h"
57 #include "gromacs/legacyheaders/types/group.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/constr.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
67 t_bin *rb;
68 int *itc0;
69 int *itc1;
70 } t_gmx_global_stat;
72 gmx_global_stat_t global_stat_init(t_inputrec *ir)
74 gmx_global_stat_t gs;
76 snew(gs, 1);
78 gs->rb = mk_bin();
79 snew(gs->itc0, ir->opts.ngtc);
80 snew(gs->itc1, ir->opts.ngtc);
82 return gs;
85 void global_stat_destroy(gmx_global_stat_t gs)
87 destroy_bin(gs->rb);
88 sfree(gs->itc0);
89 sfree(gs->itc1);
90 sfree(gs);
93 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
94 gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
96 int i, to, from;
98 from = 0;
99 to = 0;
100 for (i = 0; i < F_NRE; i++)
102 if (bToBuffer)
104 from = i;
106 else
108 to = i;
110 switch (i)
112 case F_EKIN:
113 case F_TEMP:
114 case F_DKDL:
115 if (bTemp)
117 ato[to++] = afrom[from++];
119 break;
120 case F_PRES:
121 case F_PDISPCORR:
122 if (bPres)
124 ato[to++] = afrom[from++];
126 break;
127 default:
128 if (bEner)
130 ato[to++] = afrom[from++];
132 break;
136 return to;
139 void global_stat(FILE *fplog, gmx_global_stat_t gs,
140 t_commrec *cr, gmx_enerdata_t *enerd,
141 tensor fvir, tensor svir, rvec mu_tot,
142 t_inputrec *inputrec,
143 gmx_ekindata_t *ekind, gmx_constr_t constr,
144 t_vcm *vcm,
145 int nsig, real *sig,
146 gmx_mtop_t *top_global, t_state *state_local,
147 gmx_bool bSumEkinhOld, int flags)
148 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
150 t_bin *rb;
151 int *itc0, *itc1;
152 int ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
153 int idedl = 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 int inn[egNR];
157 real copyenerd[F_NRE];
158 int nener, j;
159 real *rmsd_data = NULL;
160 double nb;
161 gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bReadEkin;
163 bVV = EI_VV(inputrec->eI);
164 bTemp = flags & CGLO_TEMPERATURE;
165 bEner = flags & CGLO_ENERGY;
166 bPres = (flags & CGLO_PRESSURE);
167 bConstrVir = (flags & CGLO_CONSTRAINT);
168 bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
169 bReadEkin = (flags & CGLO_READEKIN);
171 rb = gs->rb;
172 itc0 = gs->itc0;
173 itc1 = gs->itc1;
176 reset_bin(rb);
177 /* This routine copies all the data to be summed to one big buffer
178 * using the t_bin struct.
181 /* First, we neeed to identify which enerd->term should be
182 communicated. Temperature and pressure terms should only be
183 communicated and summed when they need to be, to avoid repeating
184 the sums and overcounting. */
186 nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
188 /* First, the data that needs to be communicated with velocity verlet every time
189 This is just the constraint virial.*/
190 if (bConstrVir)
192 isv = add_binr(rb, DIM*DIM, svir[0]);
193 where();
196 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
197 if (bTemp || !bVV)
199 if (ekind)
201 for (j = 0; (j < inputrec->opts.ngtc); j++)
203 if (bSumEkinhOld)
205 itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
207 if (bEkinAveVel && !bReadEkin)
209 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
211 else if (!bReadEkin)
213 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
216 /* these probably need to be put into one of these categories */
217 where();
218 idedl = add_binr(rb, 1, &(ekind->dekindl));
219 where();
220 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
221 where();
224 where();
226 if (bPres || !bVV)
228 ifv = add_binr(rb, DIM*DIM, fvir[0]);
232 if (bEner)
234 where();
235 ie = add_binr(rb, nener, copyenerd);
236 where();
237 if (constr)
239 rmsd_data = constr_rmsd_data(constr);
240 if (rmsd_data)
242 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
245 if (!NEED_MUTOT(*inputrec))
247 imu = add_binr(rb, DIM, mu_tot);
248 where();
251 for (j = 0; (j < egNR); j++)
253 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
255 where();
256 if (inputrec->efep != efepNO)
258 idvdll = add_bind(rb, efptNR, enerd->dvdl_lin);
259 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
260 if (enerd->n_lambda > 0)
262 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
267 if (vcm)
269 icm = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
270 where();
271 imass = add_binr(rb, vcm->nr, vcm->group_mass);
272 where();
273 if (vcm->mode == ecmANGULAR)
275 icj = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
276 where();
277 icx = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
278 where();
279 ici = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
280 where();
284 if (DOMAINDECOMP(cr))
286 nb = cr->dd->nbonded_local;
287 inb = add_bind(rb, 1, &nb);
289 where();
290 if (nsig > 0)
292 isig = add_binr(rb, nsig, sig);
295 /* Global sum it all */
296 if (debug)
298 fprintf(debug, "Summing %d energies\n", rb->maxreal);
300 sum_bin(rb, cr);
301 where();
303 /* Extract all the data locally */
305 if (bConstrVir)
307 extract_binr(rb, isv, DIM*DIM, svir[0]);
310 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
311 if (bTemp || !bVV)
313 if (ekind)
315 for (j = 0; (j < inputrec->opts.ngtc); j++)
317 if (bSumEkinhOld)
319 extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
321 if (bEkinAveVel && !bReadEkin)
323 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
325 else if (!bReadEkin)
327 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
330 extract_binr(rb, idedl, 1, &(ekind->dekindl));
331 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
332 where();
335 if (bPres || !bVV)
337 extract_binr(rb, ifv, DIM*DIM, fvir[0]);
340 if (bEner)
342 extract_binr(rb, ie, nener, copyenerd);
343 if (rmsd_data)
345 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
347 if (!NEED_MUTOT(*inputrec))
349 extract_binr(rb, imu, DIM, mu_tot);
352 for (j = 0; (j < egNR); j++)
354 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
356 if (inputrec->efep != efepNO)
358 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
359 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
360 if (enerd->n_lambda > 0)
362 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
365 if (DOMAINDECOMP(cr))
367 extract_bind(rb, inb, 1, &nb);
368 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
370 dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
373 where();
375 filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
378 if (vcm)
380 extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
381 where();
382 extract_binr(rb, imass, vcm->nr, vcm->group_mass);
383 where();
384 if (vcm->mode == ecmANGULAR)
386 extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
387 where();
388 extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
389 where();
390 extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
391 where();
395 if (nsig > 0)
397 extract_binr(rb, isig, nsig, sig);
399 where();
402 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
404 if (nstep != 0)
406 return ((step % nstep) == 0);
408 else
410 return 0;