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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
57 #include "mtop_util.h"
61 #include "mdebin_bar.h"
64 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
66 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
68 static const char *tricl_boxs_nm[] = {
69 "Box-XX", "Box-YY", "Box-ZZ",
70 "Box-YX", "Box-ZX", "Box-ZY"
73 static const char *vol_nm[] = { "Volume" };
75 static const char *dens_nm[] = {"Density" };
77 static const char *pv_nm[] = {"pV" };
79 static const char *enthalpy_nm[] = {"Enthalpy" };
81 static const char *boxvel_nm[] = {
82 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
83 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
86 #define NBOXS asize(boxs_nm)
87 #define NTRICLBOXS asize(tricl_boxs_nm)
89 t_mdebin *init_mdebin(ener_file_t fp_ene,
90 const gmx_mtop_t *mtop,
94 const char *ener_nm[F_NRE];
95 static const char *vir_nm[] = {
96 "Vir-XX", "Vir-XY", "Vir-XZ",
97 "Vir-YX", "Vir-YY", "Vir-YZ",
98 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
100 static const char *sv_nm[] = {
101 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
102 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
103 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
105 static const char *fv_nm[] = {
106 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
107 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
108 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
110 static const char *pres_nm[] = {
111 "Pres-XX","Pres-XY","Pres-XZ",
112 "Pres-YX","Pres-YY","Pres-YZ",
113 "Pres-ZX","Pres-ZY","Pres-ZZ"
115 static const char *surft_nm[] = {
118 static const char *mu_nm[] = {
119 "Mu-X", "Mu-Y", "Mu-Z"
121 static const char *vcos_nm[] = {
124 static const char *visc_nm[] = {
127 static const char *baro_nm[] = {
132 const gmx_groups_t *groups;
137 int i,j,ni,nj,n,nh,k,kk,ncon,nset;
138 gmx_bool bBHAM,bNoseHoover,b14;
147 if (EI_DYNAMICS(ir->eI))
149 md->delta_t = ir->delta_t;
156 groups = &mtop->groups;
158 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
159 b14 = (gmx_mtop_ftype_count(mtop,F_LJ14) > 0 ||
160 gmx_mtop_ftype_count(mtop,F_LJC14_Q) > 0);
162 ncon = gmx_mtop_ftype_count(mtop,F_CONSTR);
163 nset = gmx_mtop_ftype_count(mtop,F_SETTLE);
164 md->bConstr = (ncon > 0 || nset > 0);
165 md->bConstrVir = FALSE;
167 if (ncon > 0 && ir->eConstrAlg == econtLINCS) {
173 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
178 /* Energy monitoring */
185 for(i=0; i<F_NRE; i++)
187 md->bEner[i] = FALSE;
189 md->bEner[i] = !bBHAM;
190 else if (i == F_BHAM)
191 md->bEner[i] = bBHAM;
193 md->bEner[i] = ir->bQMMM;
194 else if (i == F_COUL_LR)
195 md->bEner[i] = (ir->rcoulomb > ir->rlist);
196 else if (i == F_LJ_LR)
197 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
198 else if (i == F_BHAM_LR)
199 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
200 else if (i == F_RF_EXCL)
201 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
202 else if (i == F_COUL_RECIP)
203 md->bEner[i] = EEL_FULL(ir->coulombtype);
204 else if (i == F_LJ14)
206 else if (i == F_COUL14)
208 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
209 md->bEner[i] = FALSE;
210 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
211 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
212 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
213 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
214 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
215 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
216 md->bEner[i] = (ir->efep != efepNO);
217 else if ((interaction_function[i].flags & IF_VSITE) ||
218 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
219 md->bEner[i] = FALSE;
220 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i==F_EQM))
222 else if ((i == F_GBPOL) && ir->implicit_solvent==eisGBSA)
224 else if ((i == F_NPSOLVATION) && ir->implicit_solvent==eisGBSA && (ir->sa_algorithm != esaNO))
226 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
227 md->bEner[i] = FALSE;
228 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
229 md->bEner[i] = EI_DYNAMICS(ir->eI);
231 md->bEner[i] = (EI_DYNAMICS(ir->eI) && getenv("GMX_VIRIAL_TEMPERATURE"));
232 else if (i == F_DISPCORR || i == F_PDISPCORR)
233 md->bEner[i] = (ir->eDispCorr != edispcNO);
234 else if (i == F_DISRESVIOL)
235 md->bEner[i] = (gmx_mtop_ftype_count(mtop,F_DISRES) > 0);
236 else if (i == F_ORIRESDEV)
237 md->bEner[i] = (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0);
238 else if (i == F_CONNBONDS)
239 md->bEner[i] = FALSE;
240 else if (i == F_COM_PULL)
241 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
242 else if (i == F_ECONSERVED)
243 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
244 (ir->epc == epcNO || ir->epc==epcMTTK));
246 md->bEner[i] = (gmx_mtop_ftype_count(mtop,i) > 0);
249 /* OpenMM always produces only the following 4 energy terms */
250 md->bEner[F_EPOT] = TRUE;
251 md->bEner[F_EKIN] = TRUE;
252 md->bEner[F_ETOT] = TRUE;
253 md->bEner[F_TEMP] = TRUE;
256 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
257 if (ir->bAdress && !debug) {
258 for (i = 0; i < F_NRE; i++) {
259 md->bEner[i] = FALSE;
260 if(i == F_EKIN){ md->bEner[i] = TRUE;}
261 if(i == F_TEMP){ md->bEner[i] = TRUE;}
270 for(i=0; i<F_NRE; i++)
274 ener_nm[md->f_nre]=interaction_function[i].longname;
280 md->bDiagPres = !TRICLINIC(ir->ref_p);
281 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
282 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
283 md->bDynBox = DYNAMIC_BOX(*ir);
285 md->bNHC_trotter = IR_NVT_TROTTER(ir);
286 md->bPrintNHChains = ir-> bPrintNHChains;
287 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
288 md->bMu = NEED_MUTOT(*ir);
290 md->ebin = mk_ebin();
291 /* Pass NULL for unit to let get_ebin_space determine the units
292 * for interaction_function[i].longname
294 md->ie = get_ebin_space(md->ebin,md->f_nre,ener_nm,NULL);
297 /* This should be called directly after the call for md->ie,
298 * such that md->iconrmsd follows directly in the list.
300 md->iconrmsd = get_ebin_space(md->ebin,md->nCrmsd,conrmsd_nm,"");
304 md->ib = get_ebin_space(md->ebin,
305 md->bTricl ? NTRICLBOXS : NBOXS,
306 md->bTricl ? tricl_boxs_nm : boxs_nm,
308 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
309 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
312 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
313 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
318 md->isvir = get_ebin_space(md->ebin,asize(sv_nm),sv_nm,unit_energy);
319 md->ifvir = get_ebin_space(md->ebin,asize(fv_nm),fv_nm,unit_energy);
322 md->ivir = get_ebin_space(md->ebin,asize(vir_nm),vir_nm,unit_energy);
324 md->ipres = get_ebin_space(md->ebin,asize(pres_nm),pres_nm,unit_pres_bar);
326 md->isurft = get_ebin_space(md->ebin,asize(surft_nm),surft_nm,
328 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
330 md->ipc = get_ebin_space(md->ebin,md->bTricl ? 6 : 3,
335 md->imu = get_ebin_space(md->ebin,asize(mu_nm),mu_nm,unit_dipole_D);
337 if (ir->cos_accel != 0)
339 md->ivcos = get_ebin_space(md->ebin,asize(vcos_nm),vcos_nm,unit_vel);
340 md->ivisc = get_ebin_space(md->ebin,asize(visc_nm),visc_nm,
344 /* Energy monitoring */
347 md->bEInd[i] = FALSE;
349 md->bEInd[egCOULSR] = TRUE;
350 md->bEInd[egLJSR ] = TRUE;
352 if (ir->rcoulomb > ir->rlist)
354 md->bEInd[egCOULLR] = TRUE;
358 if (ir->rvdw > ir->rlist)
360 md->bEInd[egLJLR] = TRUE;
365 md->bEInd[egLJSR] = FALSE;
366 md->bEInd[egBHAMSR] = TRUE;
367 if (ir->rvdw > ir->rlist)
369 md->bEInd[egBHAMLR] = TRUE;
374 md->bEInd[egLJ14] = TRUE;
375 md->bEInd[egCOUL14] = TRUE;
378 for(i=0; (i<egNR); i++)
386 n=groups->grps[egcENER].nr;
387 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
389 /*standard simulation*/
394 /*AdResS simulation*/
401 snew(md->igrp,md->nE);
406 for(k=0; (k<md->nEc); k++)
410 for(i=0; (i<groups->grps[egcENER].nr); i++)
412 ni=groups->grps[egcENER].nm_ind[i];
413 for(j=i; (j<groups->grps[egcENER].nr); j++)
415 nj=groups->grps[egcENER].nm_ind[j];
416 for(k=kk=0; (k<egNR); k++)
420 sprintf(gnm[kk],"%s:%s-%s",egrp_nm[k],
421 *(groups->grpname[ni]),*(groups->grpname[nj]));
425 md->igrp[n]=get_ebin_space(md->ebin,md->nEc,
426 (const char **)gnm,unit_energy);
430 for(k=0; (k<md->nEc); k++)
438 gmx_incons("Number of energy terms wrong");
442 md->nTC=groups->grps[egcTC].nr;
443 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
446 md->nTCP = 1; /* assume only one possible coupling system for barostat
453 if (md->etc == etcNOSEHOOVER)
455 if (md->bNHC_trotter)
457 md->mde_n = 2*md->nNHC*md->nTC;
461 md->mde_n = 2*md->nTC;
463 if (md->epc == epcMTTK)
465 md->mdeb_n = 2*md->nNHC*md->nTCP;
472 snew(md->tmp_r,md->mde_n);
473 snew(md->tmp_v,md->mde_n);
474 snew(md->grpnms,md->mde_n);
477 for(i=0; (i<md->nTC); i++)
479 ni=groups->grps[egcTC].nm_ind[i];
480 sprintf(buf,"T-%s",*(groups->grpname[ni]));
481 grpnms[i]=strdup(buf);
483 md->itemp=get_ebin_space(md->ebin,md->nTC,(const char **)grpnms,
486 if (md->etc == etcNOSEHOOVER)
488 if (md->bPrintNHChains)
490 if (md->bNHC_trotter)
492 for(i=0; (i<md->nTC); i++)
494 ni=groups->grps[egcTC].nm_ind[i];
495 bufi = *(groups->grpname[ni]);
496 for(j=0; (j<md->nNHC); j++)
498 sprintf(buf,"Xi-%d-%s",j,bufi);
499 grpnms[2*(i*md->nNHC+j)]=strdup(buf);
500 sprintf(buf,"vXi-%d-%s",j,bufi);
501 grpnms[2*(i*md->nNHC+j)+1]=strdup(buf);
504 md->itc=get_ebin_space(md->ebin,md->mde_n,
505 (const char **)grpnms,unit_invtime);
508 for(i=0; (i<md->nTCP); i++)
510 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
511 for(j=0; (j<md->nNHC); j++)
513 sprintf(buf,"Xi-%d-%s",j,bufi);
514 grpnms[2*(i*md->nNHC+j)]=strdup(buf);
515 sprintf(buf,"vXi-%d-%s",j,bufi);
516 grpnms[2*(i*md->nNHC+j)+1]=strdup(buf);
519 md->itcb=get_ebin_space(md->ebin,md->mdeb_n,
520 (const char **)grpnms,unit_invtime);
525 for(i=0; (i<md->nTC); i++)
527 ni=groups->grps[egcTC].nm_ind[i];
528 bufi = *(groups->grpname[ni]);
529 sprintf(buf,"Xi-%s",bufi);
530 grpnms[2*i]=strdup(buf);
531 sprintf(buf,"vXi-%s",bufi);
532 grpnms[2*i+1]=strdup(buf);
534 md->itc=get_ebin_space(md->ebin,md->mde_n,
535 (const char **)grpnms,unit_invtime);
539 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
540 md->etc == etcVRESCALE)
542 for(i=0; (i<md->nTC); i++)
544 ni=groups->grps[egcTC].nm_ind[i];
545 sprintf(buf,"Lamb-%s",*(groups->grpname[ni]));
546 grpnms[i]=strdup(buf);
548 md->itc=get_ebin_space(md->ebin,md->mde_n,(const char **)grpnms,"");
554 md->nU=groups->grps[egcACC].nr;
557 snew(grpnms,3*md->nU);
558 for(i=0; (i<md->nU); i++)
560 ni=groups->grps[egcACC].nm_ind[i];
561 sprintf(buf,"Ux-%s",*(groups->grpname[ni]));
562 grpnms[3*i+XX]=strdup(buf);
563 sprintf(buf,"Uy-%s",*(groups->grpname[ni]));
564 grpnms[3*i+YY]=strdup(buf);
565 sprintf(buf,"Uz-%s",*(groups->grpname[ni]));
566 grpnms[3*i+ZZ]=strdup(buf);
568 md->iu=get_ebin_space(md->ebin,3*md->nU,(const char **)grpnms,unit_vel);
574 do_enxnms(fp_ene,&md->ebin->nener,&md->ebin->enm);
577 md->print_grpnms=NULL;
579 /* check whether we're going to write dh histograms */
581 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO )
583 /* Currently dh histograms are only written with dynamics */
584 if (EI_DYNAMICS(ir->eI))
588 mde_delta_h_coll_init(md->dhc, ir);
594 md->fp_dhdl = fp_dhdl;
598 snew(md->temperatures,ir->fepvals->n_lambda);
599 for (i=0;i<ir->fepvals->n_lambda;i++)
601 md->temperatures[i] = ir->simtempvals->temperatures[i];
607 extern FILE *open_dhdl(const char *filename,const t_inputrec *ir,
608 const output_env_t oenv)
611 const char *dhdl="dH/d\\lambda",*deltag="\\DeltaH",*lambda="\\lambda",
612 *lambdastate="\\lambda state",*remain="remaining";
613 char title[STRLEN],label_x[STRLEN],label_y[STRLEN];
614 int i,np,nps,nsets,nsets_de,nsetsbegin;
627 if (fep->n_lambda == 0)
629 sprintf(title,"%s",dhdl);
630 sprintf(label_x,"Time (ps)");
631 sprintf(label_y,"%s (%s %s)",
632 dhdl,unit_energy,"[\\lambda]\\S-1\\N");
636 sprintf(title,"%s and %s",dhdl,deltag);
637 sprintf(label_x,"Time (ps)");
638 sprintf(label_y,"%s and %s (%s %s)",
639 dhdl,deltag,unit_energy,"[\\8l\\4]\\S-1\\N");
641 fp = gmx_fio_fopen(filename,"w+");
642 xvgr_header(fp,title,label_x,label_y,exvggtXNY,oenv);
646 bufplace = sprintf(buf,"T = %g (K) ",
649 if (ir->efep != efepSLOWGROWTH)
651 if (fep->n_lambda == 0)
653 sprintf(&(buf[bufplace]),"%s = %g",
654 lambda,fep->init_lambda);
658 sprintf(&(buf[bufplace]),"%s = %d",
659 lambdastate,fep->init_fep_state);
662 xvgr_subtitle(fp,buf,oenv);
664 for (i=0;i<efptNR;i++)
666 if (fep->separate_dvdl[i]) {nsets_dhdl++;}
669 /* count the number of delta_g states */
670 nsets_de = fep->n_lambda;
672 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
674 if (fep->n_lambda>0 && ir->bExpanded)
676 nsets += 1; /*add fep state for expanded ensemble */
679 if (fep->bPrintEnergy)
681 nsets += 1; /* add energy to the dhdl as well */
685 if ((ir->epc!=epcNO) && (fep->n_lambda>0))
687 nsetsextend += 1; /* for PV term, other terms possible if required for the reduced potential (only needed with foreign lambda) */
689 snew(setname,nsetsextend);
693 /* state for the fep_vals, if we have alchemical sampling */
694 sprintf(buf,"%s","Thermodynamic state");
695 setname[s] = strdup(buf);
699 if (fep->bPrintEnergy)
701 sprintf(buf,"%s (%s)","Energy",unit_energy);
702 setname[s] = strdup(buf);
706 for (i=0;i<efptNR;i++)
708 if (fep->separate_dvdl[i]) {
709 sprintf(buf,"%s (%s)",dhdl,efpt_names[i]);
710 setname[s] = strdup(buf);
715 if (fep->n_lambda > 0)
717 /* g_bar has to determine the lambda values used in this simulation
718 * from this xvg legend.
722 nsetsbegin = 1; /* for including the expanded ensemble */
727 if (fep->bPrintEnergy)
731 nsetsbegin += nsets_dhdl;
733 for(s=nsetsbegin; s<nsets; s++)
735 nps = sprintf(buf,"%s %s (",deltag,lambda);
736 for (i=0;i<efptNR;i++)
738 if (fep->separate_dvdl[i])
740 np = sprintf(&buf[nps],"%g,",fep->all_lambda[i][s-(nsetsbegin)]);
746 /* print the temperature for this state if doing simulated annealing */
747 sprintf(&buf[nps],"T = %g (%s))",ir->simtempvals->temperatures[s-(nsetsbegin)],unit_temp_K);
751 sprintf(&buf[nps-1],")"); /* -1 to overwrite the last comma */
753 setname[s] = strdup(buf);
755 if (ir->epc!=epcNO) {
756 np = sprintf(buf,"pV (%s)",unit_energy);
757 setname[nsetsextend-1] = strdup(buf); /* the first entry after nsets */
760 xvgr_legend(fp,nsetsextend,(const char **)setname,oenv);
762 for(s=0; s<nsetsextend; s++)
772 static void copy_energy(t_mdebin *md, real e[],real ecpy[])
776 for(i=j=0; (i<F_NRE); i++)
780 gmx_incons("Number of energy terms wrong");
783 void upd_mdebin(t_mdebin *md,
788 gmx_enerdata_t *enerd,
797 gmx_ekindata_t *ekind,
801 int i,j,k,kk,m,n,gid;
802 real crmsd[2],tmp6[6];
803 real bs[NTRICLBOXS],vol,dens,pv,enthalpy;
806 double store_dhdl[efptNR];
811 /* Do NOT use the box in the state variable, but the separate box provided
812 * as an argument. This is because we sometimes need to write the box from
813 * the last timestep to match the trajectory frames.
815 copy_energy(md, enerd->term,ecopy);
816 add_ebin(md->ebin,md->ie,md->f_nre,ecopy,bSum);
819 crmsd[0] = constr_rmsd(constr,FALSE);
822 crmsd[1] = constr_rmsd(constr,TRUE);
824 add_ebin(md->ebin,md->iconrmsd,md->nCrmsd,crmsd,FALSE);
846 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
847 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
848 add_ebin(md->ebin,md->ib ,nboxs,bs ,bSum);
849 add_ebin(md->ebin,md->ivol ,1 ,&vol ,bSum);
850 add_ebin(md->ebin,md->idens,1 ,&dens,bSum);
854 /* This is pV (in kJ/mol). The pressure is the reference pressure,
855 not the instantaneous pressure */
856 pv = vol*md->ref_p/PRESFAC;
858 add_ebin(md->ebin,md->ipv ,1 ,&pv ,bSum);
859 enthalpy = pv + enerd->term[F_ETOT];
860 add_ebin(md->ebin,md->ienthalpy ,1 ,&enthalpy ,bSum);
865 add_ebin(md->ebin,md->isvir,9,svir[0],bSum);
866 add_ebin(md->ebin,md->ifvir,9,fvir[0],bSum);
869 add_ebin(md->ebin,md->ivir,9,vir[0],bSum);
871 add_ebin(md->ebin,md->ipres,9,pres[0],bSum);
873 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
874 add_ebin(md->ebin,md->isurft,1,&tmp,bSum);
876 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
878 tmp6[0] = state->boxv[XX][XX];
879 tmp6[1] = state->boxv[YY][YY];
880 tmp6[2] = state->boxv[ZZ][ZZ];
881 tmp6[3] = state->boxv[YY][XX];
882 tmp6[4] = state->boxv[ZZ][XX];
883 tmp6[5] = state->boxv[ZZ][YY];
884 add_ebin(md->ebin,md->ipc,md->bTricl ? 6 : 3,tmp6,bSum);
888 add_ebin(md->ebin,md->imu,3,mu_tot,bSum);
890 if (ekind && ekind->cosacc.cos_accel != 0)
892 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
893 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
894 add_ebin(md->ebin,md->ivcos,1,&(ekind->cosacc.vcos),bSum);
895 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
896 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
897 *dens*vol*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
898 add_ebin(md->ebin,md->ivisc,1,&tmp,bSum);
903 for(i=0; (i<md->nEg); i++)
905 for(j=i; (j<md->nEg); j++)
907 gid=GID(i,j,md->nEg);
908 for(k=kk=0; (k<egNR); k++)
912 eee[kk++] = enerd->grpp.ener[k][gid];
915 add_ebin(md->ebin,md->igrp[n],md->nEc,eee,bSum);
923 for(i=0; (i<md->nTC); i++)
925 md->tmp_r[i] = ekind->tcstat[i].T;
927 add_ebin(md->ebin,md->itemp,md->nTC,md->tmp_r,bSum);
929 if (md->etc == etcNOSEHOOVER)
931 /* whether to print Nose-Hoover chains: */
932 if (md->bPrintNHChains)
934 if (md->bNHC_trotter)
936 for(i=0; (i<md->nTC); i++)
938 for (j=0;j<md->nNHC;j++)
941 md->tmp_r[2*k] = state->nosehoover_xi[k];
942 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
945 add_ebin(md->ebin,md->itc,md->mde_n,md->tmp_r,bSum);
948 for(i=0; (i<md->nTCP); i++)
950 for (j=0;j<md->nNHC;j++)
953 md->tmp_r[2*k] = state->nhpres_xi[k];
954 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
957 add_ebin(md->ebin,md->itcb,md->mdeb_n,md->tmp_r,bSum);
962 for(i=0; (i<md->nTC); i++)
964 md->tmp_r[2*i] = state->nosehoover_xi[i];
965 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
967 add_ebin(md->ebin,md->itc,md->mde_n,md->tmp_r,bSum);
971 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
972 md->etc == etcVRESCALE)
974 for(i=0; (i<md->nTC); i++)
976 md->tmp_r[i] = ekind->tcstat[i].lambda;
978 add_ebin(md->ebin,md->itc,md->nTC,md->tmp_r,bSum);
982 if (ekind && md->nU > 1)
984 for(i=0; (i<md->nU); i++)
986 copy_rvec(ekind->grpstat[i].u,md->tmp_v[i]);
988 add_ebin(md->ebin,md->iu,3*md->nU,md->tmp_v[0],bSum);
991 ebin_increase_count(md->ebin,bSum);
993 /* BAR + thermodynamic integration values */
994 if ((md->fp_dhdl || md->dhc) && bDoDHDL && (enerd->n_lambda > 0))
996 snew(dE,enerd->n_lambda-1);
997 for(i=0; i<enerd->n_lambda-1; i++) {
998 dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0]; /* zero for simulated tempering */
999 if (md->temperatures!=NULL)
1001 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1002 /* is this even useful to have at all? */
1003 dE[i] += (md->temperatures[i]/md->temperatures[state->fep_state]-1.0)*enerd->term[F_EKIN];
1008 if (md->fp_dhdl && bDoDHDL)
1010 fprintf(md->fp_dhdl,"%.4f",time);
1011 /* the current free energy state */
1013 /* print the current state if we are doing expanded ensemble */
1014 if (expand->elmcmove > elmcmoveNO) {
1015 fprintf(md->fp_dhdl," %4d",state->fep_state);
1017 /* total energy (for if the temperature changes */
1018 if (fep->bPrintEnergy)
1020 store_energy = enerd->term[F_ETOT];
1021 fprintf(md->fp_dhdl," %#.8g",store_energy);
1024 for (i=0;i<efptNR;i++)
1026 if (fep->separate_dvdl[i])
1028 fprintf(md->fp_dhdl," %#.8g",enerd->term[F_DVDL+i]); /* assumes F_DVDL is first */
1031 for(i=1; i<enerd->n_lambda; i++)
1033 fprintf(md->fp_dhdl," %#.8g",dE[i-1]);
1036 if ((md->epc!=epcNO) && (enerd->n_lambda > 0))
1038 fprintf(md->fp_dhdl," %#.8g",pv); /* PV term only needed when there are alternate state lambda */
1040 fprintf(md->fp_dhdl,"\n");
1041 /* and the binary free energy output */
1043 if (md->dhc && bDoDHDL)
1046 for (i=0;i<efptNR;i++)
1048 if (fep->separate_dvdl[i])
1050 store_dhdl[idhdl] = enerd->term[F_DVDL+i]; /* assumes F_DVDL is first */
1054 /* store_dh is dE */
1055 mde_delta_h_coll_add_dh(md->dhc,
1056 (double)state->fep_state,
1059 (expand->elamstats>elamstatsNO),
1060 (fep->bPrintEnergy),
1068 if ((md->fp_dhdl || md->dhc) && bDoDHDL && (enerd->n_lambda >0))
1075 void upd_mdebin_step(t_mdebin *md)
1077 ebin_increase_count(md->ebin,FALSE);
1080 static void npr(FILE *log,int n,char c)
1082 for(; (n>0); n--) fprintf(log,"%c",c);
1085 static void pprint(FILE *log,const char *s,t_mdebin *md)
1089 char buf1[22],buf2[22];
1092 fprintf(log,"\t<====== ");
1094 fprintf(log," ==>\n");
1095 fprintf(log,"\t<==== %s ====>\n",s);
1096 fprintf(log,"\t<== ");
1098 fprintf(log," ======>\n\n");
1100 fprintf(log,"\tStatistics over %s steps using %s frames\n",
1101 gmx_step_str(md->ebin->nsteps_sim,buf1),
1102 gmx_step_str(md->ebin->nsum_sim,buf2));
1106 void print_ebin_header(FILE *log,gmx_large_int_t steps,double time,real lambda)
1110 fprintf(log," %12s %12s %12s\n"
1111 " %12s %12.5f %12.5f\n\n",
1112 "Step","Time","Lambda",gmx_step_str(steps,buf),time,lambda);
1115 void print_ebin(ener_file_t fp_ene,gmx_bool bEne,gmx_bool bDR,gmx_bool bOR,
1117 gmx_large_int_t step,double time,
1118 int mode,gmx_bool bCompact,
1119 t_mdebin *md,t_fcdata *fcd,
1120 gmx_groups_t *groups,t_grpopts *opts)
1122 /*static char **grpnms=NULL;*/
1124 int i,j,n,ni,nj,ndr,nor,b;
1126 real *disre_rm3tav, *disre_rt;
1128 /* these are for the old-style blocks (1 subblock, only reals), because
1129 there can be only one per ID for these */
1134 /* temporary arrays for the lambda values to write out */
1135 double enxlambda_data[2];
1145 fr.nsteps = md->ebin->nsteps;
1146 fr.dt = md->delta_t;
1147 fr.nsum = md->ebin->nsum;
1148 fr.nre = (bEne) ? md->ebin->nener : 0;
1149 fr.ener = md->ebin->e;
1150 ndisre = bDR ? fcd->disres.npair : 0;
1151 disre_rm3tav = fcd->disres.rm3tav;
1152 disre_rt = fcd->disres.rt;
1153 /* Optional additional old-style (real-only) blocks. */
1154 for(i=0; i<enxNR; i++)
1158 if (fcd->orires.nr > 0 && bOR)
1160 diagonalize_orires_tensors(&(fcd->orires));
1161 nr[enxOR] = fcd->orires.nr;
1162 block[enxOR] = fcd->orires.otav;
1164 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1166 block[enxORI] = fcd->orires.oinsl;
1167 id[enxORI] = enxORI;
1168 nr[enxORT] = fcd->orires.nex*12;
1169 block[enxORT] = fcd->orires.eig;
1170 id[enxORT] = enxORT;
1173 /* whether we are going to wrte anything out: */
1174 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1177 /* the old-style blocks go first */
1179 for(i=0; i<enxNR; i++)
1186 add_blocks_enxframe(&fr, fr.nblock);
1187 for(b=0;b<fr.nblock;b++)
1189 add_subblocks_enxblock(&(fr.block[b]), 1);
1190 fr.block[b].id=id[b];
1191 fr.block[b].sub[0].nr = nr[b];
1193 fr.block[b].sub[0].type = xdr_datatype_float;
1194 fr.block[b].sub[0].fval = block[b];
1196 fr.block[b].sub[0].type = xdr_datatype_double;
1197 fr.block[b].sub[0].dval = block[b];
1201 /* check for disre block & fill it. */
1206 add_blocks_enxframe(&fr, fr.nblock);
1208 add_subblocks_enxblock(&(fr.block[db]), 2);
1209 fr.block[db].id=enxDISRE;
1210 fr.block[db].sub[0].nr=ndisre;
1211 fr.block[db].sub[1].nr=ndisre;
1213 fr.block[db].sub[0].type=xdr_datatype_float;
1214 fr.block[db].sub[1].type=xdr_datatype_float;
1215 fr.block[db].sub[0].fval=disre_rt;
1216 fr.block[db].sub[1].fval=disre_rm3tav;
1218 fr.block[db].sub[0].type=xdr_datatype_double;
1219 fr.block[db].sub[1].type=xdr_datatype_double;
1220 fr.block[db].sub[0].dval=disre_rt;
1221 fr.block[db].sub[1].dval=disre_rm3tav;
1224 /* here we can put new-style blocks */
1226 /* Free energy perturbation blocks */
1229 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1232 /* we can now free & reset the data in the blocks */
1235 mde_delta_h_coll_reset(md->dhc);
1238 /* do the actual I/O */
1240 gmx_fio_check_file_position(enx_file_pointer(fp_ene));
1243 /* We have stored the sums, so reset the sum history */
1244 reset_ebin_sums(md->ebin);
1252 pprint(log,"A V E R A G E S",md);
1258 pprint(log,"R M S - F L U C T U A T I O N S",md);
1262 gmx_fatal(FARGS,"Invalid print mode (%d)",mode);
1267 for(i=0;i<opts->ngtc;i++)
1269 if(opts->annealing[i]!=eannNO)
1271 fprintf(log,"Current ref_t for group %s: %8.1f\n",
1272 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1276 if (mode==eprNORMAL && fcd->orires.nr>0)
1278 print_orires_log(log,&(fcd->orires));
1280 fprintf(log," Energies (%s)\n",unit_energy);
1281 pr_ebin(log,md->ebin,md->ie,md->f_nre+md->nCrmsd,5,mode,TRUE);
1288 pr_ebin(log,md->ebin,md->ib, md->bTricl ? NTRICLBOXS : NBOXS,5,
1294 fprintf(log," Constraint Virial (%s)\n",unit_energy);
1295 pr_ebin(log,md->ebin,md->isvir,9,3,mode,FALSE);
1297 fprintf(log," Force Virial (%s)\n",unit_energy);
1298 pr_ebin(log,md->ebin,md->ifvir,9,3,mode,FALSE);
1303 fprintf(log," Total Virial (%s)\n",unit_energy);
1304 pr_ebin(log,md->ebin,md->ivir,9,3,mode,FALSE);
1309 fprintf(log," Pressure (%s)\n",unit_pres_bar);
1310 pr_ebin(log,md->ebin,md->ipres,9,3,mode,FALSE);
1315 fprintf(log," Total Dipole (%s)\n",unit_dipole_D);
1316 pr_ebin(log,md->ebin,md->imu,3,3,mode,FALSE);
1322 if (md->print_grpnms==NULL)
1324 snew(md->print_grpnms,md->nE);
1326 for(i=0; (i<md->nEg); i++)
1328 ni=groups->grps[egcENER].nm_ind[i];
1329 for(j=i; (j<md->nEg); j++)
1331 nj=groups->grps[egcENER].nm_ind[j];
1332 sprintf(buf,"%s-%s",*(groups->grpname[ni]),
1333 *(groups->grpname[nj]));
1334 md->print_grpnms[n++]=strdup(buf);
1338 sprintf(buf,"Epot (%s)",unit_energy);
1339 fprintf(log,"%15s ",buf);
1340 for(i=0; (i<egNR); i++)
1344 fprintf(log,"%12s ",egrp_nm[i]);
1348 for(i=0; (i<md->nE); i++)
1350 fprintf(log,"%15s",md->print_grpnms[i]);
1351 pr_ebin(log,md->ebin,md->igrp[i],md->nEc,md->nEc,mode,
1358 pr_ebin(log,md->ebin,md->itemp,md->nTC,4,mode,TRUE);
1363 fprintf(log,"%15s %12s %12s %12s\n",
1364 "Group","Ux","Uy","Uz");
1365 for(i=0; (i<md->nU); i++)
1367 ni=groups->grps[egcACC].nm_ind[i];
1368 fprintf(log,"%15s",*groups->grpname[ni]);
1369 pr_ebin(log,md->ebin,md->iu+3*i,3,3,mode,FALSE);
1378 void update_energyhistory(energyhistory_t * enerhist,t_mdebin * mdebin)
1382 enerhist->nsteps = mdebin->ebin->nsteps;
1383 enerhist->nsum = mdebin->ebin->nsum;
1384 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1385 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1386 enerhist->nener = mdebin->ebin->nener;
1388 if (mdebin->ebin->nsum > 0)
1390 /* Check if we need to allocate first */
1391 if(enerhist->ener_ave == NULL)
1393 snew(enerhist->ener_ave,enerhist->nener);
1394 snew(enerhist->ener_sum,enerhist->nener);
1397 for(i=0;i<enerhist->nener;i++)
1399 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1400 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1404 if (mdebin->ebin->nsum_sim > 0)
1406 /* Check if we need to allocate first */
1407 if(enerhist->ener_sum_sim == NULL)
1409 snew(enerhist->ener_sum_sim,enerhist->nener);
1412 for(i=0;i<enerhist->nener;i++)
1414 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1419 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1423 void restore_energyhistory_from_state(t_mdebin * mdebin,
1424 energyhistory_t * enerhist)
1428 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1429 mdebin->ebin->nener != enerhist->nener)
1431 gmx_fatal(FARGS,"Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1432 mdebin->ebin->nener,enerhist->nener);
1435 mdebin->ebin->nsteps = enerhist->nsteps;
1436 mdebin->ebin->nsum = enerhist->nsum;
1437 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1438 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1440 for(i=0; i<mdebin->ebin->nener; i++)
1442 mdebin->ebin->e[i].eav =
1443 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1444 mdebin->ebin->e[i].esum =
1445 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1446 mdebin->ebin->e_sim[i].esum =
1447 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1451 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);