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 SimulationGroups *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++)
452 n = groups->groups[SimulationAtomGroupType::EnergyOutput].nr;
454 md->nE = (n*(n+1))/2;
456 snew(md->igrp, md->nE);
461 for (k = 0; (k < md->nEc); k++)
463 snew(gnm[k], STRLEN);
465 for (i = 0; (i < groups->groups[SimulationAtomGroupType::EnergyOutput].nr); i++)
467 ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
468 for (j = i; (j < groups->groups[SimulationAtomGroupType::EnergyOutput].nr); j++)
470 nj = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[j];
471 for (k = kk = 0; (k < egNR); k++)
475 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
476 *(groups->groupNames[ni]), *(groups->groupNames[nj]));
480 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
485 for (k = 0; (k < md->nEc); k++)
493 gmx_incons("Number of energy terms wrong");
497 md->nTC = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].nr;
498 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
501 md->nTCP = 1; /* assume only one possible coupling system for barostat
508 if (md->etc == etcNOSEHOOVER)
510 if (md->bNHC_trotter)
512 md->mde_n = 2*md->nNHC*md->nTC;
516 md->mde_n = 2*md->nTC;
518 if (md->epc == epcMTTK)
520 md->mdeb_n = 2*md->nNHC*md->nTCP;
529 snew(md->tmp_r, md->mde_n);
530 snew(md->tmp_v, md->mde_n);
532 snew(grpnms, md->mde_n);
534 for (i = 0; (i < md->nTC); i++)
536 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
537 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
538 grpnms[i] = gmx_strdup(buf);
540 md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
543 if (md->etc == etcNOSEHOOVER)
545 if (md->bPrintNHChains)
547 if (md->bNHC_trotter)
549 for (i = 0; (i < md->nTC); i++)
551 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
552 bufi = *(groups->groupNames[ni]);
553 for (j = 0; (j < md->nNHC); j++)
555 sprintf(buf, "Xi-%d-%s", j, bufi);
556 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
557 sprintf(buf, "vXi-%d-%s", j, bufi);
558 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
561 md->itc = get_ebin_space(md->ebin, md->mde_n,
562 grpnms, unit_invtime);
565 for (i = 0; (i < md->nTCP); i++)
567 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
568 for (j = 0; (j < md->nNHC); j++)
570 sprintf(buf, "Xi-%d-%s", j, bufi);
571 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
572 sprintf(buf, "vXi-%d-%s", j, bufi);
573 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
576 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
577 grpnms, unit_invtime);
582 for (i = 0; (i < md->nTC); i++)
584 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
585 bufi = *(groups->groupNames[ni]);
586 sprintf(buf, "Xi-%s", bufi);
587 grpnms[2*i] = gmx_strdup(buf);
588 sprintf(buf, "vXi-%s", bufi);
589 grpnms[2*i+1] = gmx_strdup(buf);
591 md->itc = get_ebin_space(md->ebin, md->mde_n,
592 grpnms, unit_invtime);
596 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
597 md->etc == etcVRESCALE)
599 for (i = 0; (i < md->nTC); i++)
601 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
602 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
603 grpnms[i] = gmx_strdup(buf);
605 md->itc = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
608 for (i = 0; i < md->mde_n; i++)
614 md->nU = groups->groups[SimulationAtomGroupType::Acceleration].nr;
617 snew(grpnms, 3*md->nU);
618 for (i = 0; (i < md->nU); i++)
620 ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
621 sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
622 grpnms[3*i+XX] = gmx_strdup(buf);
623 sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
624 grpnms[3*i+YY] = gmx_strdup(buf);
625 sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
626 grpnms[3*i+ZZ] = gmx_strdup(buf);
628 md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
634 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
637 md->print_grpnms = nullptr;
639 /* check whether we're going to write dh histograms */
641 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
643 /* Currently dh histograms are only written with dynamics */
644 if (EI_DYNAMICS(ir->eI))
648 mde_delta_h_coll_init(md->dhc, ir);
650 md->fp_dhdl = nullptr;
651 snew(md->dE, ir->fepvals->n_lambda);
655 md->fp_dhdl = fp_dhdl;
656 snew(md->dE, ir->fepvals->n_lambda);
661 snew(md->temperatures, ir->fepvals->n_lambda);
662 for (i = 0; i < ir->fepvals->n_lambda; i++)
664 md->temperatures[i] = ir->simtempvals->temperatures[i];
670 //! Legacy cleanup function
671 void done_mdebin(t_mdebin *mdebin)
674 sfree(mdebin->tmp_r);
675 sfree(mdebin->tmp_v);
676 done_ebin(mdebin->ebin);
677 done_mde_delta_h_coll(mdebin->dhc);
679 sfree(mdebin->temperatures);
686 /* print a lambda vector to a string
687 fep = the inputrec's FEP input data
688 i = the index of the lambda vector
689 get_native_lambda = whether to print the native lambda
690 get_names = whether to print the names rather than the values
691 str = the pre-allocated string buffer to print to. */
692 static void print_lambda_vector(t_lambda *fep, int i,
693 gmx_bool get_native_lambda, gmx_bool get_names,
699 for (j = 0; j < efptNR; j++)
701 if (fep->separate_dvdl[j])
706 str[0] = 0; /* reset the string */
709 str += sprintf(str, "("); /* set the opening parenthesis*/
711 for (j = 0; j < efptNR; j++)
713 if (fep->separate_dvdl[j])
717 if (get_native_lambda && fep->init_lambda >= 0)
719 str += sprintf(str, "%.4f", fep->init_lambda);
723 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
728 str += sprintf(str, "%s", efpt_singular_names[j]);
730 /* print comma for the next item */
733 str += sprintf(str, ", ");
740 /* and add the closing parenthesis */
745 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
746 const gmx_output_env_t *oenv)
749 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
750 *lambdastate = "\\lambda state";
751 int i, nsets, nsets_de, nsetsbegin;
752 int n_lambda_terms = 0;
753 t_lambda *fep = ir->fepvals; /* for simplicity */
754 t_expanded *expand = ir->expandedvals;
755 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
760 gmx_bool write_pV = FALSE;
762 /* count the number of different lambda terms */
763 for (i = 0; i < efptNR; i++)
765 if (fep->separate_dvdl[i])
771 std::string title, label_x, label_y;
772 if (fep->n_lambda == 0)
774 title = gmx::formatString("%s", dhdl);
775 label_x = gmx::formatString("Time (ps)");
776 label_y = gmx::formatString("%s (%s %s)",
777 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
781 title = gmx::formatString("%s and %s", dhdl, deltag);
782 label_x = gmx::formatString("Time (ps)");
783 label_y = gmx::formatString("%s and %s (%s %s)",
784 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
786 fp = gmx_fio_fopen(filename, "w+");
787 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
792 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
794 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
796 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
798 /* compatibility output */
799 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
803 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
805 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
807 buf += gmx::formatString("%s %d: %s = %s",
808 lambdastate, fep->init_fep_state,
809 lambda_name_str, lambda_vec_str);
812 xvgr_subtitle(fp, buf.c_str(), oenv);
816 if (fep->dhdl_derivatives == edhdlderivativesYES)
818 nsets_dhdl = n_lambda_terms;
820 /* count the number of delta_g states */
821 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
823 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
825 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
827 nsets += 1; /*add fep state for expanded ensemble */
830 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
832 nsets += 1; /* add energy to the dhdl as well */
836 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
838 nsetsextend += 1; /* for PV term, other terms possible if required for
839 the reduced potential (only needed with foreign
840 lambda, and only output when init_lambda is not
841 set in order to maintain compatibility of the
845 std::vector<std::string> setname(nsetsextend);
847 if (expand->elmcmove > elmcmoveNO)
849 /* state for the fep_vals, if we have alchemical sampling */
850 setname[s++] = "Thermodynamic state";
853 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
856 switch (fep->edHdLPrintEnergy)
858 case edHdLPrintEnergyPOTENTIAL:
859 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
861 case edHdLPrintEnergyTOTAL:
862 case edHdLPrintEnergyYES:
864 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
866 setname[s++] = energy;
869 if (fep->dhdl_derivatives == edhdlderivativesYES)
871 for (i = 0; i < efptNR; i++)
873 if (fep->separate_dvdl[i])
875 std::string derivative;
876 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
878 /* compatibility output */
879 derivative = gmx::formatString("%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 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
891 setname[s++] = derivative;
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);
921 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
923 /* for compatible dhdl.xvg files */
924 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
928 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
933 /* print the temperature for this state if doing simulated annealing */
934 buf += gmx::formatString("T = %g (%s)",
935 ir->simtempvals->temperatures[s-(nsetsbegin)],
942 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
945 xvgrLegend(fp, setname, oenv);
956 //! Legacy update function
957 void upd_mdebin(t_mdebin *md,
962 gmx_enerdata_t *enerd,
971 gmx_ekindata_t *ekind,
973 const gmx::Constraints *constr)
975 int i, j, k, kk, n, gid;
976 real crmsd[2], tmp6[6];
977 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
979 double store_dhdl[efptNR];
980 real store_energy = 0;
983 /* Do NOT use the box in the state variable, but the separate box provided
984 * as an argument. This is because we sometimes need to write the box from
985 * the last timestep to match the trajectory frames.
987 add_ebin_indexed(md->ebin, md->ie, gmx::ArrayRef<gmx_bool>(md->bEner), enerd->term, bSum);
990 crmsd[0] = constr->rmsd();
991 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1000 bs[2] = box[ZZ][ZZ];
1001 bs[3] = box[YY][XX];
1002 bs[4] = box[ZZ][XX];
1003 bs[5] = box[ZZ][YY];
1008 bs[0] = box[XX][XX];
1009 bs[1] = box[YY][YY];
1010 bs[2] = box[ZZ][ZZ];
1013 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1014 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1015 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1016 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1017 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1021 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1022 not the instantaneous pressure */
1023 pv = vol*md->ref_p/PRESFAC;
1025 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1026 enthalpy = pv + enerd->term[F_ETOT];
1027 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1032 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1033 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1037 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1038 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1039 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1040 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1042 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1044 tmp6[0] = state->boxv[XX][XX];
1045 tmp6[1] = state->boxv[YY][YY];
1046 tmp6[2] = state->boxv[ZZ][ZZ];
1047 tmp6[3] = state->boxv[YY][XX];
1048 tmp6[4] = state->boxv[ZZ][XX];
1049 tmp6[5] = state->boxv[ZZ][YY];
1050 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1054 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1056 if (ekind && ekind->cosacc.cos_accel != 0)
1058 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1059 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1060 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1061 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1062 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1063 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1064 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1069 for (i = 0; (i < md->nEg); i++)
1071 for (j = i; (j < md->nEg); j++)
1073 gid = GID(i, j, md->nEg);
1074 for (k = kk = 0; (k < egNR); k++)
1078 eee[kk++] = enerd->grpp.ener[k][gid];
1081 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1089 for (i = 0; (i < md->nTC); i++)
1091 md->tmp_r[i] = ekind->tcstat[i].T;
1093 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1095 if (md->etc == etcNOSEHOOVER)
1097 /* whether to print Nose-Hoover chains: */
1098 if (md->bPrintNHChains)
1100 if (md->bNHC_trotter)
1102 for (i = 0; (i < md->nTC); i++)
1104 for (j = 0; j < md->nNHC; j++)
1107 md->tmp_r[2*k] = state->nosehoover_xi[k];
1108 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1111 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1115 for (i = 0; (i < md->nTCP); i++)
1117 for (j = 0; j < md->nNHC; j++)
1120 md->tmp_r[2*k] = state->nhpres_xi[k];
1121 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1124 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1129 for (i = 0; (i < md->nTC); i++)
1131 md->tmp_r[2*i] = state->nosehoover_xi[i];
1132 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1134 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1138 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1139 md->etc == etcVRESCALE)
1141 for (i = 0; (i < md->nTC); i++)
1143 md->tmp_r[i] = ekind->tcstat[i].lambda;
1145 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1149 if (ekind && md->nU > 1)
1151 for (i = 0; (i < md->nU); i++)
1153 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1155 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1158 ebin_increase_count(md->ebin, bSum);
1160 /* BAR + thermodynamic integration values */
1161 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1163 for (i = 0; i < enerd->n_lambda-1; i++)
1165 /* zero for simulated tempering */
1166 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1167 if (md->temperatures != nullptr)
1169 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1170 /* is this even useful to have at all? */
1171 md->dE[i] += (md->temperatures[i]/
1172 md->temperatures[state->fep_state]-1.0)*
1173 enerd->term[F_EKIN];
1179 fprintf(md->fp_dhdl, "%.4f", time);
1180 /* the current free energy state */
1182 /* print the current state if we are doing expanded ensemble */
1183 if (expand->elmcmove > elmcmoveNO)
1185 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1187 /* total energy (for if the temperature changes */
1189 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1191 switch (fep->edHdLPrintEnergy)
1193 case edHdLPrintEnergyPOTENTIAL:
1194 store_energy = enerd->term[F_EPOT];
1196 case edHdLPrintEnergyTOTAL:
1197 case edHdLPrintEnergyYES:
1199 store_energy = enerd->term[F_ETOT];
1201 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1204 if (fep->dhdl_derivatives == edhdlderivativesYES)
1206 for (i = 0; i < efptNR; i++)
1208 if (fep->separate_dvdl[i])
1210 /* assumes F_DVDL is first */
1211 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1215 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1217 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1221 (md->epc != epcNO) &&
1222 (enerd->n_lambda > 0) &&
1223 (fep->init_lambda < 0))
1225 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1226 there are alternate state
1227 lambda and we're not in
1228 compatibility mode */
1230 fprintf(md->fp_dhdl, "\n");
1231 /* and the binary free energy output */
1233 if (md->dhc && bDoDHDL)
1236 for (i = 0; i < efptNR; i++)
1238 if (fep->separate_dvdl[i])
1240 /* assumes F_DVDL is first */
1241 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1245 store_energy = enerd->term[F_ETOT];
1246 /* store_dh is dE */
1247 mde_delta_h_coll_add_dh(md->dhc,
1248 static_cast<double>(state->fep_state),
1252 md->dE + fep->lambda_start_n,
1260 void EnergyOutput::recordNonEnergyStep()
1262 ebin_increase_count(mdebin->ebin, false);
1268 //! Legacy output function
1269 void npr(FILE *log, int n, char c)
1271 for (; (n > 0); n--)
1273 fprintf(log, "%c", c);
1277 //! Legacy output function
1278 void pprint(FILE *log, const char *s, t_mdebin *md)
1282 char buf1[22], buf2[22];
1285 fprintf(log, "\t<====== ");
1286 npr(log, slen, CHAR);
1287 fprintf(log, " ==>\n");
1288 fprintf(log, "\t<==== %s ====>\n", s);
1289 fprintf(log, "\t<== ");
1290 npr(log, slen, CHAR);
1291 fprintf(log, " ======>\n\n");
1293 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1294 gmx_step_str(md->ebin->nsteps_sim, buf1),
1295 gmx_step_str(md->ebin->nsum_sim, buf2));
1301 void print_ebin_header(FILE *log, int64_t steps, double time)
1305 fprintf(log, " %12s %12s\n"
1307 "Step", "Time", gmx_step_str(steps, buf), time);
1313 // TODO It is too many responsibilities for this function to handle
1314 // both .edr and .log output for both per-time and time-average data.
1315 //! Legacy ebin output function
1316 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1318 int64_t step, double time,
1320 t_mdebin *md, t_fcdata *fcd,
1321 SimulationGroups *groups, t_grpopts *opts,
1324 /*static char **grpnms=NULL;*/
1326 int i, j, n, ni, nj, b;
1329 /* these are for the old-style blocks (1 subblock, only reals), because
1330 there can be only one per ID for these */
1337 if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1341 fprintf(log, "Not enough data recorded to report energy averages\n");
1352 fr.nsteps = md->ebin->nsteps;
1353 fr.dt = md->delta_t;
1354 fr.nsum = md->ebin->nsum;
1355 fr.nre = (bEne) ? md->ebin->nener : 0;
1356 fr.ener = md->ebin->e;
1357 ndisre = bDR ? fcd->disres.npair : 0;
1358 /* Optional additional old-style (real-only) blocks. */
1359 for (i = 0; i < enxNR; i++)
1363 if (bOR && fcd->orires.nr > 0)
1365 diagonalize_orires_tensors(&(fcd->orires));
1366 nr[enxOR] = fcd->orires.nr;
1367 block[enxOR] = fcd->orires.otav;
1369 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1371 block[enxORI] = fcd->orires.oinsl;
1372 id[enxORI] = enxORI;
1373 nr[enxORT] = fcd->orires.nex*12;
1374 block[enxORT] = fcd->orires.eig;
1375 id[enxORT] = enxORT;
1378 /* whether we are going to wrte anything out: */
1379 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1382 /* the old-style blocks go first */
1384 for (i = 0; i < enxNR; i++)
1391 add_blocks_enxframe(&fr, fr.nblock);
1392 for (b = 0; b < fr.nblock; b++)
1394 add_subblocks_enxblock(&(fr.block[b]), 1);
1395 fr.block[b].id = id[b];
1396 fr.block[b].sub[0].nr = nr[b];
1398 fr.block[b].sub[0].type = xdr_datatype_float;
1399 fr.block[b].sub[0].fval = block[b];
1401 fr.block[b].sub[0].type = xdr_datatype_double;
1402 fr.block[b].sub[0].dval = block[b];
1406 /* check for disre block & fill it. */
1411 add_blocks_enxframe(&fr, fr.nblock);
1413 add_subblocks_enxblock(&(fr.block[db]), 2);
1414 fr.block[db].id = enxDISRE;
1415 fr.block[db].sub[0].nr = ndisre;
1416 fr.block[db].sub[1].nr = ndisre;
1418 fr.block[db].sub[0].type = xdr_datatype_float;
1419 fr.block[db].sub[1].type = xdr_datatype_float;
1420 fr.block[db].sub[0].fval = fcd->disres.rt;
1421 fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1423 fr.block[db].sub[0].type = xdr_datatype_double;
1424 fr.block[db].sub[1].type = xdr_datatype_double;
1425 fr.block[db].sub[0].dval = fcd->disres.rt;
1426 fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1429 /* here we can put new-style blocks */
1431 /* Free energy perturbation blocks */
1434 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1437 /* we can now free & reset the data in the blocks */
1440 mde_delta_h_coll_reset(md->dhc);
1443 /* AWH bias blocks. */
1444 if (awh != nullptr) // TODO: add boolean in t_mdebin. See in mdebin.h.
1446 awh->writeToEnergyFrame(step, &fr);
1449 /* do the actual I/O */
1450 do_enx(fp_ene, &fr);
1453 /* We have stored the sums, so reset the sum history */
1454 reset_ebin_sums(md->ebin);
1462 pprint(log, "A V E R A G E S", md);
1468 pprint(log, "R M S - F L U C T U A T I O N S", md);
1472 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1479 for (i = 0; i < opts->ngtc; i++)
1481 if (opts->annealing[i] != eannNO)
1483 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1484 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i]]),
1489 if (mode == eprNORMAL && bOR && fcd->orires.nr > 0)
1491 print_orires_log(log, &(fcd->orires));
1493 fprintf(log, " Energies (%s)\n", unit_energy);
1494 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1497 if (mode == eprAVER)
1501 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1507 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1508 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1510 fprintf(log, " Force Virial (%s)\n", unit_energy);
1511 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1516 fprintf(log, " Total Virial (%s)\n", unit_energy);
1517 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1519 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1520 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1525 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1526 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1532 if (md->print_grpnms == nullptr)
1534 snew(md->print_grpnms, md->nE);
1536 for (i = 0; (i < md->nEg); i++)
1538 ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
1539 for (j = i; (j < md->nEg); j++)
1541 nj = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[j];
1542 sprintf(buf, "%s-%s", *(groups->groupNames[ni]),
1543 *(groups->groupNames[nj]));
1544 md->print_grpnms[n++] = gmx_strdup(buf);
1548 sprintf(buf, "Epot (%s)", unit_energy);
1549 fprintf(log, "%15s ", buf);
1550 for (i = 0; (i < egNR); i++)
1554 fprintf(log, "%12s ", egrp_nm[i]);
1558 for (i = 0; (i < md->nE); i++)
1560 fprintf(log, "%15s", md->print_grpnms[i]);
1561 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1568 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1573 fprintf(log, "%15s %12s %12s %12s\n",
1574 "Group", "Ux", "Uy", "Uz");
1575 for (i = 0; (i < md->nU); i++)
1577 ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
1578 fprintf(log, "%15s", *groups->groupNames[ni]);
1579 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1588 //! Legacy update function
1589 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1591 const t_ebin * const ebin = mdebin->ebin;
1593 enerhist->nsteps = ebin->nsteps;
1594 enerhist->nsum = ebin->nsum;
1595 enerhist->nsteps_sim = ebin->nsteps_sim;
1596 enerhist->nsum_sim = ebin->nsum_sim;
1600 /* This will only actually resize the first time */
1601 enerhist->ener_ave.resize(ebin->nener);
1602 enerhist->ener_sum.resize(ebin->nener);
1604 for (int i = 0; i < ebin->nener; i++)
1606 enerhist->ener_ave[i] = ebin->e[i].eav;
1607 enerhist->ener_sum[i] = ebin->e[i].esum;
1611 if (ebin->nsum_sim > 0)
1613 /* This will only actually resize the first time */
1614 enerhist->ener_sum_sim.resize(ebin->nener);
1616 for (int i = 0; i < ebin->nener; i++)
1618 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1623 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1627 //! Legacy restore function
1628 void restore_energyhistory_from_state(t_mdebin * mdebin,
1629 const energyhistory_t * enerhist)
1631 unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1633 GMX_RELEASE_ASSERT(enerhist, "Need valid history to restore");
1635 if ((enerhist->nsum > 0 && nener != enerhist->ener_sum.size()) ||
1636 (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1638 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1639 nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1642 mdebin->ebin->nsteps = enerhist->nsteps;
1643 mdebin->ebin->nsum = enerhist->nsum;
1644 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1645 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1647 for (int i = 0; i < mdebin->ebin->nener; i++)
1649 mdebin->ebin->e[i].eav =
1650 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1651 mdebin->ebin->e[i].esum =
1652 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1653 mdebin->ebin->e_sim[i].esum =
1654 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1658 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());
1664 EnergyOutput::EnergyOutput()
1669 void EnergyOutput::prepare(ener_file *fp_ene,
1670 const gmx_mtop_t *mtop,
1671 const t_inputrec *ir,
1675 mdebin = init_mdebin(fp_ene, mtop, ir, fp_dhdl, isRerun);
1678 EnergyOutput::~EnergyOutput()
1680 done_mdebin(mdebin);
1683 t_ebin *EnergyOutput::getEbin()
1685 return mdebin->ebin;
1688 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
1692 gmx_enerdata_t *enerd,
1701 gmx_ekindata_t *ekind,
1703 const gmx::Constraints *constr)
1705 upd_mdebin(mdebin, bDoDHDL, bSum, time, tmass, enerd, state, fep,
1706 expand, box, svir, fvir, vir, pres, ekind, mu_tot, constr);
1709 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1711 int64_t step, double time,
1714 SimulationGroups *groups, t_grpopts *opts,
1717 print_ebin(fp_ene, bEne, bDR, bOR, log, step, time, mode,
1718 mdebin, fcd, groups, opts, awh);
1721 int EnergyOutput::numEnergyTerms() const
1723 return mdebin->ebin->nener;
1726 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1728 update_energyhistory(enerhist, mdebin);
1731 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1733 restore_energyhistory_from_state(mdebin, &enerhist);