1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * GROningen Mixture of Alchemy and Childrens' Stories
44 #include "gmx_fatal.h"
67 #include "checkpoint.h"
69 #include "md_support.h"
73 typedef struct gmx_global_stat
80 gmx_global_stat_t global_stat_init(t_inputrec *ir)
87 snew(gs->itc0,ir->opts.ngtc);
88 snew(gs->itc1,ir->opts.ngtc);
93 void global_stat_destroy(gmx_global_stat_t gs)
101 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
102 gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner) {
107 for (i=0;i<F_NRE;i++)
123 ato[to++] = afrom[from++];
130 ato[to++] = afrom[from++];
136 ato[to++] = afrom[from++];
145 void global_stat(FILE *fplog,gmx_global_stat_t gs,
146 t_commrec *cr,gmx_enerdata_t *enerd,
147 tensor fvir,tensor svir,rvec mu_tot,
148 t_inputrec *inputrec,
149 gmx_ekindata_t *ekind,gmx_constr_t constr,
152 gmx_mtop_t *top_global, t_state *state_local,
153 gmx_bool bSumEkinhOld, int flags)
154 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
158 int ie=0,ifv=0,isv=0,irmsd=0,imu=0;
159 int idedl=0,idvdll=0,idvdlnl=0,iepl=0,icm=0,imass=0,ica=0,inb=0;
161 int icj=-1,ici=-1,icx=-1;
163 real copyenerd[F_NRE];
165 real *rmsd_data=NULL;
167 gmx_bool bVV,bTemp,bEner,bPres,bConstrVir,bEkinAveVel,bFirstIterate,bReadEkin;
169 bVV = EI_VV(inputrec->eI);
170 bTemp = flags & CGLO_TEMPERATURE;
171 bEner = flags & CGLO_ENERGY;
172 bPres = (flags & CGLO_PRESSURE);
173 bConstrVir = (flags & CGLO_CONSTRAINT);
174 bFirstIterate = (flags & CGLO_FIRSTITERATE);
175 bEkinAveVel = (inputrec->eI==eiVV || (inputrec->eI==eiVVAK && bPres));
176 bReadEkin = (flags & CGLO_READEKIN);
184 /* This routine copies all the data to be summed to one big buffer
185 * using the t_bin struct.
188 /* First, we neeed to identify which enerd->term should be
189 communicated. Temperature and pressure terms should only be
190 communicated and summed when they need to be, to avoid repeating
191 the sums and overcounting. */
193 nener = filter_enerdterm(enerd->term,TRUE,copyenerd,bTemp,bPres,bEner);
195 /* First, the data that needs to be communicated with velocity verlet every time
196 This is just the constraint virial.*/
198 isv = add_binr(rb,DIM*DIM,svir[0]);
202 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
207 for(j=0; (j<inputrec->opts.ngtc); j++)
211 itc0[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
213 if (bEkinAveVel && !bReadEkin)
215 itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinf[0]);
219 itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh[0]);
222 /* these probably need to be put into one of these categories */
224 idedl = add_binr(rb,1,&(ekind->dekindl));
226 ica = add_binr(rb,1,&(ekind->cosacc.mvcos));
232 if ((bPres || !bVV) && bFirstIterate)
234 ifv = add_binr(rb,DIM*DIM,fvir[0]);
243 ie = add_binr(rb,nener,copyenerd);
248 rmsd_data = constr_rmsd_data(constr);
251 irmsd = add_binr(rb,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
254 if (!NEED_MUTOT(*inputrec))
256 imu = add_binr(rb,DIM,mu_tot);
262 for(j=0; (j<egNR); j++)
264 inn[j]=add_binr(rb,enerd->grpp.nener,enerd->grpp.ener[j]);
267 if (inputrec->efep != efepNO)
269 idvdll = add_bind(rb,efptNR,enerd->dvdl_lin);
270 idvdlnl = add_bind(rb,efptNR,enerd->dvdl_nonlin);
271 if (enerd->n_lambda > 0)
273 iepl = add_bind(rb,enerd->n_lambda,enerd->enerpart_lambda);
281 icm = add_binr(rb,DIM*vcm->nr,vcm->group_p[0]);
283 imass = add_binr(rb,vcm->nr,vcm->group_mass);
285 if (vcm->mode == ecmANGULAR)
287 icj = add_binr(rb,DIM*vcm->nr,vcm->group_j[0]);
289 icx = add_binr(rb,DIM*vcm->nr,vcm->group_x[0]);
291 ici = add_binr(rb,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
296 if (DOMAINDECOMP(cr))
298 nb = cr->dd->nbonded_local;
299 inb = add_bind(rb,1,&nb);
304 isig = add_binr(rb,nsig,sig);
307 /* Global sum it all */
310 fprintf(debug,"Summing %d energies\n",rb->maxreal);
315 /* Extract all the data locally */
319 extract_binr(rb,isv ,DIM*DIM,svir[0]);
322 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
327 for(j=0; (j<inputrec->opts.ngtc); j++)
331 extract_binr(rb,itc0[j],DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
333 if (bEkinAveVel && !bReadEkin) {
334 extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinf[0]);
338 extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinh[0]);
341 extract_binr(rb,idedl,1,&(ekind->dekindl));
342 extract_binr(rb,ica,1,&(ekind->cosacc.mvcos));
346 if ((bPres || !bVV) && bFirstIterate)
348 extract_binr(rb,ifv ,DIM*DIM,fvir[0]);
355 extract_binr(rb,ie,nener,copyenerd);
358 extract_binr(rb,irmsd,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
360 if (!NEED_MUTOT(*inputrec))
362 extract_binr(rb,imu,DIM,mu_tot);
365 for(j=0; (j<egNR); j++)
367 extract_binr(rb,inn[j],enerd->grpp.nener,enerd->grpp.ener[j]);
369 if (inputrec->efep != efepNO)
371 extract_bind(rb,idvdll ,efptNR,enerd->dvdl_lin);
372 extract_bind(rb,idvdlnl,efptNR,enerd->dvdl_nonlin);
373 if (enerd->n_lambda > 0)
375 extract_bind(rb,iepl,enerd->n_lambda,enerd->enerpart_lambda);
378 if (DOMAINDECOMP(cr))
380 extract_bind(rb,inb,1,&nb);
381 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
383 dd_print_missing_interactions(fplog,cr,(int)(nb + 0.5),top_global,state_local);
388 filter_enerdterm(copyenerd,FALSE,enerd->term,bTemp,bPres,bEner);
394 extract_binr(rb,icm,DIM*vcm->nr,vcm->group_p[0]);
396 extract_binr(rb,imass,vcm->nr,vcm->group_mass);
398 if (vcm->mode == ecmANGULAR)
400 extract_binr(rb,icj,DIM*vcm->nr,vcm->group_j[0]);
402 extract_binr(rb,icx,DIM*vcm->nr,vcm->group_x[0]);
404 extract_binr(rb,ici,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
411 extract_binr(rb,isig,nsig,sig);
416 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep)
419 return ((step % nstep)==0);
424 static void moveit(t_commrec *cr,
425 int left,int right,const char *s,rvec xx[])
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)
440 gmx_bool bAppendFiles;
450 of->eIntegrator = ir->eI;
451 of->bExpanded = ir->bExpanded;
452 of->elamstats = ir->expandedvals->elamstats;
453 of->simulation_part = ir->simulation_part;
457 bAppendFiles = (mdrun_flags & MD_APPENDFILES);
459 of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT);
461 sprintf(filemode, bAppendFiles ? "a+" : "w+");
463 if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
466 !(EI_DYNAMICS(ir->eI) &&
473 of->fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
475 if (EI_DYNAMICS(ir->eI) &&
478 of->fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
479 of->xtc_prec = ir->xtcprec;
481 if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
483 of->fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
485 of->fn_cpt = opt2fn("-cpo",nfile,fnm);
487 if ((ir->efep != efepNO || ir->bSimTemp) && ir->fepvals->nstdhdl > 0 &&
488 (ir->fepvals->separate_dhdl_file == esepdhdlfileYES ) &&
493 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
497 of->fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
501 if (opt2bSet("-field",nfile,fnm) &&
502 (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
506 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field",nfile,fnm),
511 of->fp_field = xvgropen(opt2fn("-field",nfile,fnm),
512 "Applied electric field","Time (ps)",
521 void done_mdoutf(gmx_mdoutf_t *of)
523 if (of->fp_ene != NULL)
525 close_enx(of->fp_ene);
529 close_xtc(of->fp_xtc);
533 close_trn(of->fp_trn);
535 if (of->fp_dhdl != NULL)
537 gmx_fio_fclose(of->fp_dhdl);
539 if (of->fp_field != NULL)
541 gmx_fio_fclose(of->fp_field);
547 void write_traj(FILE *fplog,t_commrec *cr,
550 gmx_mtop_t *top_global,
551 gmx_large_int_t step,double t,
552 t_state *state_local,t_state *state_global,
553 rvec *f_local,rvec *f_global,
554 int *n_xtc,rvec **x_xtc)
557 gmx_groups_t *groups;
562 #define MX(xvf) moveit(cr,GMX_LEFT,GMX_RIGHT,#xvf,xvf)
564 /* MRS -- defining these variables is to manage the difference
565 * between half step and full step velocities, but there must be a better way . . . */
567 local_v = state_local->v;
568 global_v = state_global->v;
570 if (DOMAINDECOMP(cr))
572 if (mdof_flags & MDOF_CPT)
574 dd_collect_state(cr->dd,state_local,state_global);
578 if (mdof_flags & (MDOF_X | MDOF_XTC))
580 dd_collect_vec(cr->dd,state_local,state_local->x,
583 if (mdof_flags & MDOF_V)
585 dd_collect_vec(cr->dd,state_local,local_v,
589 if (mdof_flags & MDOF_F)
591 dd_collect_vec(cr->dd,state_local,f_local,f_global);
596 if (mdof_flags & MDOF_CPT)
598 /* All pointers in state_local are equal to state_global,
599 * but we need to copy the non-pointer entries.
601 state_global->lambda = state_local->lambda;
602 state_global->veta = state_local->veta;
603 state_global->vol0 = state_local->vol0;
604 copy_mat(state_local->box,state_global->box);
605 copy_mat(state_local->boxv,state_global->boxv);
606 copy_mat(state_local->svir_prev,state_global->svir_prev);
607 copy_mat(state_local->fvir_prev,state_global->fvir_prev);
608 copy_mat(state_local->pres_prev,state_global->pres_prev);
612 /* Particle decomposition, collect the data on the master node */
613 if (mdof_flags & MDOF_CPT)
615 if (state_local->flags & (1<<estX)) MX(state_global->x);
616 if (state_local->flags & (1<<estV)) MX(state_global->v);
617 if (state_local->flags & (1<<estSDX)) MX(state_global->sd_X);
618 if (state_global->nrngi > 1) {
619 if (state_local->flags & (1<<estLD_RNG)) {
621 MPI_Gather(state_local->ld_rng ,
622 state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
623 state_global->ld_rng,
624 state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
625 MASTERRANK(cr),cr->mpi_comm_mygroup);
628 if (state_local->flags & (1<<estLD_RNGI))
631 MPI_Gather(state_local->ld_rngi,
632 sizeof(state_local->ld_rngi[0]),MPI_BYTE,
633 state_global->ld_rngi,
634 sizeof(state_local->ld_rngi[0]),MPI_BYTE,
635 MASTERRANK(cr),cr->mpi_comm_mygroup);
642 if (mdof_flags & (MDOF_X | MDOF_XTC)) MX(state_global->x);
643 if (mdof_flags & MDOF_V) MX(global_v);
645 if (mdof_flags & MDOF_F) MX(f_global);
651 if (mdof_flags & MDOF_CPT)
653 write_checkpoint(of->fn_cpt,of->bKeepAndNumCPT,
654 fplog,cr,of->eIntegrator,of->simulation_part,
655 of->bExpanded,of->elamstats,step,t,state_global);
658 if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
660 fwrite_trn(of->fp_trn,step,t,state_local->lambda[efptFEP],
661 state_local->box,top_global->natoms,
662 (mdof_flags & MDOF_X) ? state_global->x : NULL,
663 (mdof_flags & MDOF_V) ? global_v : NULL,
664 (mdof_flags & MDOF_F) ? f_global : NULL);
665 if (gmx_fio_flush(of->fp_trn) != 0)
667 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
669 gmx_fio_check_file_position(of->fp_trn);
671 if (mdof_flags & MDOF_XTC) {
672 groups = &top_global->groups;
676 for(i=0; (i<top_global->natoms); i++)
678 if (ggrpnr(groups,egcXTC,i) == 0)
683 if (*n_xtc != top_global->natoms)
688 if (*n_xtc == top_global->natoms)
690 xxtc = state_global->x;
696 for(i=0; (i<top_global->natoms); i++)
698 if (ggrpnr(groups,egcXTC,i) == 0)
700 copy_rvec(state_global->x[i],xxtc[j++]);
704 if (write_xtc(of->fp_xtc,*n_xtc,step,t,
705 state_local->box,xxtc,of->xtc_prec) == 0)
707 gmx_fatal(FARGS,"XTC error - maybe you are out of disk space?");
709 gmx_fio_check_file_position(of->fp_xtc);