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.
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/mdebin.h"
45 #include "gromacs/math/units.h"
46 #include "gromacs/fileio/enxio.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/legacyheaders/disre.h"
49 #include "gromacs/legacyheaders/network.h"
50 #include "gromacs/legacyheaders/names.h"
51 #include "gromacs/legacyheaders/orires.h"
52 #include "gromacs/legacyheaders/constr.h"
53 #include "gromacs/topology/mtop_util.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/legacyheaders/macros.h"
57 #include "gromacs/legacyheaders/mdrun.h"
58 #include "mdebin_bar.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/utility/smalloc.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_LJ_RECIP)
228 md->bEner[i] = EVDW_PME(ir->vdwtype);
230 else if (i == F_LJ14)
234 else if (i == F_COUL14)
238 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
240 md->bEner[i] = FALSE;
242 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
243 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
244 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
245 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
246 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
247 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
249 md->bEner[i] = (ir->efep != efepNO);
251 else if ((interaction_function[i].flags & IF_VSITE) ||
252 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
254 md->bEner[i] = FALSE;
256 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
260 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
264 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
268 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
270 md->bEner[i] = FALSE;
272 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
274 md->bEner[i] = EI_DYNAMICS(ir->eI);
276 else if (i == F_DISPCORR || i == F_PDISPCORR)
278 md->bEner[i] = (ir->eDispCorr != edispcNO);
280 else if (i == F_DISRESVIOL)
282 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
284 else if (i == F_ORIRESDEV)
286 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
288 else if (i == F_CONNBONDS)
290 md->bEner[i] = FALSE;
292 else if (i == F_COM_PULL)
294 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
296 else if (i == F_ECONSERVED)
298 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
299 (ir->epc == epcNO || ir->epc == epcMTTK));
303 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
307 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
308 if (ir->bAdress && !debug)
310 for (i = 0; i < F_NRE; i++)
312 md->bEner[i] = FALSE;
329 for (i = 0; i < F_NRE; i++)
333 ener_nm[md->f_nre] = interaction_function[i].longname;
339 md->bDiagPres = !TRICLINIC(ir->ref_p);
340 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
341 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
342 md->bDynBox = DYNAMIC_BOX(*ir);
344 md->bNHC_trotter = IR_NVT_TROTTER(ir);
345 md->bPrintNHChains = ir->bPrintNHChains;
346 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
347 md->bMu = NEED_MUTOT(*ir);
349 md->ebin = mk_ebin();
350 /* Pass NULL for unit to let get_ebin_space determine the units
351 * for interaction_function[i].longname
353 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
356 /* This should be called directly after the call for md->ie,
357 * such that md->iconrmsd follows directly in the list.
359 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
363 md->ib = get_ebin_space(md->ebin,
364 md->bTricl ? NTRICLBOXS : NBOXS,
365 md->bTricl ? tricl_boxs_nm : boxs_nm,
367 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
368 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
371 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
372 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
377 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
378 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
382 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
386 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
390 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
393 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
395 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
396 boxvel_nm, unit_vel);
400 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
402 if (ir->cos_accel != 0)
404 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
405 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
409 /* Energy monitoring */
410 for (i = 0; i < egNR; i++)
412 md->bEInd[i] = FALSE;
414 md->bEInd[egCOULSR] = TRUE;
415 md->bEInd[egLJSR ] = TRUE;
417 if (ir->rcoulomb > ir->rlist)
419 md->bEInd[egCOULLR] = TRUE;
423 if (ir->rvdw > ir->rlist)
425 md->bEInd[egLJLR] = TRUE;
430 md->bEInd[egLJSR] = FALSE;
431 md->bEInd[egBHAMSR] = TRUE;
432 if (ir->rvdw > ir->rlist)
434 md->bEInd[egBHAMLR] = TRUE;
439 md->bEInd[egLJ14] = TRUE;
440 md->bEInd[egCOUL14] = TRUE;
443 for (i = 0; (i < egNR); i++)
451 n = groups->grps[egcENER].nr;
452 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
455 /*standard simulation*/
457 md->nE = (n*(n+1))/2;
461 /*AdResS simulation*/
468 snew(md->igrp, md->nE);
473 for (k = 0; (k < md->nEc); k++)
475 snew(gnm[k], STRLEN);
477 for (i = 0; (i < groups->grps[egcENER].nr); i++)
479 ni = groups->grps[egcENER].nm_ind[i];
480 for (j = i; (j < groups->grps[egcENER].nr); j++)
482 nj = groups->grps[egcENER].nm_ind[j];
483 for (k = kk = 0; (k < egNR); k++)
487 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
488 *(groups->grpname[ni]), *(groups->grpname[nj]));
492 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
493 (const char **)gnm, unit_energy);
497 for (k = 0; (k < md->nEc); k++)
505 gmx_incons("Number of energy terms wrong");
509 md->nTC = groups->grps[egcTC].nr;
510 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
513 md->nTCP = 1; /* assume only one possible coupling system for barostat
520 if (md->etc == etcNOSEHOOVER)
522 if (md->bNHC_trotter)
524 md->mde_n = 2*md->nNHC*md->nTC;
528 md->mde_n = 2*md->nTC;
530 if (md->epc == epcMTTK)
532 md->mdeb_n = 2*md->nNHC*md->nTCP;
541 snew(md->tmp_r, md->mde_n);
542 snew(md->tmp_v, md->mde_n);
543 snew(md->grpnms, md->mde_n);
546 for (i = 0; (i < md->nTC); i++)
548 ni = groups->grps[egcTC].nm_ind[i];
549 sprintf(buf, "T-%s", *(groups->grpname[ni]));
550 grpnms[i] = gmx_strdup(buf);
552 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
555 if (md->etc == etcNOSEHOOVER)
557 if (md->bPrintNHChains)
559 if (md->bNHC_trotter)
561 for (i = 0; (i < md->nTC); i++)
563 ni = groups->grps[egcTC].nm_ind[i];
564 bufi = *(groups->grpname[ni]);
565 for (j = 0; (j < md->nNHC); j++)
567 sprintf(buf, "Xi-%d-%s", j, bufi);
568 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
569 sprintf(buf, "vXi-%d-%s", j, bufi);
570 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
573 md->itc = get_ebin_space(md->ebin, md->mde_n,
574 (const char **)grpnms, unit_invtime);
577 for (i = 0; (i < md->nTCP); i++)
579 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
580 for (j = 0; (j < md->nNHC); j++)
582 sprintf(buf, "Xi-%d-%s", j, bufi);
583 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
584 sprintf(buf, "vXi-%d-%s", j, bufi);
585 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
588 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
589 (const char **)grpnms, unit_invtime);
594 for (i = 0; (i < md->nTC); i++)
596 ni = groups->grps[egcTC].nm_ind[i];
597 bufi = *(groups->grpname[ni]);
598 sprintf(buf, "Xi-%s", bufi);
599 grpnms[2*i] = gmx_strdup(buf);
600 sprintf(buf, "vXi-%s", bufi);
601 grpnms[2*i+1] = gmx_strdup(buf);
603 md->itc = get_ebin_space(md->ebin, md->mde_n,
604 (const char **)grpnms, unit_invtime);
608 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
609 md->etc == etcVRESCALE)
611 for (i = 0; (i < md->nTC); i++)
613 ni = groups->grps[egcTC].nm_ind[i];
614 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
615 grpnms[i] = gmx_strdup(buf);
617 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
623 md->nU = groups->grps[egcACC].nr;
626 snew(grpnms, 3*md->nU);
627 for (i = 0; (i < md->nU); i++)
629 ni = groups->grps[egcACC].nm_ind[i];
630 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
631 grpnms[3*i+XX] = gmx_strdup(buf);
632 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
633 grpnms[3*i+YY] = gmx_strdup(buf);
634 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
635 grpnms[3*i+ZZ] = gmx_strdup(buf);
637 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
643 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
646 md->print_grpnms = NULL;
648 /* check whether we're going to write dh histograms */
650 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
652 /* Currently dh histograms are only written with dynamics */
653 if (EI_DYNAMICS(ir->eI))
657 mde_delta_h_coll_init(md->dhc, ir);
660 snew(md->dE, ir->fepvals->n_lambda);
664 md->fp_dhdl = fp_dhdl;
665 snew(md->dE, ir->fepvals->n_lambda);
670 snew(md->temperatures, ir->fepvals->n_lambda);
671 for (i = 0; i < ir->fepvals->n_lambda; i++)
673 md->temperatures[i] = ir->simtempvals->temperatures[i];
679 /* print a lambda vector to a string
680 fep = the inputrec's FEP input data
681 i = the index of the lambda vector
682 get_native_lambda = whether to print the native lambda
683 get_names = whether to print the names rather than the values
684 str = the pre-allocated string buffer to print to. */
685 static void print_lambda_vector(t_lambda *fep, int i,
686 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])
712 if (get_native_lambda && fep->init_lambda >= 0)
714 str += sprintf(str, "%.4f", fep->init_lambda);
718 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
723 str += sprintf(str, "%s", efpt_singular_names[j]);
725 /* print comma for the next item */
728 str += sprintf(str, ", ");
735 /* and add the closing parenthesis */
741 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
742 const output_env_t oenv)
745 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
746 *lambdastate = "\\lambda state", *remain = "remaining";
747 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
748 int i, np, nps, nsets, nsets_de, nsetsbegin;
749 int n_lambda_terms = 0;
750 t_lambda *fep = ir->fepvals; /* for simplicity */
751 t_expanded *expand = ir->expandedvals;
753 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
759 gmx_bool write_pV = FALSE;
761 /* count the number of different lambda terms */
762 for (i = 0; i < efptNR; i++)
764 if (fep->separate_dvdl[i])
770 if (fep->n_lambda == 0)
772 sprintf(title, "%s", dhdl);
773 sprintf(label_x, "Time (ps)");
774 sprintf(label_y, "%s (%s %s)",
775 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
779 sprintf(title, "%s and %s", dhdl, deltag);
780 sprintf(label_x, "Time (ps)");
781 sprintf(label_y, "%s and %s (%s %s)",
782 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
784 fp = gmx_fio_fopen(filename, "w+");
785 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
789 bufplace = sprintf(buf, "T = %g (K) ",
792 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
794 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
796 /* compatibility output */
797 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
801 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
803 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
805 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
806 lambdastate, fep->init_fep_state,
807 lambda_name_str, lambda_vec_str);
810 xvgr_subtitle(fp, buf, oenv);
814 if (fep->dhdl_derivatives == edhdlderivativesYES)
816 nsets_dhdl = n_lambda_terms;
818 /* count the number of delta_g states */
819 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
821 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
823 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
825 nsets += 1; /*add fep state for expanded ensemble */
828 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
830 nsets += 1; /* add energy to the dhdl as well */
834 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
836 nsetsextend += 1; /* for PV term, other terms possible if required for
837 the reduced potential (only needed with foreign
838 lambda, and only output when init_lambda is not
839 set in order to maintain compatibility of the
843 snew(setname, nsetsextend);
845 if (expand->elmcmove > elmcmoveNO)
847 /* state for the fep_vals, if we have alchemical sampling */
848 sprintf(buf, "%s", "Thermodynamic state");
849 setname[s] = gmx_strdup(buf);
853 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
855 switch (fep->edHdLPrintEnergy)
857 case edHdLPrintEnergyPOTENTIAL:
858 sprintf(buf, "%s (%s)", "Potential Energy", unit_energy);
860 case edHdLPrintEnergyTOTAL:
861 case edHdLPrintEnergyYES:
863 sprintf(buf, "%s (%s)", "Total Energy", unit_energy);
865 setname[s] = gmx_strdup(buf);
869 if (fep->dhdl_derivatives == edhdlderivativesYES)
871 for (i = 0; i < efptNR; i++)
873 if (fep->separate_dvdl[i])
876 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
878 /* compatibility output */
879 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
883 double lam = fep->init_lambda;
884 if (fep->init_lambda < 0)
886 lam = fep->all_lambda[i][fep->init_fep_state];
888 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
891 setname[s] = gmx_strdup(buf);
897 if (fep->n_lambda > 0)
899 /* g_bar has to determine the lambda values used in this simulation
900 * from this xvg legend.
903 if (expand->elmcmove > elmcmoveNO)
905 nsetsbegin = 1; /* for including the expanded ensemble */
912 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
916 nsetsbegin += nsets_dhdl;
918 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
920 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
921 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
923 /* for compatible dhdl.xvg files */
924 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
928 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
933 /* print the temperature for this state if doing simulated annealing */
934 sprintf(&buf[nps], "T = %g (%s)",
935 ir->simtempvals->temperatures[s-(nsetsbegin)],
938 setname[s] = gmx_strdup(buf);
943 np = sprintf(buf, "pV (%s)", unit_energy);
944 setname[nsetsextend-1] = gmx_strdup(buf); /* the first entry after
948 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
950 for (s = 0; s < nsetsextend; s++)
960 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
964 for (i = j = 0; (i < F_NRE); i++)
973 gmx_incons("Number of energy terms wrong");
977 void upd_mdebin(t_mdebin *md,
982 gmx_enerdata_t *enerd,
991 gmx_ekindata_t *ekind,
995 int i, j, k, kk, m, n, gid;
996 real crmsd[2], tmp6[6];
997 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
1000 double store_dhdl[efptNR];
1001 real store_energy = 0;
1004 /* Do NOT use the box in the state variable, but the separate box provided
1005 * as an argument. This is because we sometimes need to write the box from
1006 * the last timestep to match the trajectory frames.
1008 copy_energy(md, enerd->term, ecopy);
1009 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1012 crmsd[0] = constr_rmsd(constr, FALSE);
1015 crmsd[1] = constr_rmsd(constr, TRUE);
1017 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1024 bs[0] = box[XX][XX];
1025 bs[1] = box[YY][YY];
1026 bs[2] = box[ZZ][ZZ];
1027 bs[3] = box[YY][XX];
1028 bs[4] = box[ZZ][XX];
1029 bs[5] = box[ZZ][YY];
1034 bs[0] = box[XX][XX];
1035 bs[1] = box[YY][YY];
1036 bs[2] = box[ZZ][ZZ];
1039 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1040 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1041 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1042 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1043 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1047 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1048 not the instantaneous pressure */
1049 pv = vol*md->ref_p/PRESFAC;
1051 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1052 enthalpy = pv + enerd->term[F_ETOT];
1053 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1058 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1059 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1063 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1067 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1071 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1072 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1074 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1076 tmp6[0] = state->boxv[XX][XX];
1077 tmp6[1] = state->boxv[YY][YY];
1078 tmp6[2] = state->boxv[ZZ][ZZ];
1079 tmp6[3] = state->boxv[YY][XX];
1080 tmp6[4] = state->boxv[ZZ][XX];
1081 tmp6[5] = state->boxv[ZZ][YY];
1082 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1086 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1088 if (ekind && ekind->cosacc.cos_accel != 0)
1090 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1091 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1092 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1093 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1094 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1095 *dens*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1096 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1101 for (i = 0; (i < md->nEg); i++)
1103 for (j = i; (j < md->nEg); j++)
1105 gid = GID(i, j, md->nEg);
1106 for (k = kk = 0; (k < egNR); k++)
1110 eee[kk++] = enerd->grpp.ener[k][gid];
1113 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1121 for (i = 0; (i < md->nTC); i++)
1123 md->tmp_r[i] = ekind->tcstat[i].T;
1125 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1127 if (md->etc == etcNOSEHOOVER)
1129 /* whether to print Nose-Hoover chains: */
1130 if (md->bPrintNHChains)
1132 if (md->bNHC_trotter)
1134 for (i = 0; (i < md->nTC); i++)
1136 for (j = 0; j < md->nNHC; j++)
1139 md->tmp_r[2*k] = state->nosehoover_xi[k];
1140 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1143 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1147 for (i = 0; (i < md->nTCP); i++)
1149 for (j = 0; j < md->nNHC; j++)
1152 md->tmp_r[2*k] = state->nhpres_xi[k];
1153 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1156 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1161 for (i = 0; (i < md->nTC); i++)
1163 md->tmp_r[2*i] = state->nosehoover_xi[i];
1164 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1166 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1170 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1171 md->etc == etcVRESCALE)
1173 for (i = 0; (i < md->nTC); i++)
1175 md->tmp_r[i] = ekind->tcstat[i].lambda;
1177 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1181 if (ekind && md->nU > 1)
1183 for (i = 0; (i < md->nU); i++)
1185 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1187 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1190 ebin_increase_count(md->ebin, bSum);
1192 /* BAR + thermodynamic integration values */
1193 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1195 for (i = 0; i < enerd->n_lambda-1; i++)
1197 /* zero for simulated tempering */
1198 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1199 if (md->temperatures != NULL)
1201 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1202 /* is this even useful to have at all? */
1203 md->dE[i] += (md->temperatures[i]/
1204 md->temperatures[state->fep_state]-1.0)*
1205 enerd->term[F_EKIN];
1211 fprintf(md->fp_dhdl, "%.4f", time);
1212 /* the current free energy state */
1214 /* print the current state if we are doing expanded ensemble */
1215 if (expand->elmcmove > elmcmoveNO)
1217 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1219 /* total energy (for if the temperature changes */
1221 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1223 switch (fep->edHdLPrintEnergy)
1225 case edHdLPrintEnergyPOTENTIAL:
1226 store_energy = enerd->term[F_EPOT];
1228 case edHdLPrintEnergyTOTAL:
1229 case edHdLPrintEnergyYES:
1231 store_energy = enerd->term[F_ETOT];
1233 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1236 if (fep->dhdl_derivatives == edhdlderivativesYES)
1238 for (i = 0; i < efptNR; i++)
1240 if (fep->separate_dvdl[i])
1242 /* assumes F_DVDL is first */
1243 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1247 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1249 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1251 if ((md->epc != epcNO) &&
1252 (enerd->n_lambda > 0) &&
1253 (fep->init_lambda < 0))
1255 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1256 there are alternate state
1257 lambda and we're not in
1258 compatibility mode */
1260 fprintf(md->fp_dhdl, "\n");
1261 /* and the binary free energy output */
1263 if (md->dhc && bDoDHDL)
1266 for (i = 0; i < efptNR; i++)
1268 if (fep->separate_dvdl[i])
1270 /* assumes F_DVDL is first */
1271 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1275 store_energy = enerd->term[F_ETOT];
1276 /* store_dh is dE */
1277 mde_delta_h_coll_add_dh(md->dhc,
1278 (double)state->fep_state,
1282 md->dE + fep->lambda_start_n,
1289 void upd_mdebin_step(t_mdebin *md)
1291 ebin_increase_count(md->ebin, FALSE);
1294 static void npr(FILE *log, int n, char c)
1296 for (; (n > 0); n--)
1298 fprintf(log, "%c", c);
1302 static void pprint(FILE *log, const char *s, t_mdebin *md)
1306 char buf1[22], buf2[22];
1309 fprintf(log, "\t<====== ");
1310 npr(log, slen, CHAR);
1311 fprintf(log, " ==>\n");
1312 fprintf(log, "\t<==== %s ====>\n", s);
1313 fprintf(log, "\t<== ");
1314 npr(log, slen, CHAR);
1315 fprintf(log, " ======>\n\n");
1317 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1318 gmx_step_str(md->ebin->nsteps_sim, buf1),
1319 gmx_step_str(md->ebin->nsum_sim, buf2));
1323 void print_ebin_header(FILE *log, gmx_int64_t steps, double time, real lambda)
1327 fprintf(log, " %12s %12s %12s\n"
1328 " %12s %12.5f %12.5f\n\n",
1329 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1332 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1334 gmx_int64_t step, double time,
1335 int mode, gmx_bool bCompact,
1336 t_mdebin *md, t_fcdata *fcd,
1337 gmx_groups_t *groups, t_grpopts *opts)
1339 /*static char **grpnms=NULL;*/
1341 int i, j, n, ni, nj, ndr, nor, b;
1343 real *disre_rm3tav, *disre_rt;
1345 /* these are for the old-style blocks (1 subblock, only reals), because
1346 there can be only one per ID for these */
1351 /* temporary arrays for the lambda values to write out */
1352 double enxlambda_data[2];
1362 fr.nsteps = md->ebin->nsteps;
1363 fr.dt = md->delta_t;
1364 fr.nsum = md->ebin->nsum;
1365 fr.nre = (bEne) ? md->ebin->nener : 0;
1366 fr.ener = md->ebin->e;
1367 ndisre = bDR ? fcd->disres.npair : 0;
1368 disre_rm3tav = fcd->disres.rm3tav;
1369 disre_rt = fcd->disres.rt;
1370 /* Optional additional old-style (real-only) blocks. */
1371 for (i = 0; i < enxNR; i++)
1375 if (fcd->orires.nr > 0 && bOR)
1377 diagonalize_orires_tensors(&(fcd->orires));
1378 nr[enxOR] = fcd->orires.nr;
1379 block[enxOR] = fcd->orires.otav;
1381 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1383 block[enxORI] = fcd->orires.oinsl;
1384 id[enxORI] = enxORI;
1385 nr[enxORT] = fcd->orires.nex*12;
1386 block[enxORT] = fcd->orires.eig;
1387 id[enxORT] = enxORT;
1390 /* whether we are going to wrte anything out: */
1391 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1394 /* the old-style blocks go first */
1396 for (i = 0; i < enxNR; i++)
1403 add_blocks_enxframe(&fr, fr.nblock);
1404 for (b = 0; b < fr.nblock; b++)
1406 add_subblocks_enxblock(&(fr.block[b]), 1);
1407 fr.block[b].id = id[b];
1408 fr.block[b].sub[0].nr = nr[b];
1410 fr.block[b].sub[0].type = xdr_datatype_float;
1411 fr.block[b].sub[0].fval = block[b];
1413 fr.block[b].sub[0].type = xdr_datatype_double;
1414 fr.block[b].sub[0].dval = block[b];
1418 /* check for disre block & fill it. */
1423 add_blocks_enxframe(&fr, fr.nblock);
1425 add_subblocks_enxblock(&(fr.block[db]), 2);
1426 fr.block[db].id = enxDISRE;
1427 fr.block[db].sub[0].nr = ndisre;
1428 fr.block[db].sub[1].nr = ndisre;
1430 fr.block[db].sub[0].type = xdr_datatype_float;
1431 fr.block[db].sub[1].type = xdr_datatype_float;
1432 fr.block[db].sub[0].fval = disre_rt;
1433 fr.block[db].sub[1].fval = disre_rm3tav;
1435 fr.block[db].sub[0].type = xdr_datatype_double;
1436 fr.block[db].sub[1].type = xdr_datatype_double;
1437 fr.block[db].sub[0].dval = disre_rt;
1438 fr.block[db].sub[1].dval = disre_rm3tav;
1441 /* here we can put new-style blocks */
1443 /* Free energy perturbation blocks */
1446 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1449 /* we can now free & reset the data in the blocks */
1452 mde_delta_h_coll_reset(md->dhc);
1455 /* do the actual I/O */
1456 do_enx(fp_ene, &fr);
1459 /* We have stored the sums, so reset the sum history */
1460 reset_ebin_sums(md->ebin);
1468 pprint(log, "A V E R A G E S", md);
1474 pprint(log, "R M S - F L U C T U A T I O N S", md);
1478 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1483 for (i = 0; i < opts->ngtc; i++)
1485 if (opts->annealing[i] != eannNO)
1487 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1488 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1492 if (mode == eprNORMAL && fcd->orires.nr > 0)
1494 print_orires_log(log, &(fcd->orires));
1496 fprintf(log, " Energies (%s)\n", unit_energy);
1497 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1504 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1510 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1511 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1513 fprintf(log, " Force Virial (%s)\n", unit_energy);
1514 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1519 fprintf(log, " Total Virial (%s)\n", unit_energy);
1520 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1525 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1526 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1531 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1532 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1538 if (md->print_grpnms == NULL)
1540 snew(md->print_grpnms, md->nE);
1542 for (i = 0; (i < md->nEg); i++)
1544 ni = groups->grps[egcENER].nm_ind[i];
1545 for (j = i; (j < md->nEg); j++)
1547 nj = groups->grps[egcENER].nm_ind[j];
1548 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1549 *(groups->grpname[nj]));
1550 md->print_grpnms[n++] = gmx_strdup(buf);
1554 sprintf(buf, "Epot (%s)", unit_energy);
1555 fprintf(log, "%15s ", buf);
1556 for (i = 0; (i < egNR); i++)
1560 fprintf(log, "%12s ", egrp_nm[i]);
1564 for (i = 0; (i < md->nE); i++)
1566 fprintf(log, "%15s", md->print_grpnms[i]);
1567 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1574 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1579 fprintf(log, "%15s %12s %12s %12s\n",
1580 "Group", "Ux", "Uy", "Uz");
1581 for (i = 0; (i < md->nU); i++)
1583 ni = groups->grps[egcACC].nm_ind[i];
1584 fprintf(log, "%15s", *groups->grpname[ni]);
1585 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1594 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1598 enerhist->nsteps = mdebin->ebin->nsteps;
1599 enerhist->nsum = mdebin->ebin->nsum;
1600 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1601 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1602 enerhist->nener = mdebin->ebin->nener;
1604 if (mdebin->ebin->nsum > 0)
1606 /* Check if we need to allocate first */
1607 if (enerhist->ener_ave == NULL)
1609 snew(enerhist->ener_ave, enerhist->nener);
1610 snew(enerhist->ener_sum, enerhist->nener);
1613 for (i = 0; i < enerhist->nener; i++)
1615 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1616 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1620 if (mdebin->ebin->nsum_sim > 0)
1622 /* Check if we need to allocate first */
1623 if (enerhist->ener_sum_sim == NULL)
1625 snew(enerhist->ener_sum_sim, enerhist->nener);
1628 for (i = 0; i < enerhist->nener; i++)
1630 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1635 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1639 void restore_energyhistory_from_state(t_mdebin * mdebin,
1640 energyhistory_t * enerhist)
1644 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1645 mdebin->ebin->nener != enerhist->nener)
1647 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1648 mdebin->ebin->nener, enerhist->nener);
1651 mdebin->ebin->nsteps = enerhist->nsteps;
1652 mdebin->ebin->nsum = enerhist->nsum;
1653 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1654 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1656 for (i = 0; i < mdebin->ebin->nener; i++)
1658 mdebin->ebin->e[i].eav =
1659 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1660 mdebin->ebin->e[i].esum =
1661 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1662 mdebin->ebin->e_sim[i].esum =
1663 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1667 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);