Merge branch 'rotation-4-5' into rotation
[gromacs/adressmacs.git] / src / mdlib / stat.c
blobddf1d3b188fca4430d56e74889b965672e70daa9
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <string.h>
41 #include <stdio.h>
42 #include "typedefs.h"
43 #include "sysstuff.h"
44 #include "gmx_fatal.h"
45 #include "network.h"
46 #include "txtdump.h"
47 #include "names.h"
48 #include "physics.h"
49 #include "vec.h"
50 #include "maths.h"
51 #include "mvdata.h"
52 #include "main.h"
53 #include "force.h"
54 #include "vcm.h"
55 #include "smalloc.h"
56 #include "futil.h"
57 #include "network.h"
58 #include "rbin.h"
59 #include "tgroup.h"
60 #include "xtcio.h"
61 #include "gmxfio.h"
62 #include "trnio.h"
63 #include "statutil.h"
64 #include "domdec.h"
65 #include "partdec.h"
66 #include "constr.h"
67 #include "checkpoint.h"
68 #include "mdrun.h"
69 #include "xvgr.h"
71 typedef struct gmx_global_stat
73 t_bin *rb;
74 int *itc0;
75 int *itc1;
76 } t_gmx_global_stat;
78 gmx_global_stat_t global_stat_init(t_inputrec *ir)
80 gmx_global_stat_t gs;
82 snew(gs,1);
84 gs->rb = mk_bin();
85 snew(gs->itc0,ir->opts.ngtc);
86 snew(gs->itc1,ir->opts.ngtc);
88 return gs;
91 void global_stat_destroy(gmx_global_stat_t gs)
93 destroy_bin(gs->rb);
94 sfree(gs->itc0);
95 sfree(gs->itc1);
96 sfree(gs);
99 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
100 gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner) {
101 int i,to,from;
103 from = 0;
104 to = 0;
105 for (i=0;i<F_NRE;i++)
107 if (bToBuffer)
109 from = i;
111 else
113 to = i;
115 switch (i) {
116 case F_EKIN:
117 case F_TEMP:
118 case F_DKDL:
119 if (bTemp)
121 ato[to++] = afrom[from++];
123 break;
124 case F_PRES:
125 case F_PDISPCORR:
126 case F_VTEMP:
127 if (bPres)
129 ato[to++] = afrom[from++];
131 break;
132 default:
133 if (bEner)
135 ato[to++] = afrom[from++];
137 break;
141 return to;
144 void global_stat(FILE *fplog,gmx_global_stat_t gs,
145 t_commrec *cr,gmx_enerdata_t *enerd,
146 tensor fvir,tensor svir,rvec mu_tot,
147 t_inputrec *inputrec,
148 gmx_ekindata_t *ekind,gmx_constr_t constr,
149 t_vcm *vcm,
150 int nsig,real *sig,
151 gmx_mtop_t *top_global, t_state *state_local,
152 gmx_bool bSumEkinhOld, int flags)
153 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
155 t_bin *rb;
156 int *itc0,*itc1;
157 int ie=0,ifv=0,isv=0,irmsd=0,imu=0;
158 int idedl=0,idvdll=0,idvdlnl=0,iepl=0,icm=0,imass=0,ica=0,inb=0;
159 int isig=-1;
160 int icj=-1,ici=-1,icx=-1;
161 int inn[egNR];
162 real copyenerd[F_NRE];
163 int nener,j;
164 real *rmsd_data=NULL;
165 double nb;
166 gmx_bool bVV,bTemp,bEner,bPres,bConstrVir,bEkinAveVel,bFirstIterate,bReadEkin;
168 bVV = EI_VV(inputrec->eI);
169 bTemp = flags & CGLO_TEMPERATURE;
170 bEner = flags & CGLO_ENERGY;
171 bPres = (flags & CGLO_PRESSURE);
172 bConstrVir = (flags & CGLO_CONSTRAINT);
173 bFirstIterate = (flags & CGLO_FIRSTITERATE);
174 bEkinAveVel = (inputrec->eI==eiVV || (inputrec->eI==eiVVAK && bPres));
175 bReadEkin = (flags & CGLO_READEKIN);
177 rb = gs->rb;
178 itc0 = gs->itc0;
179 itc1 = gs->itc1;
182 reset_bin(rb);
183 /* This routine copies all the data to be summed to one big buffer
184 * using the t_bin struct.
187 /* First, we neeed to identify which enerd->term should be
188 communicated. Temperature and pressure terms should only be
189 communicated and summed when they need to be, to avoid repeating
190 the sums and overcounting. */
192 nener = filter_enerdterm(enerd->term,TRUE,copyenerd,bTemp,bPres,bEner);
194 /* First, the data that needs to be communicated with velocity verlet every time
195 This is just the constraint virial.*/
196 if (bConstrVir) {
197 isv = add_binr(rb,DIM*DIM,svir[0]);
198 where();
201 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
202 if (bTemp || !bVV)
204 if (ekind)
206 for(j=0; (j<inputrec->opts.ngtc); j++)
208 if (bSumEkinhOld)
210 itc0[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
212 if (bEkinAveVel && !bReadEkin)
214 itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinf[0]);
216 else if (!bReadEkin)
218 itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh[0]);
221 /* these probably need to be put into one of these categories */
222 where();
223 idedl = add_binr(rb,1,&(ekind->dekindl));
224 where();
225 ica = add_binr(rb,1,&(ekind->cosacc.mvcos));
226 where();
229 where();
231 if ((bPres || !bVV) && bFirstIterate)
233 ifv = add_binr(rb,DIM*DIM,fvir[0]);
237 if (bEner)
239 where();
240 if (bFirstIterate)
242 ie = add_binr(rb,nener,copyenerd);
244 where();
245 if (constr)
247 rmsd_data = constr_rmsd_data(constr);
248 if (rmsd_data)
250 irmsd = add_binr(rb,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
253 if (!NEED_MUTOT(*inputrec))
255 imu = add_binr(rb,DIM,mu_tot);
256 where();
259 if (bFirstIterate)
261 for(j=0; (j<egNR); j++)
263 inn[j]=add_binr(rb,enerd->grpp.nener,enerd->grpp.ener[j]);
265 where();
266 if (inputrec->efep != efepNO)
268 idvdll = add_bind(rb,1,&enerd->dvdl_lin);
269 idvdlnl = add_bind(rb,1,&enerd->dvdl_nonlin);
270 if (enerd->n_lambda > 0)
272 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) {
332 extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinf[0]);
334 else if (!bReadEkin)
336 extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinh[0]);
339 extract_binr(rb,idedl,1,&(ekind->dekindl));
340 extract_binr(rb,ica,1,&(ekind->cosacc.mvcos));
341 where();
344 if ((bPres || !bVV) && bFirstIterate)
346 extract_binr(rb,ifv ,DIM*DIM,fvir[0]);
349 if (bEner)
351 if (bFirstIterate)
353 extract_binr(rb,ie,nener,copyenerd);
354 if (rmsd_data)
356 extract_binr(rb,irmsd,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
358 if (!NEED_MUTOT(*inputrec))
360 extract_binr(rb,imu,DIM,mu_tot);
363 for(j=0; (j<egNR); j++)
365 extract_binr(rb,inn[j],enerd->grpp.nener,enerd->grpp.ener[j]);
367 if (inputrec->efep != efepNO)
369 extract_bind(rb,idvdll ,1,&enerd->dvdl_lin);
370 extract_bind(rb,idvdlnl,1,&enerd->dvdl_nonlin);
371 if (enerd->n_lambda > 0)
373 extract_bind(rb,iepl,enerd->n_lambda,enerd->enerpart_lambda);
376 /* should this be here, or with ekin?*/
377 if (vcm)
379 extract_binr(rb,icm,DIM*vcm->nr,vcm->group_p[0]);
380 where();
381 extract_binr(rb,imass,vcm->nr,vcm->group_mass);
382 where();
383 if (vcm->mode == ecmANGULAR)
385 extract_binr(rb,icj,DIM*vcm->nr,vcm->group_j[0]);
386 where();
387 extract_binr(rb,icx,DIM*vcm->nr,vcm->group_x[0]);
388 where();
389 extract_binr(rb,ici,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
390 where();
393 if (DOMAINDECOMP(cr))
395 extract_bind(rb,inb,1,&nb);
396 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
398 dd_print_missing_interactions(fplog,cr,(int)(nb + 0.5),top_global,state_local);
401 where();
403 filter_enerdterm(copyenerd,FALSE,enerd->term,bTemp,bPres,bEner);
404 /* Small hack for temp only - not entirely clear if still needed?*/
405 /* enerd->term[F_TEMP] /= (cr->nnodes - cr->npmenodes); */
409 if (nsig > 0)
411 extract_binr(rb,isig,nsig,sig);
413 where();
416 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep)
418 if (nstep != 0)
419 return ((step % nstep)==0);
420 else
421 return 0;
424 static void moveit(t_commrec *cr,
425 int left,int right,const char *s,rvec xx[])
427 if (!xx)
428 return;
430 move_rvecs(cr,FALSE,FALSE,left,right,
431 xx,NULL,(cr->nnodes-cr->npmenodes)-1,NULL);
434 gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],int mdrun_flags,
435 const t_commrec *cr,const t_inputrec *ir,
436 const output_env_t oenv)
438 gmx_mdoutf_t *of;
439 char filemode[3];
440 gmx_bool bAppendFiles;
442 snew(of,1);
444 of->fp_trn = NULL;
445 of->fp_ene = NULL;
446 of->fp_xtc = NULL;
447 of->fp_dhdl = NULL;
448 of->fp_field = NULL;
450 of->eIntegrator = ir->eI;
451 of->simulation_part = ir->simulation_part;
453 if (MASTER(cr))
455 bAppendFiles = (mdrun_flags & MD_APPENDFILES);
457 of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT);
459 sprintf(filemode, bAppendFiles ? "a+" : "w+");
461 if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
462 #ifndef GMX_FAHCORE
464 !(EI_DYNAMICS(ir->eI) &&
465 ir->nstxout == 0 &&
466 ir->nstvout == 0 &&
467 ir->nstfout == 0)
468 #endif
471 of->fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
473 if (EI_DYNAMICS(ir->eI) &&
474 ir->nstxtcout > 0)
476 of->fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
477 of->xtc_prec = ir->xtcprec;
479 if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
481 of->fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
483 of->fn_cpt = opt2fn("-cpo",nfile,fnm);
485 if (ir->efep != efepNO && ir->nstdhdl > 0 &&
486 (ir->separate_dhdl_file == sepdhdlfileYES ) &&
487 EI_DYNAMICS(ir->eI))
489 if (bAppendFiles)
491 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
493 else
495 of->fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
499 if (opt2bSet("-field",nfile,fnm) &&
500 (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
502 if (bAppendFiles)
504 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field",nfile,fnm),
505 filemode);
507 else
509 of->fp_field = xvgropen(opt2fn("-field",nfile,fnm),
510 "Applied electric field","Time (ps)",
511 "E (V/nm)",oenv);
516 return of;
519 void done_mdoutf(gmx_mdoutf_t *of)
521 if (of->fp_ene != NULL)
523 close_enx(of->fp_ene);
525 if (of->fp_xtc)
527 close_xtc(of->fp_xtc);
529 if (of->fp_trn)
531 close_trn(of->fp_trn);
533 if (of->fp_dhdl != NULL)
535 gmx_fio_fclose(of->fp_dhdl);
537 if (of->fp_field != NULL)
539 gmx_fio_fclose(of->fp_field);
542 sfree(of);
545 void write_traj(FILE *fplog,t_commrec *cr,
546 gmx_mdoutf_t *of,
547 int mdof_flags,
548 gmx_mtop_t *top_global,
549 gmx_large_int_t step,double t,
550 t_state *state_local,t_state *state_global,
551 rvec *f_local,rvec *f_global,
552 int *n_xtc,rvec **x_xtc)
554 int i,j;
555 gmx_groups_t *groups;
556 rvec *xxtc;
557 rvec *local_v;
558 rvec *global_v;
560 #define MX(xvf) moveit(cr,GMX_LEFT,GMX_RIGHT,#xvf,xvf)
562 /* MRS -- defining these variables is to manage the difference
563 * between half step and full step velocities, but there must be a better way . . . */
565 local_v = state_local->v;
566 global_v = state_global->v;
568 if (DOMAINDECOMP(cr))
570 if (mdof_flags & MDOF_CPT)
572 dd_collect_state(cr->dd,state_local,state_global);
574 else
576 if (mdof_flags & (MDOF_X | MDOF_XTC))
578 dd_collect_vec(cr->dd,state_local,state_local->x,
579 state_global->x);
581 if (mdof_flags & MDOF_V)
583 dd_collect_vec(cr->dd,state_local,local_v,
584 global_v);
587 if (mdof_flags & MDOF_F)
589 dd_collect_vec(cr->dd,state_local,f_local,f_global);
592 else
594 if (mdof_flags & MDOF_CPT)
596 /* All pointers in state_local are equal to state_global,
597 * but we need to copy the non-pointer entries.
599 state_global->lambda = state_local->lambda;
600 state_global->veta = state_local->veta;
601 state_global->vol0 = state_local->vol0;
602 copy_mat(state_local->box,state_global->box);
603 copy_mat(state_local->boxv,state_global->boxv);
604 copy_mat(state_local->svir_prev,state_global->svir_prev);
605 copy_mat(state_local->fvir_prev,state_global->fvir_prev);
606 copy_mat(state_local->pres_prev,state_global->pres_prev);
608 if (cr->nnodes > 1)
610 /* Particle decomposition, collect the data on the master node */
611 if (mdof_flags & MDOF_CPT)
613 if (state_local->flags & (1<<estX)) MX(state_global->x);
614 if (state_local->flags & (1<<estV)) MX(state_global->v);
615 if (state_local->flags & (1<<estSDX)) MX(state_global->sd_X);
616 if (state_global->nrngi > 1) {
617 if (state_local->flags & (1<<estLD_RNG)) {
618 #ifdef GMX_MPI
619 MPI_Gather(state_local->ld_rng ,
620 state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
621 state_global->ld_rng,
622 state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
623 MASTERRANK(cr),cr->mpi_comm_mygroup);
624 #endif
626 if (state_local->flags & (1<<estLD_RNGI))
628 #ifdef GMX_MPI
629 MPI_Gather(state_local->ld_rngi,
630 sizeof(state_local->ld_rngi[0]),MPI_BYTE,
631 state_global->ld_rngi,
632 sizeof(state_local->ld_rngi[0]),MPI_BYTE,
633 MASTERRANK(cr),cr->mpi_comm_mygroup);
634 #endif
638 else
640 if (mdof_flags & (MDOF_X | MDOF_XTC)) MX(state_global->x);
641 if (mdof_flags & MDOF_V) MX(global_v);
643 if (mdof_flags & MDOF_F) MX(f_global);
647 if (MASTER(cr))
649 if (mdof_flags & MDOF_CPT)
651 write_checkpoint(of->fn_cpt,of->bKeepAndNumCPT,
652 fplog,cr,of->eIntegrator,
653 of->simulation_part,step,t,state_global);
656 if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
658 fwrite_trn(of->fp_trn,step,t,state_local->lambda,
659 state_local->box,top_global->natoms,
660 (mdof_flags & MDOF_X) ? state_global->x : NULL,
661 (mdof_flags & MDOF_V) ? global_v : NULL,
662 (mdof_flags & MDOF_F) ? f_global : NULL);
663 if (gmx_fio_flush(of->fp_trn) != 0)
665 gmx_file("Cannot write trajectory; maybe you are out of quota?");
667 gmx_fio_check_file_position(of->fp_trn);
669 if (mdof_flags & MDOF_XTC) {
670 groups = &top_global->groups;
671 if (*n_xtc == -1)
673 *n_xtc = 0;
674 for(i=0; (i<top_global->natoms); i++)
676 if (ggrpnr(groups,egcXTC,i) == 0)
678 (*n_xtc)++;
681 if (*n_xtc != top_global->natoms)
683 snew(*x_xtc,*n_xtc);
686 if (*n_xtc == top_global->natoms)
688 xxtc = state_global->x;
690 else
692 xxtc = *x_xtc;
693 j = 0;
694 for(i=0; (i<top_global->natoms); i++)
696 if (ggrpnr(groups,egcXTC,i) == 0)
698 copy_rvec(state_global->x[i],xxtc[j++]);
702 if (write_xtc(of->fp_xtc,*n_xtc,step,t,
703 state_local->box,xxtc,of->xtc_prec) == 0)
705 gmx_fatal(FARGS,"XTC error - maybe you are out of quota?");
707 gmx_fio_check_file_position(of->fp_xtc);