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,2016,2017,2018,2019, 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 "energyoutput.h"
47 #include "gromacs/awh/awh.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/listed_forces/disre.h"
53 #include "gromacs/listed_forces/orires.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/constr.h"
58 #include "gromacs/mdlib/ebin.h"
59 #include "gromacs/mdlib/mdebin_bar.h"
60 #include "gromacs/mdtypes/energyhistory.h"
61 #include "gromacs/mdtypes/fcdata.h"
62 #include "gromacs/mdtypes/group.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pulling/pull.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/trajectory/energyframe.h"
70 #include "gromacs/utility/arraysize.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
75 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
77 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
79 static const char *tricl_boxs_nm[] = {
80 "Box-XX", "Box-YY", "Box-ZZ",
81 "Box-YX", "Box-ZX", "Box-ZY"
84 static const char *vol_nm[] = { "Volume" };
86 static const char *dens_nm[] = {"Density" };
88 static const char *pv_nm[] = {"pV" };
90 static const char *enthalpy_nm[] = {"Enthalpy" };
92 static const char *boxvel_nm[] = {
93 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
94 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
97 #define NBOXS asize(boxs_nm)
98 #define NTRICLBOXS asize(tricl_boxs_nm)
100 const char *egrp_nm[egNR+1] = {
101 "Coul-SR", "LJ-SR", "Buck-SR",
102 "Coul-14", "LJ-14", nullptr
105 /* forward declaration */
106 typedef struct t_mde_delta_h_coll t_mde_delta_h_coll;
114 /* This is the collection of energy averages collected during mdrun, and to
115 be written out to the .edr file. */
120 int ie, iconrmsd, ib, ivol, idens, ipv, ienthalpy;
121 int isvir, ifvir, ipres, ivir, isurft, ipc, itemp, itc, itcb, iu, imu;
123 int nE, nEg, nEc, nTC, nTCP, nU, nNHC;
132 gmx_bool bNHC_trotter;
133 gmx_bool bPrintNHChains;
135 gmx_bool bMu; /* true if dipole is calculated */
143 gmx_bool bEner[F_NRE];
144 gmx_bool bEInd[egNR];
147 FILE *fp_dhdl; /* the dhdl.xvg output file */
148 double *dE; /* energy components for dhdl.xvg output */
149 t_mde_delta_h_coll *dhc; /* the delta U components (raw data + histogram) */
153 } // namespace detail
155 using detail::t_mdebin;
160 //! Legacy init function
161 t_mdebin *init_mdebin(ener_file_t fp_ene,
162 const gmx_mtop_t *mtop,
163 const t_inputrec *ir,
167 const char *ener_nm[F_NRE];
168 static const char *vir_nm[] = {
169 "Vir-XX", "Vir-XY", "Vir-XZ",
170 "Vir-YX", "Vir-YY", "Vir-YZ",
171 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
173 static const char *sv_nm[] = {
174 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
175 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
176 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
178 static const char *fv_nm[] = {
179 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
180 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
181 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
183 static const char *pres_nm[] = {
184 "Pres-XX", "Pres-XY", "Pres-XZ",
185 "Pres-YX", "Pres-YY", "Pres-YZ",
186 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
188 static const char *surft_nm[] = {
191 static const char *mu_nm[] = {
192 "Mu-X", "Mu-Y", "Mu-Z"
194 static const char *vcos_nm[] = {
197 static const char *visc_nm[] = {
200 static const char *baro_nm[] = {
204 const gmx_groups_t *groups;
209 int i, j, ni, nj, n, k, kk, ncon, nset;
214 if (EI_DYNAMICS(ir->eI))
216 md->delta_t = ir->delta_t;
223 groups = &mtop->groups;
225 bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
226 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
227 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
229 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
230 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
231 md->bConstr = (ncon > 0 || nset > 0) && !isRerun;
232 md->bConstrVir = FALSE;
235 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
239 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != nullptr);
246 /* Energy monitoring */
247 for (i = 0; i < egNR; i++)
249 md->bEInd[i] = FALSE;
252 for (i = 0; i < F_NRE; i++)
254 md->bEner[i] = FALSE;
256 (i == F_EKIN || i == F_ETOT || i == F_ECONSERVED ||
257 i == F_TEMP || i == F_PDISPCORR || i == F_PRES))
263 md->bEner[i] = !bBHAM;
265 else if (i == F_BHAM)
267 md->bEner[i] = bBHAM;
271 md->bEner[i] = ir->bQMMM;
273 else if (i == F_RF_EXCL)
275 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
277 else if (i == F_COUL_RECIP)
279 md->bEner[i] = EEL_FULL(ir->coulombtype);
281 else if (i == F_LJ_RECIP)
283 md->bEner[i] = EVDW_PME(ir->vdwtype);
285 else if (i == F_LJ14)
289 else if (i == F_COUL14)
293 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
295 md->bEner[i] = FALSE;
297 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
298 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
299 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
300 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
301 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
302 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
304 md->bEner[i] = (ir->efep != efepNO);
306 else if ((interaction_function[i].flags & IF_VSITE) ||
307 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
309 md->bEner[i] = FALSE;
311 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
315 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
317 md->bEner[i] = EI_DYNAMICS(ir->eI);
319 else if (i == F_DISPCORR || i == F_PDISPCORR)
321 md->bEner[i] = (ir->eDispCorr != edispcNO);
323 else if (i == F_DISRESVIOL)
325 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
327 else if (i == F_ORIRESDEV)
329 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
331 else if (i == F_CONNBONDS)
333 md->bEner[i] = FALSE;
335 else if (i == F_COM_PULL)
337 md->bEner[i] = ((ir->bPull && pull_have_potential(ir->pull_work)) ||
340 else if (i == F_ECONSERVED)
342 md->bEner[i] = (integratorHasConservedEnergyQuantity(ir));
346 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
351 for (i = 0; i < F_NRE; i++)
355 ener_nm[md->f_nre] = interaction_function[i].longname;
360 md->epc = isRerun ? epcNO : ir->epc;
361 md->bDiagPres = !TRICLINIC(ir->ref_p) && !isRerun;
362 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
363 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
364 md->bDynBox = inputrecDynamicBox(ir);
365 md->etc = isRerun ? etcNO : ir->etc;
366 md->bNHC_trotter = inputrecNvtTrotter(ir) && !isRerun;
367 md->bPrintNHChains = ir->bPrintNHChains && !isRerun;
368 md->bMTTK = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
369 md->bMu = inputrecNeedMutot(ir);
370 md->bPres = !isRerun;
372 md->ebin = mk_ebin();
373 /* Pass NULL for unit to let get_ebin_space determine the units
374 * for interaction_function[i].longname
376 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, nullptr);
379 /* This should be called directly after the call for md->ie,
380 * such that md->iconrmsd follows directly in the list.
382 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
386 md->ib = get_ebin_space(md->ebin,
387 md->bTricl ? NTRICLBOXS : NBOXS,
388 md->bTricl ? tricl_boxs_nm : boxs_nm,
390 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
391 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
394 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
395 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
400 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
401 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
405 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
406 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
407 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
410 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
412 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
413 boxvel_nm, unit_vel);
417 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
419 if (ir->cos_accel != 0)
421 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
422 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
426 /* Energy monitoring */
427 for (i = 0; i < egNR; i++)
429 md->bEInd[i] = FALSE;
431 md->bEInd[egCOULSR] = TRUE;
432 md->bEInd[egLJSR ] = TRUE;
436 md->bEInd[egLJSR] = FALSE;
437 md->bEInd[egBHAMSR] = TRUE;
441 md->bEInd[egLJ14] = TRUE;
442 md->bEInd[egCOUL14] = TRUE;
445 for (i = 0; (i < egNR); i++)
453 n = groups->grps[egcENER].nr;
455 md->nE = (n*(n+1))/2;
457 snew(md->igrp, md->nE);
462 for (k = 0; (k < md->nEc); k++)
464 snew(gnm[k], STRLEN);
466 for (i = 0; (i < groups->grps[egcENER].nr); i++)
468 ni = groups->grps[egcENER].nm_ind[i];
469 for (j = i; (j < groups->grps[egcENER].nr); j++)
471 nj = groups->grps[egcENER].nm_ind[j];
472 for (k = kk = 0; (k < egNR); k++)
476 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
477 *(groups->grpname[ni]), *(groups->grpname[nj]));
481 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
486 for (k = 0; (k < md->nEc); k++)
494 gmx_incons("Number of energy terms wrong");
498 md->nTC = isRerun ? 0 : groups->grps[egcTC].nr;
499 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
502 md->nTCP = 1; /* assume only one possible coupling system for barostat
509 if (md->etc == etcNOSEHOOVER)
511 if (md->bNHC_trotter)
513 md->mde_n = 2*md->nNHC*md->nTC;
517 md->mde_n = 2*md->nTC;
519 if (md->epc == epcMTTK)
521 md->mdeb_n = 2*md->nNHC*md->nTCP;
530 snew(md->tmp_r, md->mde_n);
531 snew(md->tmp_v, md->mde_n);
533 snew(grpnms, md->mde_n);
535 for (i = 0; (i < md->nTC); i++)
537 ni = groups->grps[egcTC].nm_ind[i];
538 sprintf(buf, "T-%s", *(groups->grpname[ni]));
539 grpnms[i] = gmx_strdup(buf);
541 md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
544 if (md->etc == etcNOSEHOOVER)
546 if (md->bPrintNHChains)
548 if (md->bNHC_trotter)
550 for (i = 0; (i < md->nTC); i++)
552 ni = groups->grps[egcTC].nm_ind[i];
553 bufi = *(groups->grpname[ni]);
554 for (j = 0; (j < md->nNHC); j++)
556 sprintf(buf, "Xi-%d-%s", j, bufi);
557 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
558 sprintf(buf, "vXi-%d-%s", j, bufi);
559 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
562 md->itc = get_ebin_space(md->ebin, md->mde_n,
563 grpnms, unit_invtime);
566 for (i = 0; (i < md->nTCP); i++)
568 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
569 for (j = 0; (j < md->nNHC); j++)
571 sprintf(buf, "Xi-%d-%s", j, bufi);
572 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
573 sprintf(buf, "vXi-%d-%s", j, bufi);
574 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
577 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
578 grpnms, unit_invtime);
583 for (i = 0; (i < md->nTC); i++)
585 ni = groups->grps[egcTC].nm_ind[i];
586 bufi = *(groups->grpname[ni]);
587 sprintf(buf, "Xi-%s", bufi);
588 grpnms[2*i] = gmx_strdup(buf);
589 sprintf(buf, "vXi-%s", bufi);
590 grpnms[2*i+1] = gmx_strdup(buf);
592 md->itc = get_ebin_space(md->ebin, md->mde_n,
593 grpnms, unit_invtime);
597 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
598 md->etc == etcVRESCALE)
600 for (i = 0; (i < md->nTC); i++)
602 ni = groups->grps[egcTC].nm_ind[i];
603 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
604 grpnms[i] = gmx_strdup(buf);
606 md->itc = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
609 for (i = 0; i < md->mde_n; i++)
615 md->nU = groups->grps[egcACC].nr;
618 snew(grpnms, 3*md->nU);
619 for (i = 0; (i < md->nU); i++)
621 ni = groups->grps[egcACC].nm_ind[i];
622 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
623 grpnms[3*i+XX] = gmx_strdup(buf);
624 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
625 grpnms[3*i+YY] = gmx_strdup(buf);
626 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
627 grpnms[3*i+ZZ] = gmx_strdup(buf);
629 md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
635 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
638 md->print_grpnms = nullptr;
640 /* check whether we're going to write dh histograms */
642 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
644 /* Currently dh histograms are only written with dynamics */
645 if (EI_DYNAMICS(ir->eI))
649 mde_delta_h_coll_init(md->dhc, ir);
651 md->fp_dhdl = nullptr;
652 snew(md->dE, ir->fepvals->n_lambda);
656 md->fp_dhdl = fp_dhdl;
657 snew(md->dE, ir->fepvals->n_lambda);
662 snew(md->temperatures, ir->fepvals->n_lambda);
663 for (i = 0; i < ir->fepvals->n_lambda; i++)
665 md->temperatures[i] = ir->simtempvals->temperatures[i];
671 //! Legacy cleanup function
672 void done_mdebin(t_mdebin *mdebin)
675 sfree(mdebin->tmp_r);
676 sfree(mdebin->tmp_v);
677 done_ebin(mdebin->ebin);
678 done_mde_delta_h_coll(mdebin->dhc);
680 sfree(mdebin->temperatures);
687 /* print a lambda vector to a string
688 fep = the inputrec's FEP input data
689 i = the index of the lambda vector
690 get_native_lambda = whether to print the native lambda
691 get_names = whether to print the names rather than the values
692 str = the pre-allocated string buffer to print to. */
693 static void print_lambda_vector(t_lambda *fep, int i,
694 gmx_bool get_native_lambda, gmx_bool get_names,
700 for (j = 0; j < efptNR; j++)
702 if (fep->separate_dvdl[j])
707 str[0] = 0; /* reset the string */
710 str += sprintf(str, "("); /* set the opening parenthesis*/
712 for (j = 0; j < efptNR; j++)
714 if (fep->separate_dvdl[j])
718 if (get_native_lambda && fep->init_lambda >= 0)
720 str += sprintf(str, "%.4f", fep->init_lambda);
724 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
729 str += sprintf(str, "%s", efpt_singular_names[j]);
731 /* print comma for the next item */
734 str += sprintf(str, ", ");
741 /* and add the closing parenthesis */
746 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
747 const gmx_output_env_t *oenv)
750 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
751 *lambdastate = "\\lambda state";
752 int i, nsets, nsets_de, nsetsbegin;
753 int n_lambda_terms = 0;
754 t_lambda *fep = ir->fepvals; /* for simplicity */
755 t_expanded *expand = ir->expandedvals;
756 char 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 std::string title, label_x, label_y;
773 if (fep->n_lambda == 0)
775 title = gmx::formatString("%s", dhdl);
776 label_x = gmx::formatString("Time (ps)");
777 label_y = gmx::formatString("%s (%s %s)",
778 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
782 title = gmx::formatString("%s and %s", dhdl, deltag);
783 label_x = gmx::formatString("Time (ps)");
784 label_y = gmx::formatString("%s and %s (%s %s)",
785 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
787 fp = gmx_fio_fopen(filename, "w+");
788 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
793 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
795 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
797 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
799 /* compatibility output */
800 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
804 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
806 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
808 buf += gmx::formatString("%s %d: %s = %s",
809 lambdastate, fep->init_fep_state,
810 lambda_name_str, lambda_vec_str);
813 xvgr_subtitle(fp, buf.c_str(), oenv);
817 if (fep->dhdl_derivatives == edhdlderivativesYES)
819 nsets_dhdl = n_lambda_terms;
821 /* count the number of delta_g states */
822 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
824 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
826 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
828 nsets += 1; /*add fep state for expanded ensemble */
831 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
833 nsets += 1; /* add energy to the dhdl as well */
837 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
839 nsetsextend += 1; /* for PV term, other terms possible if required for
840 the reduced potential (only needed with foreign
841 lambda, and only output when init_lambda is not
842 set in order to maintain compatibility of the
846 std::vector<std::string> setname(nsetsextend);
848 if (expand->elmcmove > elmcmoveNO)
850 /* state for the fep_vals, if we have alchemical sampling */
851 setname[s++] = "Thermodynamic state";
854 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
857 switch (fep->edHdLPrintEnergy)
859 case edHdLPrintEnergyPOTENTIAL:
860 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
862 case edHdLPrintEnergyTOTAL:
863 case edHdLPrintEnergyYES:
865 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
867 setname[s++] = energy;
870 if (fep->dhdl_derivatives == edhdlderivativesYES)
872 for (i = 0; i < efptNR; i++)
874 if (fep->separate_dvdl[i])
876 std::string derivative;
877 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
879 /* compatibility output */
880 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
884 double lam = fep->init_lambda;
885 if (fep->init_lambda < 0)
887 lam = fep->all_lambda[i][fep->init_fep_state];
889 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
892 setname[s++] = derivative;
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);
922 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
924 /* for compatible dhdl.xvg files */
925 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
929 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
934 /* print the temperature for this state if doing simulated annealing */
935 buf += gmx::formatString("T = %g (%s)",
936 ir->simtempvals->temperatures[s-(nsetsbegin)],
943 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
946 xvgrLegend(fp, setname, oenv);
957 //! Legacy update function
958 void upd_mdebin(t_mdebin *md,
963 gmx_enerdata_t *enerd,
972 gmx_ekindata_t *ekind,
974 const gmx::Constraints *constr)
976 int i, j, k, kk, n, gid;
977 real crmsd[2], tmp6[6];
978 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
980 double store_dhdl[efptNR];
981 real store_energy = 0;
984 /* Do NOT use the box in the state variable, but the separate box provided
985 * as an argument. This is because we sometimes need to write the box from
986 * the last timestep to match the trajectory frames.
988 add_ebin_indexed(md->ebin, md->ie, gmx::ArrayRef<gmx_bool>(md->bEner), enerd->term, bSum);
991 crmsd[0] = constr->rmsd();
992 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1000 bs[1] = box[YY][YY];
1001 bs[2] = box[ZZ][ZZ];
1002 bs[3] = box[YY][XX];
1003 bs[4] = box[ZZ][XX];
1004 bs[5] = box[ZZ][YY];
1009 bs[0] = box[XX][XX];
1010 bs[1] = box[YY][YY];
1011 bs[2] = box[ZZ][ZZ];
1014 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1015 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1016 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1017 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1018 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1022 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1023 not the instantaneous pressure */
1024 pv = vol*md->ref_p/PRESFAC;
1026 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1027 enthalpy = pv + enerd->term[F_ETOT];
1028 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1033 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1034 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1038 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1039 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1040 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1041 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1043 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1045 tmp6[0] = state->boxv[XX][XX];
1046 tmp6[1] = state->boxv[YY][YY];
1047 tmp6[2] = state->boxv[ZZ][ZZ];
1048 tmp6[3] = state->boxv[YY][XX];
1049 tmp6[4] = state->boxv[ZZ][XX];
1050 tmp6[5] = state->boxv[ZZ][YY];
1051 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1055 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1057 if (ekind && ekind->cosacc.cos_accel != 0)
1059 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1060 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1061 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1062 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1063 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1064 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1065 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1070 for (i = 0; (i < md->nEg); i++)
1072 for (j = i; (j < md->nEg); j++)
1074 gid = GID(i, j, md->nEg);
1075 for (k = kk = 0; (k < egNR); k++)
1079 eee[kk++] = enerd->grpp.ener[k][gid];
1082 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1090 for (i = 0; (i < md->nTC); i++)
1092 md->tmp_r[i] = ekind->tcstat[i].T;
1094 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1096 if (md->etc == etcNOSEHOOVER)
1098 /* whether to print Nose-Hoover chains: */
1099 if (md->bPrintNHChains)
1101 if (md->bNHC_trotter)
1103 for (i = 0; (i < md->nTC); i++)
1105 for (j = 0; j < md->nNHC; j++)
1108 md->tmp_r[2*k] = state->nosehoover_xi[k];
1109 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1112 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1116 for (i = 0; (i < md->nTCP); i++)
1118 for (j = 0; j < md->nNHC; j++)
1121 md->tmp_r[2*k] = state->nhpres_xi[k];
1122 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1125 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1130 for (i = 0; (i < md->nTC); i++)
1132 md->tmp_r[2*i] = state->nosehoover_xi[i];
1133 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1135 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1139 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1140 md->etc == etcVRESCALE)
1142 for (i = 0; (i < md->nTC); i++)
1144 md->tmp_r[i] = ekind->tcstat[i].lambda;
1146 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1150 if (ekind && md->nU > 1)
1152 for (i = 0; (i < md->nU); i++)
1154 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1156 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1159 ebin_increase_count(md->ebin, bSum);
1161 /* BAR + thermodynamic integration values */
1162 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1164 for (i = 0; i < enerd->n_lambda-1; i++)
1166 /* zero for simulated tempering */
1167 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1168 if (md->temperatures != nullptr)
1170 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1171 /* is this even useful to have at all? */
1172 md->dE[i] += (md->temperatures[i]/
1173 md->temperatures[state->fep_state]-1.0)*
1174 enerd->term[F_EKIN];
1180 fprintf(md->fp_dhdl, "%.4f", time);
1181 /* the current free energy state */
1183 /* print the current state if we are doing expanded ensemble */
1184 if (expand->elmcmove > elmcmoveNO)
1186 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1188 /* total energy (for if the temperature changes */
1190 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1192 switch (fep->edHdLPrintEnergy)
1194 case edHdLPrintEnergyPOTENTIAL:
1195 store_energy = enerd->term[F_EPOT];
1197 case edHdLPrintEnergyTOTAL:
1198 case edHdLPrintEnergyYES:
1200 store_energy = enerd->term[F_ETOT];
1202 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1205 if (fep->dhdl_derivatives == edhdlderivativesYES)
1207 for (i = 0; i < efptNR; i++)
1209 if (fep->separate_dvdl[i])
1211 /* assumes F_DVDL is first */
1212 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1216 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1218 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1222 (md->epc != epcNO) &&
1223 (enerd->n_lambda > 0) &&
1224 (fep->init_lambda < 0))
1226 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1227 there are alternate state
1228 lambda and we're not in
1229 compatibility mode */
1231 fprintf(md->fp_dhdl, "\n");
1232 /* and the binary free energy output */
1234 if (md->dhc && bDoDHDL)
1237 for (i = 0; i < efptNR; i++)
1239 if (fep->separate_dvdl[i])
1241 /* assumes F_DVDL is first */
1242 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1246 store_energy = enerd->term[F_ETOT];
1247 /* store_dh is dE */
1248 mde_delta_h_coll_add_dh(md->dhc,
1249 static_cast<double>(state->fep_state),
1253 md->dE + fep->lambda_start_n,
1261 void EnergyOutput::recordNonEnergyStep()
1263 ebin_increase_count(mdebin->ebin, false);
1269 //! Legacy output function
1270 void npr(FILE *log, int n, char c)
1272 for (; (n > 0); n--)
1274 fprintf(log, "%c", c);
1278 //! Legacy output function
1279 void pprint(FILE *log, const char *s, t_mdebin *md)
1283 char buf1[22], buf2[22];
1286 fprintf(log, "\t<====== ");
1287 npr(log, slen, CHAR);
1288 fprintf(log, " ==>\n");
1289 fprintf(log, "\t<==== %s ====>\n", s);
1290 fprintf(log, "\t<== ");
1291 npr(log, slen, CHAR);
1292 fprintf(log, " ======>\n\n");
1294 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1295 gmx_step_str(md->ebin->nsteps_sim, buf1),
1296 gmx_step_str(md->ebin->nsum_sim, buf2));
1302 void print_ebin_header(FILE *log, int64_t steps, double time)
1306 fprintf(log, " %12s %12s\n"
1308 "Step", "Time", gmx_step_str(steps, buf), time);
1314 // TODO It is too many responsibilities for this function to handle
1315 // both .edr and .log output for both per-time and time-average data.
1316 //! Legacy ebin output function
1317 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1319 int64_t step, double time,
1321 t_mdebin *md, t_fcdata *fcd,
1322 gmx_groups_t *groups, t_grpopts *opts,
1325 /*static char **grpnms=NULL;*/
1327 int i, j, n, ni, nj, b;
1330 /* these are for the old-style blocks (1 subblock, only reals), because
1331 there can be only one per ID for these */
1338 if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1342 fprintf(log, "Not enough data recorded to report energy averages\n");
1353 fr.nsteps = md->ebin->nsteps;
1354 fr.dt = md->delta_t;
1355 fr.nsum = md->ebin->nsum;
1356 fr.nre = (bEne) ? md->ebin->nener : 0;
1357 fr.ener = md->ebin->e;
1358 ndisre = bDR ? fcd->disres.npair : 0;
1359 /* Optional additional old-style (real-only) blocks. */
1360 for (i = 0; i < enxNR; i++)
1364 if (bOR && fcd->orires.nr > 0)
1366 diagonalize_orires_tensors(&(fcd->orires));
1367 nr[enxOR] = fcd->orires.nr;
1368 block[enxOR] = fcd->orires.otav;
1370 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1372 block[enxORI] = fcd->orires.oinsl;
1373 id[enxORI] = enxORI;
1374 nr[enxORT] = fcd->orires.nex*12;
1375 block[enxORT] = fcd->orires.eig;
1376 id[enxORT] = enxORT;
1379 /* whether we are going to wrte anything out: */
1380 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1383 /* the old-style blocks go first */
1385 for (i = 0; i < enxNR; i++)
1392 add_blocks_enxframe(&fr, fr.nblock);
1393 for (b = 0; b < fr.nblock; b++)
1395 add_subblocks_enxblock(&(fr.block[b]), 1);
1396 fr.block[b].id = id[b];
1397 fr.block[b].sub[0].nr = nr[b];
1399 fr.block[b].sub[0].type = xdr_datatype_float;
1400 fr.block[b].sub[0].fval = block[b];
1402 fr.block[b].sub[0].type = xdr_datatype_double;
1403 fr.block[b].sub[0].dval = block[b];
1407 /* check for disre block & fill it. */
1412 add_blocks_enxframe(&fr, fr.nblock);
1414 add_subblocks_enxblock(&(fr.block[db]), 2);
1415 fr.block[db].id = enxDISRE;
1416 fr.block[db].sub[0].nr = ndisre;
1417 fr.block[db].sub[1].nr = ndisre;
1419 fr.block[db].sub[0].type = xdr_datatype_float;
1420 fr.block[db].sub[1].type = xdr_datatype_float;
1421 fr.block[db].sub[0].fval = fcd->disres.rt;
1422 fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1424 fr.block[db].sub[0].type = xdr_datatype_double;
1425 fr.block[db].sub[1].type = xdr_datatype_double;
1426 fr.block[db].sub[0].dval = fcd->disres.rt;
1427 fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1430 /* here we can put new-style blocks */
1432 /* Free energy perturbation blocks */
1435 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1438 /* we can now free & reset the data in the blocks */
1441 mde_delta_h_coll_reset(md->dhc);
1444 /* AWH bias blocks. */
1445 if (awh != nullptr) // TODO: add boolean in t_mdebin. See in mdebin.h.
1447 awh->writeToEnergyFrame(step, &fr);
1450 /* do the actual I/O */
1451 do_enx(fp_ene, &fr);
1454 /* We have stored the sums, so reset the sum history */
1455 reset_ebin_sums(md->ebin);
1463 pprint(log, "A V E R A G E S", md);
1469 pprint(log, "R M S - F L U C T U A T I O N S", md);
1473 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1480 for (i = 0; i < opts->ngtc; i++)
1482 if (opts->annealing[i] != eannNO)
1484 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1485 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1490 if (mode == eprNORMAL && bOR && 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);
1498 if (mode == eprAVER)
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);
1520 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1521 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1526 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1527 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1533 if (md->print_grpnms == nullptr)
1535 snew(md->print_grpnms, md->nE);
1537 for (i = 0; (i < md->nEg); i++)
1539 ni = groups->grps[egcENER].nm_ind[i];
1540 for (j = i; (j < md->nEg); j++)
1542 nj = groups->grps[egcENER].nm_ind[j];
1543 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1544 *(groups->grpname[nj]));
1545 md->print_grpnms[n++] = gmx_strdup(buf);
1549 sprintf(buf, "Epot (%s)", unit_energy);
1550 fprintf(log, "%15s ", buf);
1551 for (i = 0; (i < egNR); i++)
1555 fprintf(log, "%12s ", egrp_nm[i]);
1559 for (i = 0; (i < md->nE); i++)
1561 fprintf(log, "%15s", md->print_grpnms[i]);
1562 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1569 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1574 fprintf(log, "%15s %12s %12s %12s\n",
1575 "Group", "Ux", "Uy", "Uz");
1576 for (i = 0; (i < md->nU); i++)
1578 ni = groups->grps[egcACC].nm_ind[i];
1579 fprintf(log, "%15s", *groups->grpname[ni]);
1580 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1589 //! Legacy update function
1590 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1592 const t_ebin * const ebin = mdebin->ebin;
1594 enerhist->nsteps = ebin->nsteps;
1595 enerhist->nsum = ebin->nsum;
1596 enerhist->nsteps_sim = ebin->nsteps_sim;
1597 enerhist->nsum_sim = ebin->nsum_sim;
1601 /* This will only actually resize the first time */
1602 enerhist->ener_ave.resize(ebin->nener);
1603 enerhist->ener_sum.resize(ebin->nener);
1605 for (int i = 0; i < ebin->nener; i++)
1607 enerhist->ener_ave[i] = ebin->e[i].eav;
1608 enerhist->ener_sum[i] = ebin->e[i].esum;
1612 if (ebin->nsum_sim > 0)
1614 /* This will only actually resize the first time */
1615 enerhist->ener_sum_sim.resize(ebin->nener);
1617 for (int i = 0; i < ebin->nener; i++)
1619 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1624 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1628 //! Legacy restore function
1629 void restore_energyhistory_from_state(t_mdebin * mdebin,
1630 const energyhistory_t * enerhist)
1632 unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1634 GMX_RELEASE_ASSERT(enerhist, "Need valid history to restore");
1636 if ((enerhist->nsum > 0 && nener != enerhist->ener_sum.size()) ||
1637 (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1639 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1640 nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1643 mdebin->ebin->nsteps = enerhist->nsteps;
1644 mdebin->ebin->nsum = enerhist->nsum;
1645 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1646 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1648 for (int i = 0; i < mdebin->ebin->nener; i++)
1650 mdebin->ebin->e[i].eav =
1651 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1652 mdebin->ebin->e[i].esum =
1653 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1654 mdebin->ebin->e_sim[i].esum =
1655 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1659 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());
1665 EnergyOutput::EnergyOutput()
1670 void EnergyOutput::prepare(ener_file *fp_ene,
1671 const gmx_mtop_t *mtop,
1672 const t_inputrec *ir,
1676 mdebin = init_mdebin(fp_ene, mtop, ir, fp_dhdl, isRerun);
1679 EnergyOutput::~EnergyOutput()
1681 done_mdebin(mdebin);
1684 t_ebin *EnergyOutput::getEbin()
1686 return mdebin->ebin;
1689 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
1693 gmx_enerdata_t *enerd,
1702 gmx_ekindata_t *ekind,
1704 const gmx::Constraints *constr)
1706 upd_mdebin(mdebin, bDoDHDL, bSum, time, tmass, enerd, state, fep,
1707 expand, box, svir, fvir, vir, pres, ekind, mu_tot, constr);
1710 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1712 int64_t step, double time,
1715 gmx_groups_t *groups, t_grpopts *opts,
1718 print_ebin(fp_ene, bEne, bDR, bOR, log, step, time, mode,
1719 mdebin, fcd, groups, opts, awh);
1722 int EnergyOutput::numEnergyTerms() const
1724 return mdebin->ebin->nener;
1727 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1729 update_energyhistory(enerhist, mdebin);
1732 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1734 restore_energyhistory_from_state(mdebin, &enerhist);