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,2015, 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.
39 #include "gromacs/legacyheaders/mdebin.h"
45 #include "gromacs/fileio/enxio.h"
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/legacyheaders/constr.h"
49 #include "gromacs/legacyheaders/disre.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/mdrun.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/network.h"
54 #include "gromacs/legacyheaders/orires.h"
55 #include "gromacs/legacyheaders/typedefs.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/mdebin_bar.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/pulling/pull.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/utility/smalloc.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, k, kk, ncon, nset;
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;
168 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
179 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
186 /* Energy monitoring */
187 for (i = 0; i < egNR; i++)
189 md->bEInd[i] = FALSE;
192 for (i = 0; i < F_NRE; i++)
194 md->bEner[i] = FALSE;
197 md->bEner[i] = !bBHAM;
199 else if (i == F_BHAM)
201 md->bEner[i] = bBHAM;
205 md->bEner[i] = ir->bQMMM;
207 else if (i == F_COUL_LR)
209 md->bEner[i] = (ir->rcoulomb > ir->rlist);
211 else if (i == F_LJ_LR)
213 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
215 else if (i == F_BHAM_LR)
217 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
219 else if (i == F_RF_EXCL)
221 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
223 else if (i == F_COUL_RECIP)
225 md->bEner[i] = EEL_FULL(ir->coulombtype);
227 else if (i == F_LJ_RECIP)
229 md->bEner[i] = EVDW_PME(ir->vdwtype);
231 else if (i == F_LJ14)
235 else if (i == F_COUL14)
239 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
241 md->bEner[i] = FALSE;
243 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
244 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
245 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
246 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
247 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
248 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
250 md->bEner[i] = (ir->efep != efepNO);
252 else if ((interaction_function[i].flags & IF_VSITE) ||
253 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
255 md->bEner[i] = FALSE;
257 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
261 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
265 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
269 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
271 md->bEner[i] = FALSE;
273 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
275 md->bEner[i] = EI_DYNAMICS(ir->eI);
277 else if (i == F_DISPCORR || i == F_PDISPCORR)
279 md->bEner[i] = (ir->eDispCorr != edispcNO);
281 else if (i == F_DISRESVIOL)
283 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
285 else if (i == F_ORIRESDEV)
287 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
289 else if (i == F_CONNBONDS)
291 md->bEner[i] = FALSE;
293 else if (i == F_COM_PULL)
295 md->bEner[i] = (ir->bPull && pull_have_potential(ir->pull_work));
297 else if (i == F_ECONSERVED)
299 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
300 (ir->epc == epcNO || ir->epc == epcMTTK));
304 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
308 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
309 if (ir->bAdress && !debug)
311 for (i = 0; i < F_NRE; i++)
313 md->bEner[i] = FALSE;
330 for (i = 0; i < F_NRE; i++)
334 ener_nm[md->f_nre] = interaction_function[i].longname;
340 md->bDiagPres = !TRICLINIC(ir->ref_p);
341 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
342 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
343 md->bDynBox = DYNAMIC_BOX(*ir);
345 md->bNHC_trotter = IR_NVT_TROTTER(ir);
346 md->bPrintNHChains = ir->bPrintNHChains;
347 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
348 md->bMu = NEED_MUTOT(*ir);
350 md->ebin = mk_ebin();
351 /* Pass NULL for unit to let get_ebin_space determine the units
352 * for interaction_function[i].longname
354 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
357 /* This should be called directly after the call for md->ie,
358 * such that md->iconrmsd follows directly in the list.
360 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
364 md->ib = get_ebin_space(md->ebin,
365 md->bTricl ? NTRICLBOXS : NBOXS,
366 md->bTricl ? tricl_boxs_nm : boxs_nm,
368 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
369 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
372 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
373 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
378 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
379 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
383 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
387 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
391 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
394 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
396 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
397 boxvel_nm, unit_vel);
401 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
403 if (ir->cos_accel != 0)
405 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
406 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
410 /* Energy monitoring */
411 for (i = 0; i < egNR; i++)
413 md->bEInd[i] = FALSE;
415 md->bEInd[egCOULSR] = TRUE;
416 md->bEInd[egLJSR ] = TRUE;
418 if (ir->rcoulomb > ir->rlist)
420 md->bEInd[egCOULLR] = TRUE;
424 if (ir->rvdw > ir->rlist)
426 md->bEInd[egLJLR] = TRUE;
431 md->bEInd[egLJSR] = FALSE;
432 md->bEInd[egBHAMSR] = TRUE;
433 if (ir->rvdw > ir->rlist)
435 md->bEInd[egBHAMLR] = TRUE;
440 md->bEInd[egLJ14] = TRUE;
441 md->bEInd[egCOUL14] = TRUE;
444 for (i = 0; (i < egNR); i++)
452 n = groups->grps[egcENER].nr;
453 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
456 /*standard simulation*/
458 md->nE = (n*(n+1))/2;
462 /*AdResS simulation*/
469 snew(md->igrp, md->nE);
474 for (k = 0; (k < md->nEc); k++)
476 snew(gnm[k], STRLEN);
478 for (i = 0; (i < groups->grps[egcENER].nr); i++)
480 ni = groups->grps[egcENER].nm_ind[i];
481 for (j = i; (j < groups->grps[egcENER].nr); j++)
483 nj = groups->grps[egcENER].nm_ind[j];
484 for (k = kk = 0; (k < egNR); k++)
488 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
489 *(groups->grpname[ni]), *(groups->grpname[nj]));
493 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
494 (const char **)gnm, unit_energy);
498 for (k = 0; (k < md->nEc); k++)
506 gmx_incons("Number of energy terms wrong");
510 md->nTC = groups->grps[egcTC].nr;
511 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
514 md->nTCP = 1; /* assume only one possible coupling system for barostat
521 if (md->etc == etcNOSEHOOVER)
523 if (md->bNHC_trotter)
525 md->mde_n = 2*md->nNHC*md->nTC;
529 md->mde_n = 2*md->nTC;
531 if (md->epc == epcMTTK)
533 md->mdeb_n = 2*md->nNHC*md->nTCP;
542 snew(md->tmp_r, md->mde_n);
543 snew(md->tmp_v, md->mde_n);
544 snew(md->grpnms, md->mde_n);
547 for (i = 0; (i < md->nTC); i++)
549 ni = groups->grps[egcTC].nm_ind[i];
550 sprintf(buf, "T-%s", *(groups->grpname[ni]));
551 grpnms[i] = gmx_strdup(buf);
553 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
556 if (md->etc == etcNOSEHOOVER)
558 if (md->bPrintNHChains)
560 if (md->bNHC_trotter)
562 for (i = 0; (i < md->nTC); i++)
564 ni = groups->grps[egcTC].nm_ind[i];
565 bufi = *(groups->grpname[ni]);
566 for (j = 0; (j < md->nNHC); j++)
568 sprintf(buf, "Xi-%d-%s", j, bufi);
569 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
570 sprintf(buf, "vXi-%d-%s", j, bufi);
571 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
574 md->itc = get_ebin_space(md->ebin, md->mde_n,
575 (const char **)grpnms, unit_invtime);
578 for (i = 0; (i < md->nTCP); i++)
580 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
581 for (j = 0; (j < md->nNHC); j++)
583 sprintf(buf, "Xi-%d-%s", j, bufi);
584 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
585 sprintf(buf, "vXi-%d-%s", j, bufi);
586 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
589 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
590 (const char **)grpnms, unit_invtime);
595 for (i = 0; (i < md->nTC); i++)
597 ni = groups->grps[egcTC].nm_ind[i];
598 bufi = *(groups->grpname[ni]);
599 sprintf(buf, "Xi-%s", bufi);
600 grpnms[2*i] = gmx_strdup(buf);
601 sprintf(buf, "vXi-%s", bufi);
602 grpnms[2*i+1] = gmx_strdup(buf);
604 md->itc = get_ebin_space(md->ebin, md->mde_n,
605 (const char **)grpnms, unit_invtime);
609 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
610 md->etc == etcVRESCALE)
612 for (i = 0; (i < md->nTC); i++)
614 ni = groups->grps[egcTC].nm_ind[i];
615 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
616 grpnms[i] = gmx_strdup(buf);
618 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
624 md->nU = groups->grps[egcACC].nr;
627 snew(grpnms, 3*md->nU);
628 for (i = 0; (i < md->nU); i++)
630 ni = groups->grps[egcACC].nm_ind[i];
631 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
632 grpnms[3*i+XX] = gmx_strdup(buf);
633 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
634 grpnms[3*i+YY] = gmx_strdup(buf);
635 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
636 grpnms[3*i+ZZ] = gmx_strdup(buf);
638 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
644 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
647 md->print_grpnms = NULL;
649 /* check whether we're going to write dh histograms */
651 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
653 /* Currently dh histograms are only written with dynamics */
654 if (EI_DYNAMICS(ir->eI))
658 mde_delta_h_coll_init(md->dhc, ir);
661 snew(md->dE, ir->fepvals->n_lambda);
665 md->fp_dhdl = fp_dhdl;
666 snew(md->dE, ir->fepvals->n_lambda);
671 snew(md->temperatures, ir->fepvals->n_lambda);
672 for (i = 0; i < ir->fepvals->n_lambda; i++)
674 md->temperatures[i] = ir->simtempvals->temperatures[i];
680 /* print a lambda vector to a string
681 fep = the inputrec's FEP input data
682 i = the index of the lambda vector
683 get_native_lambda = whether to print the native lambda
684 get_names = whether to print the names rather than the values
685 str = the pre-allocated string buffer to print to. */
686 static void print_lambda_vector(t_lambda *fep, int i,
687 gmx_bool get_native_lambda, gmx_bool get_names,
693 for (j = 0; j < efptNR; j++)
695 if (fep->separate_dvdl[j])
700 str[0] = 0; /* reset the string */
703 str += sprintf(str, "("); /* set the opening parenthesis*/
705 for (j = 0; j < efptNR; j++)
707 if (fep->separate_dvdl[j])
711 if (get_native_lambda && fep->init_lambda >= 0)
713 str += sprintf(str, "%.4f", fep->init_lambda);
717 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
722 str += sprintf(str, "%s", efpt_singular_names[j]);
724 /* print comma for the next item */
727 str += sprintf(str, ", ");
734 /* and add the closing parenthesis */
740 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
741 const output_env_t oenv)
744 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
745 *lambdastate = "\\lambda state";
746 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
747 int i, nps, nsets, nsets_de, nsetsbegin;
748 int n_lambda_terms = 0;
749 t_lambda *fep = ir->fepvals; /* for simplicity */
750 t_expanded *expand = ir->expandedvals;
752 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
758 gmx_bool write_pV = FALSE;
760 /* count the number of different lambda terms */
761 for (i = 0; i < efptNR; i++)
763 if (fep->separate_dvdl[i])
769 if (fep->n_lambda == 0)
771 sprintf(title, "%s", dhdl);
772 sprintf(label_x, "Time (ps)");
773 sprintf(label_y, "%s (%s %s)",
774 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
778 sprintf(title, "%s and %s", dhdl, deltag);
779 sprintf(label_x, "Time (ps)");
780 sprintf(label_y, "%s and %s (%s %s)",
781 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
783 fp = gmx_fio_fopen(filename, "w+");
784 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
788 bufplace = sprintf(buf, "T = %g (K) ",
791 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
793 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
795 /* compatibility output */
796 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
800 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
802 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
804 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
805 lambdastate, fep->init_fep_state,
806 lambda_name_str, lambda_vec_str);
809 xvgr_subtitle(fp, buf, oenv);
813 if (fep->dhdl_derivatives == edhdlderivativesYES)
815 nsets_dhdl = n_lambda_terms;
817 /* count the number of delta_g states */
818 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
820 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
822 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
824 nsets += 1; /*add fep state for expanded ensemble */
827 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
829 nsets += 1; /* add energy to the dhdl as well */
833 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
835 nsetsextend += 1; /* for PV term, other terms possible if required for
836 the reduced potential (only needed with foreign
837 lambda, and only output when init_lambda is not
838 set in order to maintain compatibility of the
842 snew(setname, nsetsextend);
844 if (expand->elmcmove > elmcmoveNO)
846 /* state for the fep_vals, if we have alchemical sampling */
847 sprintf(buf, "%s", "Thermodynamic state");
848 setname[s] = gmx_strdup(buf);
852 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
854 switch (fep->edHdLPrintEnergy)
856 case edHdLPrintEnergyPOTENTIAL:
857 sprintf(buf, "%s (%s)", "Potential Energy", unit_energy);
859 case edHdLPrintEnergyTOTAL:
860 case edHdLPrintEnergyYES:
862 sprintf(buf, "%s (%s)", "Total Energy", unit_energy);
864 setname[s] = gmx_strdup(buf);
868 if (fep->dhdl_derivatives == edhdlderivativesYES)
870 for (i = 0; i < efptNR; i++)
872 if (fep->separate_dvdl[i])
875 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
877 /* compatibility output */
878 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
882 double lam = fep->init_lambda;
883 if (fep->init_lambda < 0)
885 lam = fep->all_lambda[i][fep->init_fep_state];
887 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
890 setname[s] = gmx_strdup(buf);
896 if (fep->n_lambda > 0)
898 /* g_bar has to determine the lambda values used in this simulation
899 * from this xvg legend.
902 if (expand->elmcmove > elmcmoveNO)
904 nsetsbegin = 1; /* for including the expanded ensemble */
911 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
915 nsetsbegin += nsets_dhdl;
917 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
919 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
920 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
922 /* for compatible dhdl.xvg files */
923 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
927 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
932 /* print the temperature for this state if doing simulated annealing */
933 sprintf(&buf[nps], "T = %g (%s)",
934 ir->simtempvals->temperatures[s-(nsetsbegin)],
937 setname[s] = gmx_strdup(buf);
942 sprintf(buf, "pV (%s)", unit_energy);
943 setname[nsetsextend-1] = gmx_strdup(buf); /* the first entry after
947 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
949 for (s = 0; s < nsetsextend; s++)
959 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
963 for (i = j = 0; (i < F_NRE); i++)
972 gmx_incons("Number of energy terms wrong");
976 void upd_mdebin(t_mdebin *md,
981 gmx_enerdata_t *enerd,
990 gmx_ekindata_t *ekind,
994 int i, j, k, kk, n, gid;
995 real crmsd[2], tmp6[6];
996 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
999 double store_dhdl[efptNR];
1000 real store_energy = 0;
1003 /* Do NOT use the box in the state variable, but the separate box provided
1004 * as an argument. This is because we sometimes need to write the box from
1005 * the last timestep to match the trajectory frames.
1007 copy_energy(md, enerd->term, ecopy);
1008 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1011 crmsd[0] = constr_rmsd(constr, FALSE);
1014 crmsd[1] = constr_rmsd(constr, TRUE);
1016 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1023 bs[0] = box[XX][XX];
1024 bs[1] = box[YY][YY];
1025 bs[2] = box[ZZ][ZZ];
1026 bs[3] = box[YY][XX];
1027 bs[4] = box[ZZ][XX];
1028 bs[5] = box[ZZ][YY];
1033 bs[0] = box[XX][XX];
1034 bs[1] = box[YY][YY];
1035 bs[2] = box[ZZ][ZZ];
1038 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1039 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1040 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1041 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1042 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1046 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1047 not the instantaneous pressure */
1048 pv = vol*md->ref_p/PRESFAC;
1050 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1051 enthalpy = pv + enerd->term[F_ETOT];
1052 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1057 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1058 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1062 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1066 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1070 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1071 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1073 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1075 tmp6[0] = state->boxv[XX][XX];
1076 tmp6[1] = state->boxv[YY][YY];
1077 tmp6[2] = state->boxv[ZZ][ZZ];
1078 tmp6[3] = state->boxv[YY][XX];
1079 tmp6[4] = state->boxv[ZZ][XX];
1080 tmp6[5] = state->boxv[ZZ][YY];
1081 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1085 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1087 if (ekind && ekind->cosacc.cos_accel != 0)
1089 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1090 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1091 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1092 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1093 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1094 *dens*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1095 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1100 for (i = 0; (i < md->nEg); i++)
1102 for (j = i; (j < md->nEg); j++)
1104 gid = GID(i, j, md->nEg);
1105 for (k = kk = 0; (k < egNR); k++)
1109 eee[kk++] = enerd->grpp.ener[k][gid];
1112 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1120 for (i = 0; (i < md->nTC); i++)
1122 md->tmp_r[i] = ekind->tcstat[i].T;
1124 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1126 if (md->etc == etcNOSEHOOVER)
1128 /* whether to print Nose-Hoover chains: */
1129 if (md->bPrintNHChains)
1131 if (md->bNHC_trotter)
1133 for (i = 0; (i < md->nTC); i++)
1135 for (j = 0; j < md->nNHC; j++)
1138 md->tmp_r[2*k] = state->nosehoover_xi[k];
1139 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1142 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1146 for (i = 0; (i < md->nTCP); i++)
1148 for (j = 0; j < md->nNHC; j++)
1151 md->tmp_r[2*k] = state->nhpres_xi[k];
1152 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1155 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1160 for (i = 0; (i < md->nTC); i++)
1162 md->tmp_r[2*i] = state->nosehoover_xi[i];
1163 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1165 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1169 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1170 md->etc == etcVRESCALE)
1172 for (i = 0; (i < md->nTC); i++)
1174 md->tmp_r[i] = ekind->tcstat[i].lambda;
1176 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1180 if (ekind && md->nU > 1)
1182 for (i = 0; (i < md->nU); i++)
1184 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1186 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1189 ebin_increase_count(md->ebin, bSum);
1191 /* BAR + thermodynamic integration values */
1192 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1194 for (i = 0; i < enerd->n_lambda-1; i++)
1196 /* zero for simulated tempering */
1197 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1198 if (md->temperatures != NULL)
1200 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1201 /* is this even useful to have at all? */
1202 md->dE[i] += (md->temperatures[i]/
1203 md->temperatures[state->fep_state]-1.0)*
1204 enerd->term[F_EKIN];
1210 fprintf(md->fp_dhdl, "%.4f", time);
1211 /* the current free energy state */
1213 /* print the current state if we are doing expanded ensemble */
1214 if (expand->elmcmove > elmcmoveNO)
1216 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1218 /* total energy (for if the temperature changes */
1220 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1222 switch (fep->edHdLPrintEnergy)
1224 case edHdLPrintEnergyPOTENTIAL:
1225 store_energy = enerd->term[F_EPOT];
1227 case edHdLPrintEnergyTOTAL:
1228 case edHdLPrintEnergyYES:
1230 store_energy = enerd->term[F_ETOT];
1232 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1235 if (fep->dhdl_derivatives == edhdlderivativesYES)
1237 for (i = 0; i < efptNR; i++)
1239 if (fep->separate_dvdl[i])
1241 /* assumes F_DVDL is first */
1242 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1246 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1248 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1252 (md->epc != epcNO) &&
1253 (enerd->n_lambda > 0) &&
1254 (fep->init_lambda < 0))
1256 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1257 there are alternate state
1258 lambda and we're not in
1259 compatibility mode */
1261 fprintf(md->fp_dhdl, "\n");
1262 /* and the binary free energy output */
1264 if (md->dhc && bDoDHDL)
1267 for (i = 0; i < efptNR; i++)
1269 if (fep->separate_dvdl[i])
1271 /* assumes F_DVDL is first */
1272 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1276 store_energy = enerd->term[F_ETOT];
1277 /* store_dh is dE */
1278 mde_delta_h_coll_add_dh(md->dhc,
1279 (double)state->fep_state,
1283 md->dE + fep->lambda_start_n,
1290 void upd_mdebin_step(t_mdebin *md)
1292 ebin_increase_count(md->ebin, FALSE);
1295 static void npr(FILE *log, int n, char c)
1297 for (; (n > 0); n--)
1299 fprintf(log, "%c", c);
1303 static void pprint(FILE *log, const char *s, t_mdebin *md)
1307 char buf1[22], buf2[22];
1310 fprintf(log, "\t<====== ");
1311 npr(log, slen, CHAR);
1312 fprintf(log, " ==>\n");
1313 fprintf(log, "\t<==== %s ====>\n", s);
1314 fprintf(log, "\t<== ");
1315 npr(log, slen, CHAR);
1316 fprintf(log, " ======>\n\n");
1318 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1319 gmx_step_str(md->ebin->nsteps_sim, buf1),
1320 gmx_step_str(md->ebin->nsum_sim, buf2));
1324 void print_ebin_header(FILE *log, gmx_int64_t steps, double time, real lambda)
1328 fprintf(log, " %12s %12s %12s\n"
1329 " %12s %12.5f %12.5f\n\n",
1330 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1333 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1335 gmx_int64_t step, double time,
1336 int mode, gmx_bool bCompact,
1337 t_mdebin *md, t_fcdata *fcd,
1338 gmx_groups_t *groups, t_grpopts *opts)
1340 /*static char **grpnms=NULL;*/
1342 int i, j, n, ni, nj, b;
1344 real *disre_rm3tav, *disre_rt;
1346 /* these are for the old-style blocks (1 subblock, only reals), because
1347 there can be only one per ID for these */
1360 fr.nsteps = md->ebin->nsteps;
1361 fr.dt = md->delta_t;
1362 fr.nsum = md->ebin->nsum;
1363 fr.nre = (bEne) ? md->ebin->nener : 0;
1364 fr.ener = md->ebin->e;
1365 ndisre = bDR ? fcd->disres.npair : 0;
1366 disre_rm3tav = fcd->disres.rm3tav;
1367 disre_rt = fcd->disres.rt;
1368 /* Optional additional old-style (real-only) blocks. */
1369 for (i = 0; i < enxNR; i++)
1373 if (fcd->orires.nr > 0 && bOR)
1375 diagonalize_orires_tensors(&(fcd->orires));
1376 nr[enxOR] = fcd->orires.nr;
1377 block[enxOR] = fcd->orires.otav;
1379 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1381 block[enxORI] = fcd->orires.oinsl;
1382 id[enxORI] = enxORI;
1383 nr[enxORT] = fcd->orires.nex*12;
1384 block[enxORT] = fcd->orires.eig;
1385 id[enxORT] = enxORT;
1388 /* whether we are going to wrte anything out: */
1389 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1392 /* the old-style blocks go first */
1394 for (i = 0; i < enxNR; i++)
1401 add_blocks_enxframe(&fr, fr.nblock);
1402 for (b = 0; b < fr.nblock; b++)
1404 add_subblocks_enxblock(&(fr.block[b]), 1);
1405 fr.block[b].id = id[b];
1406 fr.block[b].sub[0].nr = nr[b];
1408 fr.block[b].sub[0].type = xdr_datatype_float;
1409 fr.block[b].sub[0].fval = block[b];
1411 fr.block[b].sub[0].type = xdr_datatype_double;
1412 fr.block[b].sub[0].dval = block[b];
1416 /* check for disre block & fill it. */
1421 add_blocks_enxframe(&fr, fr.nblock);
1423 add_subblocks_enxblock(&(fr.block[db]), 2);
1424 fr.block[db].id = enxDISRE;
1425 fr.block[db].sub[0].nr = ndisre;
1426 fr.block[db].sub[1].nr = ndisre;
1428 fr.block[db].sub[0].type = xdr_datatype_float;
1429 fr.block[db].sub[1].type = xdr_datatype_float;
1430 fr.block[db].sub[0].fval = disre_rt;
1431 fr.block[db].sub[1].fval = disre_rm3tav;
1433 fr.block[db].sub[0].type = xdr_datatype_double;
1434 fr.block[db].sub[1].type = xdr_datatype_double;
1435 fr.block[db].sub[0].dval = disre_rt;
1436 fr.block[db].sub[1].dval = disre_rm3tav;
1439 /* here we can put new-style blocks */
1441 /* Free energy perturbation blocks */
1444 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1447 /* we can now free & reset the data in the blocks */
1450 mde_delta_h_coll_reset(md->dhc);
1453 /* do the actual I/O */
1454 do_enx(fp_ene, &fr);
1457 /* We have stored the sums, so reset the sum history */
1458 reset_ebin_sums(md->ebin);
1466 pprint(log, "A V E R A G E S", md);
1472 pprint(log, "R M S - F L U C T U A T I O N S", md);
1476 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1481 for (i = 0; i < opts->ngtc; i++)
1483 if (opts->annealing[i] != eannNO)
1485 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1486 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1490 if (mode == eprNORMAL && fcd->orires.nr > 0)
1492 print_orires_log(log, &(fcd->orires));
1494 fprintf(log, " Energies (%s)\n", unit_energy);
1495 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1502 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1508 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1509 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1511 fprintf(log, " Force Virial (%s)\n", unit_energy);
1512 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1517 fprintf(log, " Total Virial (%s)\n", unit_energy);
1518 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1523 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1524 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1529 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1530 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1536 if (md->print_grpnms == NULL)
1538 snew(md->print_grpnms, md->nE);
1540 for (i = 0; (i < md->nEg); i++)
1542 ni = groups->grps[egcENER].nm_ind[i];
1543 for (j = i; (j < md->nEg); j++)
1545 nj = groups->grps[egcENER].nm_ind[j];
1546 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1547 *(groups->grpname[nj]));
1548 md->print_grpnms[n++] = gmx_strdup(buf);
1552 sprintf(buf, "Epot (%s)", unit_energy);
1553 fprintf(log, "%15s ", buf);
1554 for (i = 0; (i < egNR); i++)
1558 fprintf(log, "%12s ", egrp_nm[i]);
1562 for (i = 0; (i < md->nE); i++)
1564 fprintf(log, "%15s", md->print_grpnms[i]);
1565 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1572 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1577 fprintf(log, "%15s %12s %12s %12s\n",
1578 "Group", "Ux", "Uy", "Uz");
1579 for (i = 0; (i < md->nU); i++)
1581 ni = groups->grps[egcACC].nm_ind[i];
1582 fprintf(log, "%15s", *groups->grpname[ni]);
1583 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1592 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1596 enerhist->nsteps = mdebin->ebin->nsteps;
1597 enerhist->nsum = mdebin->ebin->nsum;
1598 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1599 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1600 enerhist->nener = mdebin->ebin->nener;
1602 if (mdebin->ebin->nsum > 0)
1604 /* Check if we need to allocate first */
1605 if (enerhist->ener_ave == NULL)
1607 snew(enerhist->ener_ave, enerhist->nener);
1608 snew(enerhist->ener_sum, enerhist->nener);
1611 for (i = 0; i < enerhist->nener; i++)
1613 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1614 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1618 if (mdebin->ebin->nsum_sim > 0)
1620 /* Check if we need to allocate first */
1621 if (enerhist->ener_sum_sim == NULL)
1623 snew(enerhist->ener_sum_sim, enerhist->nener);
1626 for (i = 0; i < enerhist->nener; i++)
1628 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1633 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1637 void restore_energyhistory_from_state(t_mdebin * mdebin,
1638 energyhistory_t * enerhist)
1642 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1643 mdebin->ebin->nener != enerhist->nener)
1645 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1646 mdebin->ebin->nener, enerhist->nener);
1649 mdebin->ebin->nsteps = enerhist->nsteps;
1650 mdebin->ebin->nsum = enerhist->nsum;
1651 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1652 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1654 for (i = 0; i < mdebin->ebin->nener; i++)
1656 mdebin->ebin->e[i].eav =
1657 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1658 mdebin->ebin->e[i].esum =
1659 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1660 mdebin->ebin->e_sim[i].esum =
1661 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1665 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);