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)
108 for (i = 0; i < F_NRE; i++)
125 ato[to++] = afrom[from++];
132 ato[to++] = afrom[from++];
138 ato[to++] = afrom[from++];
147 void global_stat(FILE *fplog, gmx_global_stat_t gs,
148 t_commrec *cr, gmx_enerdata_t *enerd,
149 tensor fvir, tensor svir, rvec mu_tot,
150 t_inputrec *inputrec,
151 gmx_ekindata_t *ekind, gmx_constr_t constr,
154 gmx_mtop_t *top_global, t_state *state_local,
155 gmx_bool bSumEkinhOld, int flags)
156 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
160 int ie = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
161 int idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
163 int icj = -1, ici = -1, icx = -1;
165 real copyenerd[F_NRE];
167 real *rmsd_data = NULL;
169 gmx_bool bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bFirstIterate, bReadEkin;
171 bVV = EI_VV(inputrec->eI);
172 bTemp = flags & CGLO_TEMPERATURE;
173 bEner = flags & CGLO_ENERGY;
174 bPres = (flags & CGLO_PRESSURE);
175 bConstrVir = (flags & CGLO_CONSTRAINT);
176 bFirstIterate = (flags & CGLO_FIRSTITERATE);
177 bEkinAveVel = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
178 bReadEkin = (flags & CGLO_READEKIN);
186 /* This routine copies all the data to be summed to one big buffer
187 * using the t_bin struct.
190 /* First, we neeed to identify which enerd->term should be
191 communicated. Temperature and pressure terms should only be
192 communicated and summed when they need to be, to avoid repeating
193 the sums and overcounting. */
195 nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
197 /* First, the data that needs to be communicated with velocity verlet every time
198 This is just the constraint virial.*/
201 isv = add_binr(rb, DIM*DIM, svir[0]);
205 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
210 for (j = 0; (j < inputrec->opts.ngtc); j++)
214 itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
216 if (bEkinAveVel && !bReadEkin)
218 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
222 itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
225 /* these probably need to be put into one of these categories */
227 idedl = add_binr(rb, 1, &(ekind->dekindl));
229 ica = add_binr(rb, 1, &(ekind->cosacc.mvcos));
235 if ((bPres || !bVV) && bFirstIterate)
237 ifv = add_binr(rb, DIM*DIM, fvir[0]);
246 ie = add_binr(rb, nener, copyenerd);
251 rmsd_data = constr_rmsd_data(constr);
254 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
257 if (!NEED_MUTOT(*inputrec))
259 imu = add_binr(rb, DIM, mu_tot);
265 for (j = 0; (j < egNR); j++)
267 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
270 if (inputrec->efep != efepNO)
272 idvdll = add_bind(rb, efptNR, enerd->dvdl_lin);
273 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
274 if (enerd->n_lambda > 0)
276 iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
284 icm = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
286 imass = add_binr(rb, vcm->nr, vcm->group_mass);
288 if (vcm->mode == ecmANGULAR)
290 icj = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
292 icx = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
294 ici = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
299 if (DOMAINDECOMP(cr))
301 nb = cr->dd->nbonded_local;
302 inb = add_bind(rb, 1, &nb);
307 isig = add_binr(rb, nsig, sig);
310 /* Global sum it all */
313 fprintf(debug, "Summing %d energies\n", rb->maxreal);
318 /* Extract all the data locally */
322 extract_binr(rb, isv, DIM*DIM, svir[0]);
325 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
330 for (j = 0; (j < inputrec->opts.ngtc); j++)
334 extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
336 if (bEkinAveVel && !bReadEkin)
338 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
342 extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
345 extract_binr(rb, idedl, 1, &(ekind->dekindl));
346 extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
350 if ((bPres || !bVV) && bFirstIterate)
352 extract_binr(rb, ifv, DIM*DIM, fvir[0]);
359 extract_binr(rb, ie, nener, copyenerd);
362 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
364 if (!NEED_MUTOT(*inputrec))
366 extract_binr(rb, imu, DIM, mu_tot);
369 for (j = 0; (j < egNR); j++)
371 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
373 if (inputrec->efep != efepNO)
375 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
376 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
377 if (enerd->n_lambda > 0)
379 extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
382 if (DOMAINDECOMP(cr))
384 extract_bind(rb, inb, 1, &nb);
385 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
387 dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
392 filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
398 extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
400 extract_binr(rb, imass, vcm->nr, vcm->group_mass);
402 if (vcm->mode == ecmANGULAR)
404 extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
406 extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
408 extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
415 extract_binr(rb, isig, nsig, sig);
420 int do_per_step(gmx_large_int_t step, gmx_large_int_t nstep)
424 return ((step % nstep) == 0);
432 static void moveit(t_commrec *cr,
433 int left, int right, const char *s, rvec xx[])
440 move_rvecs(cr, FALSE, FALSE, left, right,
441 xx, NULL, (cr->nnodes-cr->npmenodes)-1, NULL);
444 gmx_mdoutf_t *init_mdoutf(int nfile, const t_filenm fnm[], int mdrun_flags,
445 const t_commrec *cr, const t_inputrec *ir,
446 const output_env_t oenv)
450 gmx_bool bAppendFiles;
460 of->eIntegrator = ir->eI;
461 of->bExpanded = ir->bExpanded;
462 of->elamstats = ir->expandedvals->elamstats;
463 of->simulation_part = ir->simulation_part;
467 bAppendFiles = (mdrun_flags & MD_APPENDFILES);
469 of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT);
471 sprintf(filemode, bAppendFiles ? "a+" : "w+");
473 if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
476 !(EI_DYNAMICS(ir->eI) &&
483 of->fp_trn = open_trn(ftp2fn(efTRN, nfile, fnm), filemode);
485 if (EI_DYNAMICS(ir->eI) &&
488 of->fp_xtc = open_xtc(ftp2fn(efXTC, nfile, fnm), filemode);
489 of->xtc_prec = ir->xtcprec;
491 if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
493 of->fp_ene = open_enx(ftp2fn(efEDR, nfile, fnm), filemode);
495 of->fn_cpt = opt2fn("-cpo", nfile, fnm);
497 if ((ir->efep != efepNO || ir->bSimTemp) && ir->fepvals->nstdhdl > 0 &&
498 (ir->fepvals->separate_dhdl_file == esepdhdlfileYES ) &&
503 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl", nfile, fnm), filemode);
507 of->fp_dhdl = open_dhdl(opt2fn("-dhdl", nfile, fnm), ir, oenv);
511 if (opt2bSet("-field", nfile, fnm) &&
512 (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
516 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field", nfile, fnm),
521 of->fp_field = xvgropen(opt2fn("-field", nfile, fnm),
522 "Applied electric field", "Time (ps)",
531 void done_mdoutf(gmx_mdoutf_t *of)
533 if (of->fp_ene != NULL)
535 close_enx(of->fp_ene);
539 close_xtc(of->fp_xtc);
543 close_trn(of->fp_trn);
545 if (of->fp_dhdl != NULL)
547 gmx_fio_fclose(of->fp_dhdl);
549 if (of->fp_field != NULL)
551 gmx_fio_fclose(of->fp_field);
557 void write_traj(FILE *fplog, t_commrec *cr,
560 gmx_mtop_t *top_global,
561 gmx_large_int_t step, double t,
562 t_state *state_local, t_state *state_global,
563 rvec *f_local, rvec *f_global,
564 int *n_xtc, rvec **x_xtc)
567 gmx_groups_t *groups;
572 #define MX(xvf) moveit(cr, GMX_LEFT, GMX_RIGHT,#xvf, xvf)
574 /* MRS -- defining these variables is to manage the difference
575 * between half step and full step velocities, but there must be a better way . . . */
577 local_v = state_local->v;
578 global_v = state_global->v;
580 if (DOMAINDECOMP(cr))
582 if (mdof_flags & MDOF_CPT)
584 dd_collect_state(cr->dd, state_local, state_global);
588 if (mdof_flags & (MDOF_X | MDOF_XTC))
590 dd_collect_vec(cr->dd, state_local, state_local->x,
593 if (mdof_flags & MDOF_V)
595 dd_collect_vec(cr->dd, state_local, local_v,
599 if (mdof_flags & MDOF_F)
601 dd_collect_vec(cr->dd, state_local, f_local, f_global);
606 if (mdof_flags & MDOF_CPT)
608 /* All pointers in state_local are equal to state_global,
609 * but we need to copy the non-pointer entries.
611 state_global->lambda = state_local->lambda;
612 state_global->veta = state_local->veta;
613 state_global->vol0 = state_local->vol0;
614 copy_mat(state_local->box, state_global->box);
615 copy_mat(state_local->boxv, state_global->boxv);
616 copy_mat(state_local->svir_prev, state_global->svir_prev);
617 copy_mat(state_local->fvir_prev, state_global->fvir_prev);
618 copy_mat(state_local->pres_prev, state_global->pres_prev);
622 /* Particle decomposition, collect the data on the master node */
623 if (mdof_flags & MDOF_CPT)
625 if (state_local->flags & (1<<estX))
629 if (state_local->flags & (1<<estV))
633 if (state_local->flags & (1<<estSDX))
635 MX(state_global->sd_X);
637 if (state_global->nrngi > 1)
639 if (state_local->flags & (1<<estLD_RNG))
642 MPI_Gather(state_local->ld_rng,
643 state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
644 state_global->ld_rng,
645 state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
646 MASTERRANK(cr), cr->mpi_comm_mygroup);
649 if (state_local->flags & (1<<estLD_RNGI))
652 MPI_Gather(state_local->ld_rngi,
653 sizeof(state_local->ld_rngi[0]), MPI_BYTE,
654 state_global->ld_rngi,
655 sizeof(state_local->ld_rngi[0]), MPI_BYTE,
656 MASTERRANK(cr), cr->mpi_comm_mygroup);
663 if (mdof_flags & (MDOF_X | MDOF_XTC))
667 if (mdof_flags & MDOF_V)
672 if (mdof_flags & MDOF_F)
681 if (mdof_flags & MDOF_CPT)
683 write_checkpoint(of->fn_cpt, of->bKeepAndNumCPT,
684 fplog, cr, of->eIntegrator, of->simulation_part,
685 of->bExpanded, of->elamstats, step, t, state_global);
688 if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
690 fwrite_trn(of->fp_trn, step, t, state_local->lambda[efptFEP],
691 state_local->box, top_global->natoms,
692 (mdof_flags & MDOF_X) ? state_global->x : NULL,
693 (mdof_flags & MDOF_V) ? global_v : NULL,
694 (mdof_flags & MDOF_F) ? f_global : NULL);
695 if (gmx_fio_flush(of->fp_trn) != 0)
697 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
699 gmx_fio_check_file_position(of->fp_trn);
701 if (mdof_flags & MDOF_XTC)
703 groups = &top_global->groups;
707 for (i = 0; (i < top_global->natoms); i++)
709 if (ggrpnr(groups, egcXTC, i) == 0)
714 if (*n_xtc != top_global->natoms)
716 snew(*x_xtc, *n_xtc);
719 if (*n_xtc == top_global->natoms)
721 xxtc = state_global->x;
727 for (i = 0; (i < top_global->natoms); i++)
729 if (ggrpnr(groups, egcXTC, i) == 0)
731 copy_rvec(state_global->x[i], xxtc[j++]);
735 if (write_xtc(of->fp_xtc, *n_xtc, step, t,
736 state_local->box, xxtc, of->xtc_prec) == 0)
738 gmx_fatal(FARGS, "XTC error - maybe you are out of disk space?");
740 gmx_fio_check_file_position(of->fp_xtc);