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 * GROwing Monsters And Cloning Shrimps
47 #include "gromacs/fileio/enxio.h"
55 #include "mtop_util.h"
57 #include "gromacs/fileio/gmxfio.h"
60 #include "mdebin_bar.h"
63 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
65 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
67 static const char *tricl_boxs_nm[] = {
68 "Box-XX", "Box-YY", "Box-ZZ",
69 "Box-YX", "Box-ZX", "Box-ZY"
72 static const char *vol_nm[] = { "Volume" };
74 static const char *dens_nm[] = {"Density" };
76 static const char *pv_nm[] = {"pV" };
78 static const char *enthalpy_nm[] = {"Enthalpy" };
80 static const char *boxvel_nm[] = {
81 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
82 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
85 #define NBOXS asize(boxs_nm)
86 #define NTRICLBOXS asize(tricl_boxs_nm)
88 t_mdebin *init_mdebin(ener_file_t fp_ene,
89 const gmx_mtop_t *mtop,
93 const char *ener_nm[F_NRE];
94 static const char *vir_nm[] = {
95 "Vir-XX", "Vir-XY", "Vir-XZ",
96 "Vir-YX", "Vir-YY", "Vir-YZ",
97 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
99 static const char *sv_nm[] = {
100 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
101 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
102 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
104 static const char *fv_nm[] = {
105 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
106 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
107 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
109 static const char *pres_nm[] = {
110 "Pres-XX", "Pres-XY", "Pres-XZ",
111 "Pres-YX", "Pres-YY", "Pres-YZ",
112 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
114 static const char *surft_nm[] = {
117 static const char *mu_nm[] = {
118 "Mu-X", "Mu-Y", "Mu-Z"
120 static const char *vcos_nm[] = {
123 static const char *visc_nm[] = {
126 static const char *baro_nm[] = {
131 const gmx_groups_t *groups;
136 int i, j, ni, nj, n, nh, k, kk, ncon, nset;
137 gmx_bool bBHAM, bNoseHoover, b14;
146 if (EI_DYNAMICS(ir->eI))
148 md->delta_t = ir->delta_t;
155 groups = &mtop->groups;
157 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
158 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
159 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
161 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
162 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
163 md->bConstr = (ncon > 0 || nset > 0);
164 md->bConstrVir = FALSE;
167 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
178 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
185 /* Energy monitoring */
186 for (i = 0; i < egNR; i++)
188 md->bEInd[i] = FALSE;
191 for (i = 0; i < F_NRE; i++)
193 md->bEner[i] = FALSE;
196 md->bEner[i] = !bBHAM;
198 else if (i == F_BHAM)
200 md->bEner[i] = bBHAM;
204 md->bEner[i] = ir->bQMMM;
206 else if (i == F_COUL_LR)
208 md->bEner[i] = (ir->rcoulomb > ir->rlist);
210 else if (i == F_LJ_LR)
212 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
214 else if (i == F_BHAM_LR)
216 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
218 else if (i == F_RF_EXCL)
220 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
222 else if (i == F_COUL_RECIP)
224 md->bEner[i] = EEL_FULL(ir->coulombtype);
226 else if (i == F_LJ14)
230 else if (i == F_COUL14)
234 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
236 md->bEner[i] = FALSE;
238 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
239 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
240 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
241 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
242 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
243 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
245 md->bEner[i] = (ir->efep != efepNO);
247 else if ((interaction_function[i].flags & IF_VSITE) ||
248 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
250 md->bEner[i] = FALSE;
252 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
256 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
260 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
264 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
266 md->bEner[i] = FALSE;
268 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
270 md->bEner[i] = EI_DYNAMICS(ir->eI);
272 else if (i == F_DISPCORR || i == F_PDISPCORR)
274 md->bEner[i] = (ir->eDispCorr != edispcNO);
276 else if (i == F_DISRESVIOL)
278 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
280 else if (i == F_ORIRESDEV)
282 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
284 else if (i == F_CONNBONDS)
286 md->bEner[i] = FALSE;
288 else if (i == F_COM_PULL)
290 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
292 else if (i == F_ECONSERVED)
294 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
295 (ir->epc == epcNO || ir->epc == epcMTTK));
299 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
303 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
304 if (ir->bAdress && !debug)
306 for (i = 0; i < F_NRE; i++)
308 md->bEner[i] = FALSE;
325 for (i = 0; i < F_NRE; i++)
329 ener_nm[md->f_nre] = interaction_function[i].longname;
335 md->bDiagPres = !TRICLINIC(ir->ref_p);
336 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
337 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
338 md->bDynBox = DYNAMIC_BOX(*ir);
340 md->bNHC_trotter = IR_NVT_TROTTER(ir);
341 md->bPrintNHChains = ir->bPrintNHChains;
342 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
343 md->bMu = NEED_MUTOT(*ir);
345 md->ebin = mk_ebin();
346 /* Pass NULL for unit to let get_ebin_space determine the units
347 * for interaction_function[i].longname
349 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
352 /* This should be called directly after the call for md->ie,
353 * such that md->iconrmsd follows directly in the list.
355 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
359 md->ib = get_ebin_space(md->ebin,
360 md->bTricl ? NTRICLBOXS : NBOXS,
361 md->bTricl ? tricl_boxs_nm : boxs_nm,
363 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
364 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
367 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
368 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
373 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
374 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
378 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
382 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
386 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
389 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
391 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
392 boxvel_nm, unit_vel);
396 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
398 if (ir->cos_accel != 0)
400 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
401 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
405 /* Energy monitoring */
406 for (i = 0; i < egNR; i++)
408 md->bEInd[i] = FALSE;
410 md->bEInd[egCOULSR] = TRUE;
411 md->bEInd[egLJSR ] = TRUE;
413 if (ir->rcoulomb > ir->rlist)
415 md->bEInd[egCOULLR] = TRUE;
419 if (ir->rvdw > ir->rlist)
421 md->bEInd[egLJLR] = TRUE;
426 md->bEInd[egLJSR] = FALSE;
427 md->bEInd[egBHAMSR] = TRUE;
428 if (ir->rvdw > ir->rlist)
430 md->bEInd[egBHAMLR] = TRUE;
435 md->bEInd[egLJ14] = TRUE;
436 md->bEInd[egCOUL14] = TRUE;
439 for (i = 0; (i < egNR); i++)
447 n = groups->grps[egcENER].nr;
448 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
451 /*standard simulation*/
453 md->nE = (n*(n+1))/2;
457 /*AdResS simulation*/
464 snew(md->igrp, md->nE);
469 for (k = 0; (k < md->nEc); k++)
471 snew(gnm[k], STRLEN);
473 for (i = 0; (i < groups->grps[egcENER].nr); i++)
475 ni = groups->grps[egcENER].nm_ind[i];
476 for (j = i; (j < groups->grps[egcENER].nr); j++)
478 nj = groups->grps[egcENER].nm_ind[j];
479 for (k = kk = 0; (k < egNR); k++)
483 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
484 *(groups->grpname[ni]), *(groups->grpname[nj]));
488 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
489 (const char **)gnm, unit_energy);
493 for (k = 0; (k < md->nEc); k++)
501 gmx_incons("Number of energy terms wrong");
505 md->nTC = groups->grps[egcTC].nr;
506 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
509 md->nTCP = 1; /* assume only one possible coupling system for barostat
516 if (md->etc == etcNOSEHOOVER)
518 if (md->bNHC_trotter)
520 md->mde_n = 2*md->nNHC*md->nTC;
524 md->mde_n = 2*md->nTC;
526 if (md->epc == epcMTTK)
528 md->mdeb_n = 2*md->nNHC*md->nTCP;
537 snew(md->tmp_r, md->mde_n);
538 snew(md->tmp_v, md->mde_n);
539 snew(md->grpnms, md->mde_n);
542 for (i = 0; (i < md->nTC); i++)
544 ni = groups->grps[egcTC].nm_ind[i];
545 sprintf(buf, "T-%s", *(groups->grpname[ni]));
546 grpnms[i] = strdup(buf);
548 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
551 if (md->etc == etcNOSEHOOVER)
553 if (md->bPrintNHChains)
555 if (md->bNHC_trotter)
557 for (i = 0; (i < md->nTC); i++)
559 ni = groups->grps[egcTC].nm_ind[i];
560 bufi = *(groups->grpname[ni]);
561 for (j = 0; (j < md->nNHC); j++)
563 sprintf(buf, "Xi-%d-%s", j, bufi);
564 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
565 sprintf(buf, "vXi-%d-%s", j, bufi);
566 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
569 md->itc = get_ebin_space(md->ebin, md->mde_n,
570 (const char **)grpnms, unit_invtime);
573 for (i = 0; (i < md->nTCP); i++)
575 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
576 for (j = 0; (j < md->nNHC); j++)
578 sprintf(buf, "Xi-%d-%s", j, bufi);
579 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
580 sprintf(buf, "vXi-%d-%s", j, bufi);
581 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
584 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
585 (const char **)grpnms, unit_invtime);
590 for (i = 0; (i < md->nTC); i++)
592 ni = groups->grps[egcTC].nm_ind[i];
593 bufi = *(groups->grpname[ni]);
594 sprintf(buf, "Xi-%s", bufi);
595 grpnms[2*i] = strdup(buf);
596 sprintf(buf, "vXi-%s", bufi);
597 grpnms[2*i+1] = strdup(buf);
599 md->itc = get_ebin_space(md->ebin, md->mde_n,
600 (const char **)grpnms, unit_invtime);
604 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
605 md->etc == etcVRESCALE)
607 for (i = 0; (i < md->nTC); i++)
609 ni = groups->grps[egcTC].nm_ind[i];
610 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
611 grpnms[i] = strdup(buf);
613 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
619 md->nU = groups->grps[egcACC].nr;
622 snew(grpnms, 3*md->nU);
623 for (i = 0; (i < md->nU); i++)
625 ni = groups->grps[egcACC].nm_ind[i];
626 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
627 grpnms[3*i+XX] = strdup(buf);
628 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
629 grpnms[3*i+YY] = strdup(buf);
630 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
631 grpnms[3*i+ZZ] = strdup(buf);
633 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
639 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
642 md->print_grpnms = NULL;
644 /* check whether we're going to write dh histograms */
646 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
648 /* Currently dh histograms are only written with dynamics */
649 if (EI_DYNAMICS(ir->eI))
653 mde_delta_h_coll_init(md->dhc, ir);
656 snew(md->dE, ir->fepvals->n_lambda);
660 md->fp_dhdl = fp_dhdl;
661 snew(md->dE, ir->fepvals->n_lambda);
666 snew(md->temperatures, ir->fepvals->n_lambda);
667 for (i = 0; i < ir->fepvals->n_lambda; i++)
669 md->temperatures[i] = ir->simtempvals->temperatures[i];
675 /* print a lambda vector to a string
676 fep = the inputrec's FEP input data
677 i = the index of the lambda vector
678 get_native_lambda = whether to print the native lambda
679 get_names = whether to print the names rather than the values
680 str = the pre-allocated string buffer to print to. */
681 static void print_lambda_vector(t_lambda *fep, int i,
682 gmx_bool get_native_lambda, gmx_bool get_names,
689 for (j = 0; j < efptNR; j++)
691 if (fep->separate_dvdl[j])
696 str[0] = 0; /* reset the string */
699 str += sprintf(str, "("); /* set the opening parenthesis*/
701 for (j = 0; j < efptNR; j++)
703 if (fep->separate_dvdl[j])
708 if (get_native_lambda && fep->init_lambda >= 0)
710 str += sprintf(str, "%.4f", fep->init_lambda);
714 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
719 str += sprintf(str, "%s", efpt_singular_names[j]);
721 /* print comma for the next item */
724 str += sprintf(str, ", ");
731 /* and add the closing parenthesis */
732 str += sprintf(str, ")");
737 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
738 const output_env_t oenv)
741 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
742 *lambdastate = "\\lambda state", *remain = "remaining";
743 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
744 int i, np, nps, nsets, nsets_de, nsetsbegin;
745 int n_lambda_terms = 0;
746 t_lambda *fep = ir->fepvals; /* for simplicity */
747 t_expanded *expand = ir->expandedvals;
749 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
755 gmx_bool write_pV = FALSE;
757 /* count the number of different lambda terms */
758 for (i = 0; i < efptNR; i++)
760 if (fep->separate_dvdl[i])
766 if (fep->n_lambda == 0)
768 sprintf(title, "%s", dhdl);
769 sprintf(label_x, "Time (ps)");
770 sprintf(label_y, "%s (%s %s)",
771 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
775 sprintf(title, "%s and %s", dhdl, deltag);
776 sprintf(label_x, "Time (ps)");
777 sprintf(label_y, "%s and %s (%s %s)",
778 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
780 fp = gmx_fio_fopen(filename, "w+");
781 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
785 bufplace = sprintf(buf, "T = %g (K) ",
788 if (ir->efep != efepSLOWGROWTH)
790 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
792 /* compatibility output */
793 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
797 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
799 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
801 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
802 lambdastate, fep->init_fep_state,
803 lambda_name_str, lambda_vec_str);
806 xvgr_subtitle(fp, buf, oenv);
810 if (fep->dhdl_derivatives == edhdlderivativesYES)
812 nsets_dhdl = n_lambda_terms;
814 /* count the number of delta_g states */
815 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
817 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
819 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
821 nsets += 1; /*add fep state for expanded ensemble */
824 if (fep->bPrintEnergy)
826 nsets += 1; /* add energy to the dhdl as well */
830 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
832 nsetsextend += 1; /* for PV term, other terms possible if required for
833 the reduced potential (only needed with foreign
834 lambda, and only output when init_lambda is not
835 set in order to maintain compatibility of the
839 snew(setname, nsetsextend);
841 if (expand->elmcmove > elmcmoveNO)
843 /* state for the fep_vals, if we have alchemical sampling */
844 sprintf(buf, "%s", "Thermodynamic state");
845 setname[s] = strdup(buf);
849 if (fep->bPrintEnergy)
851 sprintf(buf, "%s (%s)", "Energy", unit_energy);
852 setname[s] = strdup(buf);
856 if (fep->dhdl_derivatives == edhdlderivativesYES)
858 for (i = 0; i < efptNR; i++)
860 if (fep->separate_dvdl[i])
863 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
865 /* compatibility output */
866 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
870 double lam = fep->init_lambda;
871 if (fep->init_lambda < 0)
873 lam = fep->all_lambda[i][fep->init_fep_state];
875 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
878 setname[s] = strdup(buf);
884 if (fep->n_lambda > 0)
886 /* g_bar has to determine the lambda values used in this simulation
887 * from this xvg legend.
890 if (expand->elmcmove > elmcmoveNO)
892 nsetsbegin = 1; /* for including the expanded ensemble */
899 if (fep->bPrintEnergy)
903 nsetsbegin += nsets_dhdl;
905 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
907 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
908 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
910 /* for compatible dhdl.xvg files */
911 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
915 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
920 /* print the temperature for this state if doing simulated annealing */
921 sprintf(&buf[nps], "T = %g (%s)",
922 ir->simtempvals->temperatures[s-(nsetsbegin)],
925 setname[s] = strdup(buf);
930 np = sprintf(buf, "pV (%s)", unit_energy);
931 setname[nsetsextend-1] = strdup(buf); /* the first entry after
935 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
937 for (s = 0; s < nsetsextend; s++)
947 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
951 for (i = j = 0; (i < F_NRE); i++)
960 gmx_incons("Number of energy terms wrong");
964 void upd_mdebin(t_mdebin *md,
969 gmx_enerdata_t *enerd,
978 gmx_ekindata_t *ekind,
982 int i, j, k, kk, m, n, gid;
983 real crmsd[2], tmp6[6];
984 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
987 double store_dhdl[efptNR];
988 real store_energy = 0;
991 /* Do NOT use the box in the state variable, but the separate box provided
992 * as an argument. This is because we sometimes need to write the box from
993 * the last timestep to match the trajectory frames.
995 copy_energy(md, enerd->term, ecopy);
996 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
999 crmsd[0] = constr_rmsd(constr, FALSE);
1002 crmsd[1] = constr_rmsd(constr, TRUE);
1004 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1011 bs[0] = box[XX][XX];
1012 bs[1] = box[YY][YY];
1013 bs[2] = box[ZZ][ZZ];
1014 bs[3] = box[YY][XX];
1015 bs[4] = box[ZZ][XX];
1016 bs[5] = box[ZZ][YY];
1021 bs[0] = box[XX][XX];
1022 bs[1] = box[YY][YY];
1023 bs[2] = box[ZZ][ZZ];
1026 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1027 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1028 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1029 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1030 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1034 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1035 not the instantaneous pressure */
1036 pv = vol*md->ref_p/PRESFAC;
1038 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1039 enthalpy = pv + enerd->term[F_ETOT];
1040 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1045 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1046 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1050 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1054 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1058 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1059 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1061 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1063 tmp6[0] = state->boxv[XX][XX];
1064 tmp6[1] = state->boxv[YY][YY];
1065 tmp6[2] = state->boxv[ZZ][ZZ];
1066 tmp6[3] = state->boxv[YY][XX];
1067 tmp6[4] = state->boxv[ZZ][XX];
1068 tmp6[5] = state->boxv[ZZ][YY];
1069 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1073 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1075 if (ekind && ekind->cosacc.cos_accel != 0)
1077 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1078 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1079 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1080 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1081 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1082 *dens*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1083 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1088 for (i = 0; (i < md->nEg); i++)
1090 for (j = i; (j < md->nEg); j++)
1092 gid = GID(i, j, md->nEg);
1093 for (k = kk = 0; (k < egNR); k++)
1097 eee[kk++] = enerd->grpp.ener[k][gid];
1100 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1108 for (i = 0; (i < md->nTC); i++)
1110 md->tmp_r[i] = ekind->tcstat[i].T;
1112 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1114 if (md->etc == etcNOSEHOOVER)
1116 /* whether to print Nose-Hoover chains: */
1117 if (md->bPrintNHChains)
1119 if (md->bNHC_trotter)
1121 for (i = 0; (i < md->nTC); i++)
1123 for (j = 0; j < md->nNHC; j++)
1126 md->tmp_r[2*k] = state->nosehoover_xi[k];
1127 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1130 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1134 for (i = 0; (i < md->nTCP); i++)
1136 for (j = 0; j < md->nNHC; j++)
1139 md->tmp_r[2*k] = state->nhpres_xi[k];
1140 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1143 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1148 for (i = 0; (i < md->nTC); i++)
1150 md->tmp_r[2*i] = state->nosehoover_xi[i];
1151 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1153 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1157 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1158 md->etc == etcVRESCALE)
1160 for (i = 0; (i < md->nTC); i++)
1162 md->tmp_r[i] = ekind->tcstat[i].lambda;
1164 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1168 if (ekind && md->nU > 1)
1170 for (i = 0; (i < md->nU); i++)
1172 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1174 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1177 ebin_increase_count(md->ebin, bSum);
1179 /* BAR + thermodynamic integration values */
1180 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1182 for (i = 0; i < enerd->n_lambda-1; i++)
1184 /* zero for simulated tempering */
1185 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1186 if (md->temperatures != NULL)
1188 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1189 /* is this even useful to have at all? */
1190 md->dE[i] += (md->temperatures[i]/
1191 md->temperatures[state->fep_state]-1.0)*
1192 enerd->term[F_EKIN];
1198 fprintf(md->fp_dhdl, "%.4f", time);
1199 /* the current free energy state */
1201 /* print the current state if we are doing expanded ensemble */
1202 if (expand->elmcmove > elmcmoveNO)
1204 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1206 /* total energy (for if the temperature changes */
1207 if (fep->bPrintEnergy)
1209 store_energy = enerd->term[F_ETOT];
1210 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1213 if (fep->dhdl_derivatives == edhdlderivativesYES)
1215 for (i = 0; i < efptNR; i++)
1217 if (fep->separate_dvdl[i])
1219 /* assumes F_DVDL is first */
1220 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1224 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1226 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1228 if ((md->epc != epcNO) &&
1229 (enerd->n_lambda > 0) &&
1230 (fep->init_lambda < 0))
1232 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1233 there are alternate state
1234 lambda and we're not in
1235 compatibility mode */
1237 fprintf(md->fp_dhdl, "\n");
1238 /* and the binary free energy output */
1240 if (md->dhc && bDoDHDL)
1243 for (i = 0; i < efptNR; i++)
1245 if (fep->separate_dvdl[i])
1247 /* assumes F_DVDL is first */
1248 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1252 store_energy = enerd->term[F_ETOT];
1253 /* store_dh is dE */
1254 mde_delta_h_coll_add_dh(md->dhc,
1255 (double)state->fep_state,
1259 md->dE + fep->lambda_start_n,
1266 void upd_mdebin_step(t_mdebin *md)
1268 ebin_increase_count(md->ebin, FALSE);
1271 static void npr(FILE *log, int n, char c)
1273 for (; (n > 0); n--)
1275 fprintf(log, "%c", c);
1279 static void pprint(FILE *log, const char *s, t_mdebin *md)
1283 char buf1[22], buf2[22];
1286 fprintf(log, "\t<====== ");
1287 npr(log, slen, CHAR);
1288 fprintf(log, " ==>\n");
1289 fprintf(log, "\t<==== %s ====>\n", s);
1290 fprintf(log, "\t<== ");
1291 npr(log, slen, CHAR);
1292 fprintf(log, " ======>\n\n");
1294 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1295 gmx_step_str(md->ebin->nsteps_sim, buf1),
1296 gmx_step_str(md->ebin->nsum_sim, buf2));
1300 void print_ebin_header(FILE *log, gmx_large_int_t steps, double time, real lambda)
1304 fprintf(log, " %12s %12s %12s\n"
1305 " %12s %12.5f %12.5f\n\n",
1306 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1309 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1311 gmx_large_int_t step, double time,
1312 int mode, gmx_bool bCompact,
1313 t_mdebin *md, t_fcdata *fcd,
1314 gmx_groups_t *groups, t_grpopts *opts)
1316 /*static char **grpnms=NULL;*/
1318 int i, j, n, ni, nj, ndr, nor, b;
1320 real *disre_rm3tav, *disre_rt;
1322 /* these are for the old-style blocks (1 subblock, only reals), because
1323 there can be only one per ID for these */
1328 /* temporary arrays for the lambda values to write out */
1329 double enxlambda_data[2];
1339 fr.nsteps = md->ebin->nsteps;
1340 fr.dt = md->delta_t;
1341 fr.nsum = md->ebin->nsum;
1342 fr.nre = (bEne) ? md->ebin->nener : 0;
1343 fr.ener = md->ebin->e;
1344 ndisre = bDR ? fcd->disres.npair : 0;
1345 disre_rm3tav = fcd->disres.rm3tav;
1346 disre_rt = fcd->disres.rt;
1347 /* Optional additional old-style (real-only) blocks. */
1348 for (i = 0; i < enxNR; i++)
1352 if (fcd->orires.nr > 0 && bOR)
1354 diagonalize_orires_tensors(&(fcd->orires));
1355 nr[enxOR] = fcd->orires.nr;
1356 block[enxOR] = fcd->orires.otav;
1358 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1360 block[enxORI] = fcd->orires.oinsl;
1361 id[enxORI] = enxORI;
1362 nr[enxORT] = fcd->orires.nex*12;
1363 block[enxORT] = fcd->orires.eig;
1364 id[enxORT] = enxORT;
1367 /* whether we are going to wrte anything out: */
1368 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1371 /* the old-style blocks go first */
1373 for (i = 0; i < enxNR; i++)
1380 add_blocks_enxframe(&fr, fr.nblock);
1381 for (b = 0; b < fr.nblock; b++)
1383 add_subblocks_enxblock(&(fr.block[b]), 1);
1384 fr.block[b].id = id[b];
1385 fr.block[b].sub[0].nr = nr[b];
1387 fr.block[b].sub[0].type = xdr_datatype_float;
1388 fr.block[b].sub[0].fval = block[b];
1390 fr.block[b].sub[0].type = xdr_datatype_double;
1391 fr.block[b].sub[0].dval = block[b];
1395 /* check for disre block & fill it. */
1400 add_blocks_enxframe(&fr, fr.nblock);
1402 add_subblocks_enxblock(&(fr.block[db]), 2);
1403 fr.block[db].id = enxDISRE;
1404 fr.block[db].sub[0].nr = ndisre;
1405 fr.block[db].sub[1].nr = ndisre;
1407 fr.block[db].sub[0].type = xdr_datatype_float;
1408 fr.block[db].sub[1].type = xdr_datatype_float;
1409 fr.block[db].sub[0].fval = disre_rt;
1410 fr.block[db].sub[1].fval = disre_rm3tav;
1412 fr.block[db].sub[0].type = xdr_datatype_double;
1413 fr.block[db].sub[1].type = xdr_datatype_double;
1414 fr.block[db].sub[0].dval = disre_rt;
1415 fr.block[db].sub[1].dval = disre_rm3tav;
1418 /* here we can put new-style blocks */
1420 /* Free energy perturbation blocks */
1423 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1426 /* we can now free & reset the data in the blocks */
1429 mde_delta_h_coll_reset(md->dhc);
1432 /* do the actual I/O */
1433 do_enx(fp_ene, &fr);
1434 gmx_fio_check_file_position(enx_file_pointer(fp_ene));
1437 /* We have stored the sums, so reset the sum history */
1438 reset_ebin_sums(md->ebin);
1446 pprint(log, "A V E R A G E S", md);
1452 pprint(log, "R M S - F L U C T U A T I O N S", md);
1456 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1461 for (i = 0; i < opts->ngtc; i++)
1463 if (opts->annealing[i] != eannNO)
1465 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1466 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1470 if (mode == eprNORMAL && fcd->orires.nr > 0)
1472 print_orires_log(log, &(fcd->orires));
1474 fprintf(log, " Energies (%s)\n", unit_energy);
1475 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1482 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1488 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1489 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1491 fprintf(log, " Force Virial (%s)\n", unit_energy);
1492 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1497 fprintf(log, " Total Virial (%s)\n", unit_energy);
1498 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1503 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1504 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1509 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1510 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1516 if (md->print_grpnms == NULL)
1518 snew(md->print_grpnms, md->nE);
1520 for (i = 0; (i < md->nEg); i++)
1522 ni = groups->grps[egcENER].nm_ind[i];
1523 for (j = i; (j < md->nEg); j++)
1525 nj = groups->grps[egcENER].nm_ind[j];
1526 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1527 *(groups->grpname[nj]));
1528 md->print_grpnms[n++] = strdup(buf);
1532 sprintf(buf, "Epot (%s)", unit_energy);
1533 fprintf(log, "%15s ", buf);
1534 for (i = 0; (i < egNR); i++)
1538 fprintf(log, "%12s ", egrp_nm[i]);
1542 for (i = 0; (i < md->nE); i++)
1544 fprintf(log, "%15s", md->print_grpnms[i]);
1545 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1552 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1557 fprintf(log, "%15s %12s %12s %12s\n",
1558 "Group", "Ux", "Uy", "Uz");
1559 for (i = 0; (i < md->nU); i++)
1561 ni = groups->grps[egcACC].nm_ind[i];
1562 fprintf(log, "%15s", *groups->grpname[ni]);
1563 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1572 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1576 enerhist->nsteps = mdebin->ebin->nsteps;
1577 enerhist->nsum = mdebin->ebin->nsum;
1578 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1579 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1580 enerhist->nener = mdebin->ebin->nener;
1582 if (mdebin->ebin->nsum > 0)
1584 /* Check if we need to allocate first */
1585 if (enerhist->ener_ave == NULL)
1587 snew(enerhist->ener_ave, enerhist->nener);
1588 snew(enerhist->ener_sum, enerhist->nener);
1591 for (i = 0; i < enerhist->nener; i++)
1593 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1594 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1598 if (mdebin->ebin->nsum_sim > 0)
1600 /* Check if we need to allocate first */
1601 if (enerhist->ener_sum_sim == NULL)
1603 snew(enerhist->ener_sum_sim, enerhist->nener);
1606 for (i = 0; i < enerhist->nener; i++)
1608 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1613 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1617 void restore_energyhistory_from_state(t_mdebin * mdebin,
1618 energyhistory_t * enerhist)
1622 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1623 mdebin->ebin->nener != enerhist->nener)
1625 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1626 mdebin->ebin->nener, enerhist->nener);
1629 mdebin->ebin->nsteps = enerhist->nsteps;
1630 mdebin->ebin->nsum = enerhist->nsum;
1631 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1632 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1634 for (i = 0; i < mdebin->ebin->nener; i++)
1636 mdebin->ebin->e[i].eav =
1637 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1638 mdebin->ebin->e[i].esum =
1639 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1640 mdebin->ebin->e_sim[i].esum =
1641 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1645 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);