Move physics.* to math/units.*
[gromacs.git] / src / gromacs / mdlib / stat.c
blob272271e31f9e72c31b45426af5c353795d553761
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include <string.h>
42 #include <stdio.h>
43 #include "typedefs.h"
44 #include "types/commrec.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "network.h"
47 #include "txtdump.h"
48 #include "names.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/math/utilities.h"
51 #include "force.h"
52 #include "vcm.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/utility/futil.h"
55 #include "network.h"
56 #include "rbin.h"
57 #include "tgroup.h"
58 #include "gromacs/fileio/xtcio.h"
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/fileio/trnio.h"
61 #include "domdec.h"
62 #include "constr.h"
63 #include "checkpoint.h"
64 #include "md_support.h"
65 #include "mdrun.h"
66 #include "sim_util.h"
68 typedef struct gmx_global_stat
70 t_bin *rb;
71 int *itc0;
72 int *itc1;
73 } t_gmx_global_stat;
75 gmx_global_stat_t global_stat_init(t_inputrec *ir)
77 gmx_global_stat_t gs;
79 snew(gs, 1);
81 gs->rb = mk_bin();
82 snew(gs->itc0, ir->opts.ngtc);
83 snew(gs->itc1, ir->opts.ngtc);
85 return gs;
88 void global_stat_destroy(gmx_global_stat_t gs)
90 destroy_bin(gs->rb);
91 sfree(gs->itc0);
92 sfree(gs->itc1);
93 sfree(gs);
96 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
97 gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
99 int i, to, from;
101 from = 0;
102 to = 0;
103 for (i = 0; i < F_NRE; i++)
105 if (bToBuffer)
107 from = i;
109 else
111 to = i;
113 switch (i)
115 case F_EKIN:
116 case F_TEMP:
117 case F_DKDL:
118 if (bTemp)
120 ato[to++] = afrom[from++];
122 break;
123 case F_PRES:
124 case F_PDISPCORR:
125 if (bPres)
127 ato[to++] = afrom[from++];
129 break;
130 default:
131 if (bEner)
133 ato[to++] = afrom[from++];
135 break;
139 return to;
142 void global_stat(FILE *fplog, gmx_global_stat_t gs,
143 t_commrec *cr, gmx_enerdata_t *enerd,
144 tensor fvir, tensor svir, rvec mu_tot,
145 t_inputrec *inputrec,
146 gmx_ekindata_t *ekind, gmx_constr_t constr,
147 t_vcm *vcm,
148 int nsig, real *sig,
149 gmx_mtop_t *top_global, t_state *state_local,
150 gmx_bool bSumEkinhOld, int flags)
151 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
153 t_bin *rb;
154 int *itc0, *itc1;
155 int ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
156 int idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
157 int isig = -1;
158 int icj = -1, ici = -1, icx = -1;
159 int inn[egNR];
160 real copyenerd[F_NRE];
161 int nener, j;
162 real *rmsd_data = NULL;
163 double nb;
164 gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bFirstIterate, bReadEkin;
166 bVV = EI_VV(inputrec->eI);
167 bTemp = flags & CGLO_TEMPERATURE;
168 bEner = flags & CGLO_ENERGY;
169 bPres = (flags & CGLO_PRESSURE);
170 bConstrVir = (flags & CGLO_CONSTRAINT);
171 bFirstIterate = (flags & CGLO_FIRSTITERATE);
172 bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
173 bReadEkin = (flags & CGLO_READEKIN);
175 rb = gs->rb;
176 itc0 = gs->itc0;
177 itc1 = gs->itc1;
180 reset_bin(rb);
181 /* This routine copies all the data to be summed to one big buffer
182 * using the t_bin struct.
185 /* First, we neeed to identify which enerd->term should be
186 communicated. Temperature and pressure terms should only be
187 communicated and summed when they need to be, to avoid repeating
188 the sums and overcounting. */
190 nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
192 /* First, the data that needs to be communicated with velocity verlet every time
193 This is just the constraint virial.*/
194 if (bConstrVir)
196 isv = add_binr(rb, DIM*DIM, svir[0]);
197 where();
200 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
201 if (bTemp || !bVV)
203 if (ekind)
205 for (j = 0; (j < inputrec->opts.ngtc); j++)
207 if (bSumEkinhOld)
209 itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
211 if (bEkinAveVel && !bReadEkin)
213 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
215 else if (!bReadEkin)
217 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
220 /* these probably need to be put into one of these categories */
221 where();
222 idedl = add_binr(rb, 1, &(ekind->dekindl));
223 where();
224 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
225 where();
228 where();
230 if ((bPres || !bVV) && bFirstIterate)
232 ifv = add_binr(rb, DIM*DIM, fvir[0]);
236 if (bEner)
238 where();
239 if (bFirstIterate)
241 ie = add_binr(rb, nener, copyenerd);
243 where();
244 if (constr)
246 rmsd_data = constr_rmsd_data(constr);
247 if (rmsd_data)
249 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
252 if (!NEED_MUTOT(*inputrec))
254 imu = add_binr(rb, DIM, mu_tot);
255 where();
258 if (bFirstIterate)
260 for (j = 0; (j < egNR); j++)
262 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
264 where();
265 if (inputrec->efep != efepNO)
267 idvdll = add_bind(rb, efptNR, enerd->dvdl_lin);
268 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
269 if (enerd->n_lambda > 0)
271 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
277 if (vcm)
279 icm = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
280 where();
281 imass = add_binr(rb, vcm->nr, vcm->group_mass);
282 where();
283 if (vcm->mode == ecmANGULAR)
285 icj = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
286 where();
287 icx = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
288 where();
289 ici = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
290 where();
294 if (DOMAINDECOMP(cr))
296 nb = cr->dd->nbonded_local;
297 inb = add_bind(rb, 1, &nb);
299 where();
300 if (nsig > 0)
302 isig = add_binr(rb, nsig, sig);
305 /* Global sum it all */
306 if (debug)
308 fprintf(debug, "Summing %d energies\n", rb->maxreal);
310 sum_bin(rb, cr);
311 where();
313 /* Extract all the data locally */
315 if (bConstrVir)
317 extract_binr(rb, isv, DIM*DIM, svir[0]);
320 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
321 if (bTemp || !bVV)
323 if (ekind)
325 for (j = 0; (j < inputrec->opts.ngtc); j++)
327 if (bSumEkinhOld)
329 extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
331 if (bEkinAveVel && !bReadEkin)
333 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
335 else if (!bReadEkin)
337 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
340 extract_binr(rb, idedl, 1, &(ekind->dekindl));
341 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
342 where();
345 if ((bPres || !bVV) && bFirstIterate)
347 extract_binr(rb, ifv, DIM*DIM, fvir[0]);
350 if (bEner)
352 if (bFirstIterate)
354 extract_binr(rb, ie, nener, copyenerd);
355 if (rmsd_data)
357 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
359 if (!NEED_MUTOT(*inputrec))
361 extract_binr(rb, imu, DIM, mu_tot);
364 for (j = 0; (j < egNR); j++)
366 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
368 if (inputrec->efep != efepNO)
370 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
371 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
372 if (enerd->n_lambda > 0)
374 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
377 if (DOMAINDECOMP(cr))
379 extract_bind(rb, inb, 1, &nb);
380 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
382 dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
385 where();
387 filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
391 if (vcm)
393 extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
394 where();
395 extract_binr(rb, imass, vcm->nr, vcm->group_mass);
396 where();
397 if (vcm->mode == ecmANGULAR)
399 extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
400 where();
401 extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
402 where();
403 extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
404 where();
408 if (nsig > 0)
410 extract_binr(rb, isig, nsig, sig);
412 where();
415 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
417 if (nstep != 0)
419 return ((step % nstep) == 0);
421 else
423 return 0;