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 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/legacyheaders/mdebin.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/legacyheaders/disre.h"
51 #include "gromacs/legacyheaders/network.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/orires.h"
54 #include "gromacs/legacyheaders/constr.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/legacyheaders/macros.h"
59 #include "gromacs/legacyheaders/mdrun.h"
60 #include "mdebin_bar.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/utility/smalloc.h"
65 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
67 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
69 static const char *tricl_boxs_nm[] = {
70 "Box-XX", "Box-YY", "Box-ZZ",
71 "Box-YX", "Box-ZX", "Box-ZY"
74 static const char *vol_nm[] = { "Volume" };
76 static const char *dens_nm[] = {"Density" };
78 static const char *pv_nm[] = {"pV" };
80 static const char *enthalpy_nm[] = {"Enthalpy" };
82 static const char *boxvel_nm[] = {
83 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
84 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
87 #define NBOXS asize(boxs_nm)
88 #define NTRICLBOXS asize(tricl_boxs_nm)
90 t_mdebin *init_mdebin(ener_file_t fp_ene,
91 const gmx_mtop_t *mtop,
95 const char *ener_nm[F_NRE];
96 static const char *vir_nm[] = {
97 "Vir-XX", "Vir-XY", "Vir-XZ",
98 "Vir-YX", "Vir-YY", "Vir-YZ",
99 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
101 static const char *sv_nm[] = {
102 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
103 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
104 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
106 static const char *fv_nm[] = {
107 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
108 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
109 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
111 static const char *pres_nm[] = {
112 "Pres-XX", "Pres-XY", "Pres-XZ",
113 "Pres-YX", "Pres-YY", "Pres-YZ",
114 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
116 static const char *surft_nm[] = {
119 static const char *mu_nm[] = {
120 "Mu-X", "Mu-Y", "Mu-Z"
122 static const char *vcos_nm[] = {
125 static const char *visc_nm[] = {
128 static const char *baro_nm[] = {
133 const gmx_groups_t *groups;
138 int i, j, ni, nj, n, nh, k, kk, ncon, nset;
139 gmx_bool bBHAM, bNoseHoover, b14;
148 if (EI_DYNAMICS(ir->eI))
150 md->delta_t = ir->delta_t;
157 groups = &mtop->groups;
159 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
160 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
161 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
163 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
164 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
165 md->bConstr = (ncon > 0 || nset > 0);
166 md->bConstrVir = FALSE;
169 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
180 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
187 /* Energy monitoring */
188 for (i = 0; i < egNR; i++)
190 md->bEInd[i] = FALSE;
193 for (i = 0; i < F_NRE; i++)
195 md->bEner[i] = FALSE;
198 md->bEner[i] = !bBHAM;
200 else if (i == F_BHAM)
202 md->bEner[i] = bBHAM;
206 md->bEner[i] = ir->bQMMM;
208 else if (i == F_COUL_LR)
210 md->bEner[i] = (ir->rcoulomb > ir->rlist);
212 else if (i == F_LJ_LR)
214 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
216 else if (i == F_BHAM_LR)
218 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
220 else if (i == F_RF_EXCL)
222 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
224 else if (i == F_COUL_RECIP)
226 md->bEner[i] = EEL_FULL(ir->coulombtype);
228 else if (i == F_LJ_RECIP)
230 md->bEner[i] = EVDW_PME(ir->vdwtype);
232 else if (i == F_LJ14)
236 else if (i == F_COUL14)
240 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
242 md->bEner[i] = FALSE;
244 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
245 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
246 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
247 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
248 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
249 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
251 md->bEner[i] = (ir->efep != efepNO);
253 else if ((interaction_function[i].flags & IF_VSITE) ||
254 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
256 md->bEner[i] = FALSE;
258 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
262 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
266 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
270 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
272 md->bEner[i] = FALSE;
274 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
276 md->bEner[i] = EI_DYNAMICS(ir->eI);
278 else if (i == F_DISPCORR || i == F_PDISPCORR)
280 md->bEner[i] = (ir->eDispCorr != edispcNO);
282 else if (i == F_DISRESVIOL)
284 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
286 else if (i == F_ORIRESDEV)
288 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
290 else if (i == F_CONNBONDS)
292 md->bEner[i] = FALSE;
294 else if (i == F_COM_PULL)
296 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
298 else if (i == F_ECONSERVED)
300 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
301 (ir->epc == epcNO || ir->epc == epcMTTK));
305 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
309 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
310 if (ir->bAdress && !debug)
312 for (i = 0; i < F_NRE; i++)
314 md->bEner[i] = FALSE;
331 for (i = 0; i < F_NRE; i++)
335 ener_nm[md->f_nre] = interaction_function[i].longname;
341 md->bDiagPres = !TRICLINIC(ir->ref_p);
342 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
343 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
344 md->bDynBox = DYNAMIC_BOX(*ir);
346 md->bNHC_trotter = IR_NVT_TROTTER(ir);
347 md->bPrintNHChains = ir->bPrintNHChains;
348 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
349 md->bMu = NEED_MUTOT(*ir);
351 md->ebin = mk_ebin();
352 /* Pass NULL for unit to let get_ebin_space determine the units
353 * for interaction_function[i].longname
355 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
358 /* This should be called directly after the call for md->ie,
359 * such that md->iconrmsd follows directly in the list.
361 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
365 md->ib = get_ebin_space(md->ebin,
366 md->bTricl ? NTRICLBOXS : NBOXS,
367 md->bTricl ? tricl_boxs_nm : boxs_nm,
369 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
370 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
373 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
374 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
379 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
380 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
384 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
388 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
392 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
395 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
397 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
398 boxvel_nm, unit_vel);
402 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
404 if (ir->cos_accel != 0)
406 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
407 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
411 /* Energy monitoring */
412 for (i = 0; i < egNR; i++)
414 md->bEInd[i] = FALSE;
416 md->bEInd[egCOULSR] = TRUE;
417 md->bEInd[egLJSR ] = TRUE;
419 if (ir->rcoulomb > ir->rlist)
421 md->bEInd[egCOULLR] = TRUE;
425 if (ir->rvdw > ir->rlist)
427 md->bEInd[egLJLR] = TRUE;
432 md->bEInd[egLJSR] = FALSE;
433 md->bEInd[egBHAMSR] = TRUE;
434 if (ir->rvdw > ir->rlist)
436 md->bEInd[egBHAMLR] = TRUE;
441 md->bEInd[egLJ14] = TRUE;
442 md->bEInd[egCOUL14] = TRUE;
445 for (i = 0; (i < egNR); i++)
453 n = groups->grps[egcENER].nr;
454 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
457 /*standard simulation*/
459 md->nE = (n*(n+1))/2;
463 /*AdResS simulation*/
470 snew(md->igrp, md->nE);
475 for (k = 0; (k < md->nEc); k++)
477 snew(gnm[k], STRLEN);
479 for (i = 0; (i < groups->grps[egcENER].nr); i++)
481 ni = groups->grps[egcENER].nm_ind[i];
482 for (j = i; (j < groups->grps[egcENER].nr); j++)
484 nj = groups->grps[egcENER].nm_ind[j];
485 for (k = kk = 0; (k < egNR); k++)
489 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
490 *(groups->grpname[ni]), *(groups->grpname[nj]));
494 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
495 (const char **)gnm, unit_energy);
499 for (k = 0; (k < md->nEc); k++)
507 gmx_incons("Number of energy terms wrong");
511 md->nTC = groups->grps[egcTC].nr;
512 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
515 md->nTCP = 1; /* assume only one possible coupling system for barostat
522 if (md->etc == etcNOSEHOOVER)
524 if (md->bNHC_trotter)
526 md->mde_n = 2*md->nNHC*md->nTC;
530 md->mde_n = 2*md->nTC;
532 if (md->epc == epcMTTK)
534 md->mdeb_n = 2*md->nNHC*md->nTCP;
543 snew(md->tmp_r, md->mde_n);
544 snew(md->tmp_v, md->mde_n);
545 snew(md->grpnms, md->mde_n);
548 for (i = 0; (i < md->nTC); i++)
550 ni = groups->grps[egcTC].nm_ind[i];
551 sprintf(buf, "T-%s", *(groups->grpname[ni]));
552 grpnms[i] = gmx_strdup(buf);
554 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
557 if (md->etc == etcNOSEHOOVER)
559 if (md->bPrintNHChains)
561 if (md->bNHC_trotter)
563 for (i = 0; (i < md->nTC); i++)
565 ni = groups->grps[egcTC].nm_ind[i];
566 bufi = *(groups->grpname[ni]);
567 for (j = 0; (j < md->nNHC); j++)
569 sprintf(buf, "Xi-%d-%s", j, bufi);
570 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
571 sprintf(buf, "vXi-%d-%s", j, bufi);
572 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
575 md->itc = get_ebin_space(md->ebin, md->mde_n,
576 (const char **)grpnms, unit_invtime);
579 for (i = 0; (i < md->nTCP); i++)
581 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
582 for (j = 0; (j < md->nNHC); j++)
584 sprintf(buf, "Xi-%d-%s", j, bufi);
585 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
586 sprintf(buf, "vXi-%d-%s", j, bufi);
587 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
590 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
591 (const char **)grpnms, unit_invtime);
596 for (i = 0; (i < md->nTC); i++)
598 ni = groups->grps[egcTC].nm_ind[i];
599 bufi = *(groups->grpname[ni]);
600 sprintf(buf, "Xi-%s", bufi);
601 grpnms[2*i] = gmx_strdup(buf);
602 sprintf(buf, "vXi-%s", bufi);
603 grpnms[2*i+1] = gmx_strdup(buf);
605 md->itc = get_ebin_space(md->ebin, md->mde_n,
606 (const char **)grpnms, unit_invtime);
610 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
611 md->etc == etcVRESCALE)
613 for (i = 0; (i < md->nTC); i++)
615 ni = groups->grps[egcTC].nm_ind[i];
616 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
617 grpnms[i] = gmx_strdup(buf);
619 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
625 md->nU = groups->grps[egcACC].nr;
628 snew(grpnms, 3*md->nU);
629 for (i = 0; (i < md->nU); i++)
631 ni = groups->grps[egcACC].nm_ind[i];
632 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
633 grpnms[3*i+XX] = gmx_strdup(buf);
634 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
635 grpnms[3*i+YY] = gmx_strdup(buf);
636 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
637 grpnms[3*i+ZZ] = gmx_strdup(buf);
639 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
645 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
648 md->print_grpnms = NULL;
650 /* check whether we're going to write dh histograms */
652 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
654 /* Currently dh histograms are only written with dynamics */
655 if (EI_DYNAMICS(ir->eI))
659 mde_delta_h_coll_init(md->dhc, ir);
662 snew(md->dE, ir->fepvals->n_lambda);
666 md->fp_dhdl = fp_dhdl;
667 snew(md->dE, ir->fepvals->n_lambda);
672 snew(md->temperatures, ir->fepvals->n_lambda);
673 for (i = 0; i < ir->fepvals->n_lambda; i++)
675 md->temperatures[i] = ir->simtempvals->temperatures[i];
681 /* print a lambda vector to a string
682 fep = the inputrec's FEP input data
683 i = the index of the lambda vector
684 get_native_lambda = whether to print the native lambda
685 get_names = whether to print the names rather than the values
686 str = the pre-allocated string buffer to print to. */
687 static void print_lambda_vector(t_lambda *fep, int i,
688 gmx_bool get_native_lambda, gmx_bool get_names,
695 for (j = 0; j < efptNR; j++)
697 if (fep->separate_dvdl[j])
702 str[0] = 0; /* reset the string */
705 str += sprintf(str, "("); /* set the opening parenthesis*/
707 for (j = 0; j < efptNR; j++)
709 if (fep->separate_dvdl[j])
714 if (get_native_lambda && fep->init_lambda >= 0)
716 str += sprintf(str, "%.4f", fep->init_lambda);
720 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
725 str += sprintf(str, "%s", efpt_singular_names[j]);
727 /* print comma for the next item */
730 str += sprintf(str, ", ");
737 /* and add the closing parenthesis */
743 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
744 const output_env_t oenv)
747 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
748 *lambdastate = "\\lambda state", *remain = "remaining";
749 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
750 int i, np, nps, nsets, nsets_de, nsetsbegin;
751 int n_lambda_terms = 0;
752 t_lambda *fep = ir->fepvals; /* for simplicity */
753 t_expanded *expand = ir->expandedvals;
755 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
761 gmx_bool write_pV = FALSE;
763 /* count the number of different lambda terms */
764 for (i = 0; i < efptNR; i++)
766 if (fep->separate_dvdl[i])
772 if (fep->n_lambda == 0)
774 sprintf(title, "%s", dhdl);
775 sprintf(label_x, "Time (ps)");
776 sprintf(label_y, "%s (%s %s)",
777 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
781 sprintf(title, "%s and %s", dhdl, deltag);
782 sprintf(label_x, "Time (ps)");
783 sprintf(label_y, "%s and %s (%s %s)",
784 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
786 fp = gmx_fio_fopen(filename, "w+");
787 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
791 bufplace = sprintf(buf, "T = %g (K) ",
794 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
796 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
798 /* compatibility output */
799 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
803 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
805 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
807 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
808 lambdastate, fep->init_fep_state,
809 lambda_name_str, lambda_vec_str);
812 xvgr_subtitle(fp, buf, oenv);
816 if (fep->dhdl_derivatives == edhdlderivativesYES)
818 nsets_dhdl = n_lambda_terms;
820 /* count the number of delta_g states */
821 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
823 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
825 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
827 nsets += 1; /*add fep state for expanded ensemble */
830 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
832 nsets += 1; /* add energy to the dhdl as well */
836 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
838 nsetsextend += 1; /* for PV term, other terms possible if required for
839 the reduced potential (only needed with foreign
840 lambda, and only output when init_lambda is not
841 set in order to maintain compatibility of the
845 snew(setname, nsetsextend);
847 if (expand->elmcmove > elmcmoveNO)
849 /* state for the fep_vals, if we have alchemical sampling */
850 sprintf(buf, "%s", "Thermodynamic state");
851 setname[s] = gmx_strdup(buf);
855 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
857 switch (fep->edHdLPrintEnergy)
859 case edHdLPrintEnergyPOTENTIAL:
860 sprintf(buf, "%s (%s)", "Potential Energy", unit_energy);
862 case edHdLPrintEnergyTOTAL:
863 case edHdLPrintEnergyYES:
865 sprintf(buf, "%s (%s)", "Total Energy", unit_energy);
867 setname[s] = gmx_strdup(buf);
871 if (fep->dhdl_derivatives == edhdlderivativesYES)
873 for (i = 0; i < efptNR; i++)
875 if (fep->separate_dvdl[i])
878 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
880 /* compatibility output */
881 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
885 double lam = fep->init_lambda;
886 if (fep->init_lambda < 0)
888 lam = fep->all_lambda[i][fep->init_fep_state];
890 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
893 setname[s] = gmx_strdup(buf);
899 if (fep->n_lambda > 0)
901 /* g_bar has to determine the lambda values used in this simulation
902 * from this xvg legend.
905 if (expand->elmcmove > elmcmoveNO)
907 nsetsbegin = 1; /* for including the expanded ensemble */
914 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
918 nsetsbegin += nsets_dhdl;
920 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
922 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
923 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
925 /* for compatible dhdl.xvg files */
926 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
930 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
935 /* print the temperature for this state if doing simulated annealing */
936 sprintf(&buf[nps], "T = %g (%s)",
937 ir->simtempvals->temperatures[s-(nsetsbegin)],
940 setname[s] = gmx_strdup(buf);
945 np = sprintf(buf, "pV (%s)", unit_energy);
946 setname[nsetsextend-1] = gmx_strdup(buf); /* the first entry after
950 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
952 for (s = 0; s < nsetsextend; s++)
962 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
966 for (i = j = 0; (i < F_NRE); i++)
975 gmx_incons("Number of energy terms wrong");
979 void upd_mdebin(t_mdebin *md,
984 gmx_enerdata_t *enerd,
993 gmx_ekindata_t *ekind,
997 int i, j, k, kk, m, n, gid;
998 real crmsd[2], tmp6[6];
999 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
1002 double store_dhdl[efptNR];
1003 real store_energy = 0;
1006 /* Do NOT use the box in the state variable, but the separate box provided
1007 * as an argument. This is because we sometimes need to write the box from
1008 * the last timestep to match the trajectory frames.
1010 copy_energy(md, enerd->term, ecopy);
1011 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1014 crmsd[0] = constr_rmsd(constr, FALSE);
1017 crmsd[1] = constr_rmsd(constr, TRUE);
1019 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1026 bs[0] = box[XX][XX];
1027 bs[1] = box[YY][YY];
1028 bs[2] = box[ZZ][ZZ];
1029 bs[3] = box[YY][XX];
1030 bs[4] = box[ZZ][XX];
1031 bs[5] = box[ZZ][YY];
1036 bs[0] = box[XX][XX];
1037 bs[1] = box[YY][YY];
1038 bs[2] = box[ZZ][ZZ];
1041 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1042 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1043 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1044 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1045 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1049 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1050 not the instantaneous pressure */
1051 pv = vol*md->ref_p/PRESFAC;
1053 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1054 enthalpy = pv + enerd->term[F_ETOT];
1055 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1060 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1061 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1065 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1069 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1073 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1074 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1076 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1078 tmp6[0] = state->boxv[XX][XX];
1079 tmp6[1] = state->boxv[YY][YY];
1080 tmp6[2] = state->boxv[ZZ][ZZ];
1081 tmp6[3] = state->boxv[YY][XX];
1082 tmp6[4] = state->boxv[ZZ][XX];
1083 tmp6[5] = state->boxv[ZZ][YY];
1084 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1088 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1090 if (ekind && ekind->cosacc.cos_accel != 0)
1092 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1093 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1094 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1095 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1096 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1097 *dens*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1098 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1103 for (i = 0; (i < md->nEg); i++)
1105 for (j = i; (j < md->nEg); j++)
1107 gid = GID(i, j, md->nEg);
1108 for (k = kk = 0; (k < egNR); k++)
1112 eee[kk++] = enerd->grpp.ener[k][gid];
1115 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1123 for (i = 0; (i < md->nTC); i++)
1125 md->tmp_r[i] = ekind->tcstat[i].T;
1127 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1129 if (md->etc == etcNOSEHOOVER)
1131 /* whether to print Nose-Hoover chains: */
1132 if (md->bPrintNHChains)
1134 if (md->bNHC_trotter)
1136 for (i = 0; (i < md->nTC); i++)
1138 for (j = 0; j < md->nNHC; j++)
1141 md->tmp_r[2*k] = state->nosehoover_xi[k];
1142 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1145 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1149 for (i = 0; (i < md->nTCP); i++)
1151 for (j = 0; j < md->nNHC; j++)
1154 md->tmp_r[2*k] = state->nhpres_xi[k];
1155 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1158 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1163 for (i = 0; (i < md->nTC); i++)
1165 md->tmp_r[2*i] = state->nosehoover_xi[i];
1166 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1168 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1172 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1173 md->etc == etcVRESCALE)
1175 for (i = 0; (i < md->nTC); i++)
1177 md->tmp_r[i] = ekind->tcstat[i].lambda;
1179 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1183 if (ekind && md->nU > 1)
1185 for (i = 0; (i < md->nU); i++)
1187 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1189 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1192 ebin_increase_count(md->ebin, bSum);
1194 /* BAR + thermodynamic integration values */
1195 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1197 for (i = 0; i < enerd->n_lambda-1; i++)
1199 /* zero for simulated tempering */
1200 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1201 if (md->temperatures != NULL)
1203 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1204 /* is this even useful to have at all? */
1205 md->dE[i] += (md->temperatures[i]/
1206 md->temperatures[state->fep_state]-1.0)*
1207 enerd->term[F_EKIN];
1213 fprintf(md->fp_dhdl, "%.4f", time);
1214 /* the current free energy state */
1216 /* print the current state if we are doing expanded ensemble */
1217 if (expand->elmcmove > elmcmoveNO)
1219 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1221 /* total energy (for if the temperature changes */
1223 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1225 switch (fep->edHdLPrintEnergy)
1227 case edHdLPrintEnergyPOTENTIAL:
1228 store_energy = enerd->term[F_EPOT];
1230 case edHdLPrintEnergyTOTAL:
1231 case edHdLPrintEnergyYES:
1233 store_energy = enerd->term[F_ETOT];
1235 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1238 if (fep->dhdl_derivatives == edhdlderivativesYES)
1240 for (i = 0; i < efptNR; i++)
1242 if (fep->separate_dvdl[i])
1244 /* assumes F_DVDL is first */
1245 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1249 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1251 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1253 if ((md->epc != epcNO) &&
1254 (enerd->n_lambda > 0) &&
1255 (fep->init_lambda < 0))
1257 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1258 there are alternate state
1259 lambda and we're not in
1260 compatibility mode */
1262 fprintf(md->fp_dhdl, "\n");
1263 /* and the binary free energy output */
1265 if (md->dhc && bDoDHDL)
1268 for (i = 0; i < efptNR; i++)
1270 if (fep->separate_dvdl[i])
1272 /* assumes F_DVDL is first */
1273 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1277 store_energy = enerd->term[F_ETOT];
1278 /* store_dh is dE */
1279 mde_delta_h_coll_add_dh(md->dhc,
1280 (double)state->fep_state,
1284 md->dE + fep->lambda_start_n,
1291 void upd_mdebin_step(t_mdebin *md)
1293 ebin_increase_count(md->ebin, FALSE);
1296 static void npr(FILE *log, int n, char c)
1298 for (; (n > 0); n--)
1300 fprintf(log, "%c", c);
1304 static void pprint(FILE *log, const char *s, t_mdebin *md)
1308 char buf1[22], buf2[22];
1311 fprintf(log, "\t<====== ");
1312 npr(log, slen, CHAR);
1313 fprintf(log, " ==>\n");
1314 fprintf(log, "\t<==== %s ====>\n", s);
1315 fprintf(log, "\t<== ");
1316 npr(log, slen, CHAR);
1317 fprintf(log, " ======>\n\n");
1319 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1320 gmx_step_str(md->ebin->nsteps_sim, buf1),
1321 gmx_step_str(md->ebin->nsum_sim, buf2));
1325 void print_ebin_header(FILE *log, gmx_int64_t steps, double time, real lambda)
1329 fprintf(log, " %12s %12s %12s\n"
1330 " %12s %12.5f %12.5f\n\n",
1331 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1334 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1336 gmx_int64_t step, double time,
1337 int mode, gmx_bool bCompact,
1338 t_mdebin *md, t_fcdata *fcd,
1339 gmx_groups_t *groups, t_grpopts *opts)
1341 /*static char **grpnms=NULL;*/
1343 int i, j, n, ni, nj, ndr, nor, b;
1345 real *disre_rm3tav, *disre_rt;
1347 /* these are for the old-style blocks (1 subblock, only reals), because
1348 there can be only one per ID for these */
1353 /* temporary arrays for the lambda values to write out */
1354 double enxlambda_data[2];
1364 fr.nsteps = md->ebin->nsteps;
1365 fr.dt = md->delta_t;
1366 fr.nsum = md->ebin->nsum;
1367 fr.nre = (bEne) ? md->ebin->nener : 0;
1368 fr.ener = md->ebin->e;
1369 ndisre = bDR ? fcd->disres.npair : 0;
1370 disre_rm3tav = fcd->disres.rm3tav;
1371 disre_rt = fcd->disres.rt;
1372 /* Optional additional old-style (real-only) blocks. */
1373 for (i = 0; i < enxNR; i++)
1377 if (fcd->orires.nr > 0 && bOR)
1379 diagonalize_orires_tensors(&(fcd->orires));
1380 nr[enxOR] = fcd->orires.nr;
1381 block[enxOR] = fcd->orires.otav;
1383 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1385 block[enxORI] = fcd->orires.oinsl;
1386 id[enxORI] = enxORI;
1387 nr[enxORT] = fcd->orires.nex*12;
1388 block[enxORT] = fcd->orires.eig;
1389 id[enxORT] = enxORT;
1392 /* whether we are going to wrte anything out: */
1393 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1396 /* the old-style blocks go first */
1398 for (i = 0; i < enxNR; i++)
1405 add_blocks_enxframe(&fr, fr.nblock);
1406 for (b = 0; b < fr.nblock; b++)
1408 add_subblocks_enxblock(&(fr.block[b]), 1);
1409 fr.block[b].id = id[b];
1410 fr.block[b].sub[0].nr = nr[b];
1412 fr.block[b].sub[0].type = xdr_datatype_float;
1413 fr.block[b].sub[0].fval = block[b];
1415 fr.block[b].sub[0].type = xdr_datatype_double;
1416 fr.block[b].sub[0].dval = block[b];
1420 /* check for disre block & fill it. */
1425 add_blocks_enxframe(&fr, fr.nblock);
1427 add_subblocks_enxblock(&(fr.block[db]), 2);
1428 fr.block[db].id = enxDISRE;
1429 fr.block[db].sub[0].nr = ndisre;
1430 fr.block[db].sub[1].nr = ndisre;
1432 fr.block[db].sub[0].type = xdr_datatype_float;
1433 fr.block[db].sub[1].type = xdr_datatype_float;
1434 fr.block[db].sub[0].fval = disre_rt;
1435 fr.block[db].sub[1].fval = disre_rm3tav;
1437 fr.block[db].sub[0].type = xdr_datatype_double;
1438 fr.block[db].sub[1].type = xdr_datatype_double;
1439 fr.block[db].sub[0].dval = disre_rt;
1440 fr.block[db].sub[1].dval = disre_rm3tav;
1443 /* here we can put new-style blocks */
1445 /* Free energy perturbation blocks */
1448 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1451 /* we can now free & reset the data in the blocks */
1454 mde_delta_h_coll_reset(md->dhc);
1457 /* do the actual I/O */
1458 do_enx(fp_ene, &fr);
1461 /* We have stored the sums, so reset the sum history */
1462 reset_ebin_sums(md->ebin);
1470 pprint(log, "A V E R A G E S", md);
1476 pprint(log, "R M S - F L U C T U A T I O N S", md);
1480 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1485 for (i = 0; i < opts->ngtc; i++)
1487 if (opts->annealing[i] != eannNO)
1489 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1490 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1494 if (mode == eprNORMAL && fcd->orires.nr > 0)
1496 print_orires_log(log, &(fcd->orires));
1498 fprintf(log, " Energies (%s)\n", unit_energy);
1499 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1506 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1512 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1513 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1515 fprintf(log, " Force Virial (%s)\n", unit_energy);
1516 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1521 fprintf(log, " Total Virial (%s)\n", unit_energy);
1522 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1527 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1528 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1533 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1534 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1540 if (md->print_grpnms == NULL)
1542 snew(md->print_grpnms, md->nE);
1544 for (i = 0; (i < md->nEg); i++)
1546 ni = groups->grps[egcENER].nm_ind[i];
1547 for (j = i; (j < md->nEg); j++)
1549 nj = groups->grps[egcENER].nm_ind[j];
1550 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1551 *(groups->grpname[nj]));
1552 md->print_grpnms[n++] = gmx_strdup(buf);
1556 sprintf(buf, "Epot (%s)", unit_energy);
1557 fprintf(log, "%15s ", buf);
1558 for (i = 0; (i < egNR); i++)
1562 fprintf(log, "%12s ", egrp_nm[i]);
1566 for (i = 0; (i < md->nE); i++)
1568 fprintf(log, "%15s", md->print_grpnms[i]);
1569 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1576 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1581 fprintf(log, "%15s %12s %12s %12s\n",
1582 "Group", "Ux", "Uy", "Uz");
1583 for (i = 0; (i < md->nU); i++)
1585 ni = groups->grps[egcACC].nm_ind[i];
1586 fprintf(log, "%15s", *groups->grpname[ni]);
1587 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1596 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1600 enerhist->nsteps = mdebin->ebin->nsteps;
1601 enerhist->nsum = mdebin->ebin->nsum;
1602 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1603 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1604 enerhist->nener = mdebin->ebin->nener;
1606 if (mdebin->ebin->nsum > 0)
1608 /* Check if we need to allocate first */
1609 if (enerhist->ener_ave == NULL)
1611 snew(enerhist->ener_ave, enerhist->nener);
1612 snew(enerhist->ener_sum, enerhist->nener);
1615 for (i = 0; i < enerhist->nener; i++)
1617 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1618 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1622 if (mdebin->ebin->nsum_sim > 0)
1624 /* Check if we need to allocate first */
1625 if (enerhist->ener_sum_sim == NULL)
1627 snew(enerhist->ener_sum_sim, enerhist->nener);
1630 for (i = 0; i < enerhist->nener; i++)
1632 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1637 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1641 void restore_energyhistory_from_state(t_mdebin * mdebin,
1642 energyhistory_t * enerhist)
1646 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1647 mdebin->ebin->nener != enerhist->nener)
1649 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1650 mdebin->ebin->nener, enerhist->nener);
1653 mdebin->ebin->nsteps = enerhist->nsteps;
1654 mdebin->ebin->nsum = enerhist->nsum;
1655 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1656 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1658 for (i = 0; i < mdebin->ebin->nener; i++)
1660 mdebin->ebin->e[i].eav =
1661 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1662 mdebin->ebin->e[i].esum =
1663 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1664 mdebin->ebin->e_sim[i].esum =
1665 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1669 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);