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,
164 const pull_t *pull_work,
168 const char *ener_nm[F_NRE];
169 static const char *vir_nm[] = {
170 "Vir-XX", "Vir-XY", "Vir-XZ",
171 "Vir-YX", "Vir-YY", "Vir-YZ",
172 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
174 static const char *sv_nm[] = {
175 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
176 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
177 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
179 static const char *fv_nm[] = {
180 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
181 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
182 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
184 static const char *pres_nm[] = {
185 "Pres-XX", "Pres-XY", "Pres-XZ",
186 "Pres-YX", "Pres-YY", "Pres-YZ",
187 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
189 static const char *surft_nm[] = {
192 static const char *mu_nm[] = {
193 "Mu-X", "Mu-Y", "Mu-Z"
195 static const char *vcos_nm[] = {
198 static const char *visc_nm[] = {
201 static const char *baro_nm[] = {
205 const SimulationGroups *groups;
210 int i, j, ni, nj, n, k, kk, ncon, nset;
215 if (EI_DYNAMICS(ir->eI))
217 md->delta_t = ir->delta_t;
224 groups = &mtop->groups;
226 bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
227 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
228 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
230 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
231 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
232 md->bConstr = (ncon > 0 || nset > 0) && !isRerun;
233 md->bConstrVir = FALSE;
236 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
240 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != nullptr);
247 /* Energy monitoring */
248 for (i = 0; i < egNR; i++)
250 md->bEInd[i] = FALSE;
253 for (i = 0; i < F_NRE; i++)
255 md->bEner[i] = FALSE;
257 (i == F_EKIN || i == F_ETOT || i == F_ECONSERVED ||
258 i == F_TEMP || i == F_PDISPCORR || i == F_PRES))
264 md->bEner[i] = !bBHAM;
266 else if (i == F_BHAM)
268 md->bEner[i] = bBHAM;
272 md->bEner[i] = ir->bQMMM;
274 else if (i == F_RF_EXCL)
276 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
278 else if (i == F_COUL_RECIP)
280 md->bEner[i] = EEL_FULL(ir->coulombtype);
282 else if (i == F_LJ_RECIP)
284 md->bEner[i] = EVDW_PME(ir->vdwtype);
286 else if (i == F_LJ14)
290 else if (i == F_COUL14)
294 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
296 md->bEner[i] = FALSE;
298 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
299 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
300 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
301 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
302 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
303 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
305 md->bEner[i] = (ir->efep != efepNO);
307 else if ((interaction_function[i].flags & IF_VSITE) ||
308 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
310 md->bEner[i] = FALSE;
312 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
316 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
318 md->bEner[i] = EI_DYNAMICS(ir->eI);
320 else if (i == F_DISPCORR || i == F_PDISPCORR)
322 md->bEner[i] = (ir->eDispCorr != edispcNO);
324 else if (i == F_DISRESVIOL)
326 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
328 else if (i == F_ORIRESDEV)
330 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
332 else if (i == F_CONNBONDS)
334 md->bEner[i] = FALSE;
336 else if (i == F_COM_PULL)
338 md->bEner[i] = ((ir->bPull && pull_have_potential(pull_work)) ||
341 else if (i == F_ECONSERVED)
343 md->bEner[i] = (integratorHasConservedEnergyQuantity(ir));
347 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
352 for (i = 0; i < F_NRE; i++)
356 ener_nm[md->f_nre] = interaction_function[i].longname;
361 md->epc = isRerun ? epcNO : ir->epc;
362 md->bDiagPres = !TRICLINIC(ir->ref_p) && !isRerun;
363 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
364 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
365 md->bDynBox = inputrecDynamicBox(ir);
366 md->etc = isRerun ? etcNO : ir->etc;
367 md->bNHC_trotter = inputrecNvtTrotter(ir) && !isRerun;
368 md->bPrintNHChains = ir->bPrintNHChains && !isRerun;
369 md->bMTTK = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
370 md->bMu = inputrecNeedMutot(ir);
371 md->bPres = !isRerun;
373 md->ebin = mk_ebin();
374 /* Pass NULL for unit to let get_ebin_space determine the units
375 * for interaction_function[i].longname
377 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, nullptr);
380 /* This should be called directly after the call for md->ie,
381 * such that md->iconrmsd follows directly in the list.
383 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
387 md->ib = get_ebin_space(md->ebin,
388 md->bTricl ? NTRICLBOXS : NBOXS,
389 md->bTricl ? tricl_boxs_nm : boxs_nm,
391 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
392 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
395 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
396 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
401 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
402 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
406 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
407 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
408 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
411 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
413 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
414 boxvel_nm, unit_vel);
418 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
420 if (ir->cos_accel != 0)
422 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
423 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
427 /* Energy monitoring */
428 for (i = 0; i < egNR; i++)
430 md->bEInd[i] = FALSE;
432 md->bEInd[egCOULSR] = TRUE;
433 md->bEInd[egLJSR ] = TRUE;
437 md->bEInd[egLJSR] = FALSE;
438 md->bEInd[egBHAMSR] = TRUE;
442 md->bEInd[egLJ14] = TRUE;
443 md->bEInd[egCOUL14] = TRUE;
446 for (i = 0; (i < egNR); i++)
453 n = groups->groups[SimulationAtomGroupType::EnergyOutput].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->groups[SimulationAtomGroupType::EnergyOutput].nr); i++)
468 ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
469 for (j = i; (j < groups->groups[SimulationAtomGroupType::EnergyOutput].nr); j++)
471 nj = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[j];
472 for (k = kk = 0; (k < egNR); k++)
476 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
477 *(groups->groupNames[ni]), *(groups->groupNames[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->groups[SimulationAtomGroupType::TemperatureCoupling].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);
532 // TODO redo the group name memory management to make it more clear
534 snew(grpnms, std::max(md->mde_n, md->mdeb_n)); // Just in case md->mdeb_n > md->mde_n
536 for (i = 0; (i < md->nTC); i++)
538 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
539 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
540 grpnms[i] = gmx_strdup(buf);
542 md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
544 for (i = 0; i < md->nTC; i++)
550 if (md->etc == etcNOSEHOOVER)
552 if (md->bPrintNHChains)
554 if (md->bNHC_trotter)
556 for (i = 0; (i < md->nTC); i++)
558 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
559 bufi = *(groups->groupNames[ni]);
560 for (j = 0; (j < md->nNHC); j++)
562 sprintf(buf, "Xi-%d-%s", j, bufi);
563 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
564 sprintf(buf, "vXi-%d-%s", j, bufi);
565 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
568 md->itc = get_ebin_space(md->ebin, md->mde_n,
569 grpnms, unit_invtime);
570 allocated = md->mde_n;
573 for (i = 0; (i < md->nTCP); i++)
575 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
576 for (j = 0; (j < md->nNHC); j++)
578 sprintf(buf, "Xi-%d-%s", j, bufi);
579 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
580 sprintf(buf, "vXi-%d-%s", j, bufi);
581 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
584 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
585 grpnms, unit_invtime);
586 allocated = md->mdeb_n;
591 for (i = 0; (i < md->nTC); i++)
593 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
594 bufi = *(groups->groupNames[ni]);
595 sprintf(buf, "Xi-%s", bufi);
596 grpnms[2*i] = gmx_strdup(buf);
597 sprintf(buf, "vXi-%s", bufi);
598 grpnms[2*i+1] = gmx_strdup(buf);
600 md->itc = get_ebin_space(md->ebin, md->mde_n,
601 grpnms, unit_invtime);
602 allocated = md->mde_n;
606 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
607 md->etc == etcVRESCALE)
609 for (i = 0; (i < md->nTC); i++)
611 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
612 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
613 grpnms[i] = gmx_strdup(buf);
615 md->itc = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
616 allocated = md->mde_n;
619 for (i = 0; i < allocated; i++)
625 md->nU = groups->groups[SimulationAtomGroupType::Acceleration].nr;
626 snew(md->tmp_v, md->nU);
629 snew(grpnms, 3*md->nU);
630 for (i = 0; (i < md->nU); i++)
632 ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
633 sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
634 grpnms[3*i+XX] = gmx_strdup(buf);
635 sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
636 grpnms[3*i+YY] = gmx_strdup(buf);
637 sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
638 grpnms[3*i+ZZ] = gmx_strdup(buf);
640 md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
641 for (i = 0; i < 3*md->nU; i++)
650 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
653 md->print_grpnms = nullptr;
655 /* check whether we're going to write dh histograms */
657 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
659 /* Currently dh histograms are only written with dynamics */
660 if (EI_DYNAMICS(ir->eI))
664 mde_delta_h_coll_init(md->dhc, ir);
666 md->fp_dhdl = nullptr;
667 snew(md->dE, ir->fepvals->n_lambda);
671 md->fp_dhdl = fp_dhdl;
672 snew(md->dE, ir->fepvals->n_lambda);
677 snew(md->temperatures, ir->fepvals->n_lambda);
678 for (i = 0; i < ir->fepvals->n_lambda; i++)
680 md->temperatures[i] = ir->simtempvals->temperatures[i];
686 //! Legacy cleanup function
687 void done_mdebin(t_mdebin *mdebin)
690 sfree(mdebin->tmp_r);
691 sfree(mdebin->tmp_v);
692 done_ebin(mdebin->ebin);
693 done_mde_delta_h_coll(mdebin->dhc);
695 sfree(mdebin->temperatures);
696 if (mdebin->nE > 1 && mdebin->print_grpnms != nullptr)
698 for (int n = 0; n < mdebin->nE; n++)
700 sfree(mdebin->print_grpnms[n]);
702 sfree(mdebin->print_grpnms);
710 /* print a lambda vector to a string
711 fep = the inputrec's FEP input data
712 i = the index of the lambda vector
713 get_native_lambda = whether to print the native lambda
714 get_names = whether to print the names rather than the values
715 str = the pre-allocated string buffer to print to. */
716 static void print_lambda_vector(t_lambda *fep, int i,
717 gmx_bool get_native_lambda, gmx_bool get_names,
723 for (j = 0; j < efptNR; j++)
725 if (fep->separate_dvdl[j])
730 str[0] = 0; /* reset the string */
733 str += sprintf(str, "("); /* set the opening parenthesis*/
735 for (j = 0; j < efptNR; j++)
737 if (fep->separate_dvdl[j])
741 if (get_native_lambda && fep->init_lambda >= 0)
743 str += sprintf(str, "%.4f", fep->init_lambda);
747 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
752 str += sprintf(str, "%s", efpt_singular_names[j]);
754 /* print comma for the next item */
757 str += sprintf(str, ", ");
764 /* and add the closing parenthesis */
769 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
770 const gmx_output_env_t *oenv)
773 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
774 *lambdastate = "\\lambda state";
775 int i, nsets, nsets_de, nsetsbegin;
776 int n_lambda_terms = 0;
777 t_lambda *fep = ir->fepvals; /* for simplicity */
778 t_expanded *expand = ir->expandedvals;
779 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
784 gmx_bool write_pV = FALSE;
786 /* count the number of different lambda terms */
787 for (i = 0; i < efptNR; i++)
789 if (fep->separate_dvdl[i])
795 std::string title, label_x, label_y;
796 if (fep->n_lambda == 0)
798 title = gmx::formatString("%s", dhdl);
799 label_x = gmx::formatString("Time (ps)");
800 label_y = gmx::formatString("%s (%s %s)",
801 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
805 title = gmx::formatString("%s and %s", dhdl, deltag);
806 label_x = gmx::formatString("Time (ps)");
807 label_y = gmx::formatString("%s and %s (%s %s)",
808 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
810 fp = gmx_fio_fopen(filename, "w+");
811 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
816 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
818 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
820 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
822 /* compatibility output */
823 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
827 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
829 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
831 buf += gmx::formatString("%s %d: %s = %s",
832 lambdastate, fep->init_fep_state,
833 lambda_name_str, lambda_vec_str);
836 xvgr_subtitle(fp, buf.c_str(), oenv);
840 if (fep->dhdl_derivatives == edhdlderivativesYES)
842 nsets_dhdl = n_lambda_terms;
844 /* count the number of delta_g states */
845 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
847 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
849 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
851 nsets += 1; /*add fep state for expanded ensemble */
854 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
856 nsets += 1; /* add energy to the dhdl as well */
860 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
862 nsetsextend += 1; /* for PV term, other terms possible if required for
863 the reduced potential (only needed with foreign
864 lambda, and only output when init_lambda is not
865 set in order to maintain compatibility of the
869 std::vector<std::string> setname(nsetsextend);
871 if (expand->elmcmove > elmcmoveNO)
873 /* state for the fep_vals, if we have alchemical sampling */
874 setname[s++] = "Thermodynamic state";
877 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
880 switch (fep->edHdLPrintEnergy)
882 case edHdLPrintEnergyPOTENTIAL:
883 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
885 case edHdLPrintEnergyTOTAL:
886 case edHdLPrintEnergyYES:
888 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
890 setname[s++] = energy;
893 if (fep->dhdl_derivatives == edhdlderivativesYES)
895 for (i = 0; i < efptNR; i++)
897 if (fep->separate_dvdl[i])
899 std::string derivative;
900 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
902 /* compatibility output */
903 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
907 double lam = fep->init_lambda;
908 if (fep->init_lambda < 0)
910 lam = fep->all_lambda[i][fep->init_fep_state];
912 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
915 setname[s++] = derivative;
920 if (fep->n_lambda > 0)
922 /* g_bar has to determine the lambda values used in this simulation
923 * from this xvg legend.
926 if (expand->elmcmove > elmcmoveNO)
928 nsetsbegin = 1; /* for including the expanded ensemble */
935 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
939 nsetsbegin += nsets_dhdl;
941 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
943 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
945 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
947 /* for compatible dhdl.xvg files */
948 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
952 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
957 /* print the temperature for this state if doing simulated annealing */
958 buf += gmx::formatString("T = %g (%s)",
959 ir->simtempvals->temperatures[s-(nsetsbegin)],
966 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
969 xvgrLegend(fp, setname, oenv);
980 //! Legacy update function
981 void upd_mdebin(t_mdebin *md,
986 gmx_enerdata_t *enerd,
995 gmx_ekindata_t *ekind,
997 const gmx::Constraints *constr)
999 int i, j, k, kk, n, gid;
1000 real crmsd[2], tmp6[6];
1001 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
1003 double store_dhdl[efptNR];
1004 real store_energy = 0;
1007 /* Do NOT use the box in the state variable, but the separate box provided
1008 * as an argument. This is because we sometimes need to write the box from
1009 * the last timestep to match the trajectory frames.
1011 add_ebin_indexed(md->ebin, md->ie, gmx::ArrayRef<gmx_bool>(md->bEner), enerd->term, bSum);
1014 crmsd[0] = constr->rmsd();
1015 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1022 bs[0] = box[XX][XX];
1023 bs[1] = box[YY][YY];
1024 bs[2] = box[ZZ][ZZ];
1025 bs[3] = box[YY][XX];
1026 bs[4] = box[ZZ][XX];
1027 bs[5] = box[ZZ][YY];
1032 bs[0] = box[XX][XX];
1033 bs[1] = box[YY][YY];
1034 bs[2] = box[ZZ][ZZ];
1037 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1038 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1039 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1040 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1041 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1045 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1046 not the instantaneous pressure */
1047 pv = vol*md->ref_p/PRESFAC;
1049 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1050 enthalpy = pv + enerd->term[F_ETOT];
1051 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1056 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1057 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1061 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1062 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1063 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1064 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1066 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1068 tmp6[0] = state->boxv[XX][XX];
1069 tmp6[1] = state->boxv[YY][YY];
1070 tmp6[2] = state->boxv[ZZ][ZZ];
1071 tmp6[3] = state->boxv[YY][XX];
1072 tmp6[4] = state->boxv[ZZ][XX];
1073 tmp6[5] = state->boxv[ZZ][YY];
1074 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1078 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1080 if (ekind && ekind->cosacc.cos_accel != 0)
1082 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1083 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1084 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1085 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1086 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1087 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1088 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1093 for (i = 0; (i < md->nEg); i++)
1095 for (j = i; (j < md->nEg); j++)
1097 gid = GID(i, j, md->nEg);
1098 for (k = kk = 0; (k < egNR); k++)
1102 eee[kk++] = enerd->grpp.ener[k][gid];
1105 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1113 for (i = 0; (i < md->nTC); i++)
1115 md->tmp_r[i] = ekind->tcstat[i].T;
1117 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1119 if (md->etc == etcNOSEHOOVER)
1121 /* whether to print Nose-Hoover chains: */
1122 if (md->bPrintNHChains)
1124 if (md->bNHC_trotter)
1126 for (i = 0; (i < md->nTC); i++)
1128 for (j = 0; j < md->nNHC; j++)
1131 md->tmp_r[2*k] = state->nosehoover_xi[k];
1132 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1135 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1139 for (i = 0; (i < md->nTCP); i++)
1141 for (j = 0; j < md->nNHC; j++)
1144 md->tmp_r[2*k] = state->nhpres_xi[k];
1145 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1148 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1153 for (i = 0; (i < md->nTC); i++)
1155 md->tmp_r[2*i] = state->nosehoover_xi[i];
1156 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1158 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1162 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1163 md->etc == etcVRESCALE)
1165 for (i = 0; (i < md->nTC); i++)
1167 md->tmp_r[i] = ekind->tcstat[i].lambda;
1169 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1173 if (ekind && md->nU > 1)
1175 for (i = 0; (i < md->nU); i++)
1177 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1179 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1182 ebin_increase_count(md->ebin, bSum);
1184 /* BAR + thermodynamic integration values */
1185 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1187 for (i = 0; i < enerd->n_lambda-1; i++)
1189 /* zero for simulated tempering */
1190 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1191 if (md->temperatures != nullptr)
1193 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1194 /* is this even useful to have at all? */
1195 md->dE[i] += (md->temperatures[i]/
1196 md->temperatures[state->fep_state]-1.0)*
1197 enerd->term[F_EKIN];
1203 fprintf(md->fp_dhdl, "%.4f", time);
1204 /* the current free energy state */
1206 /* print the current state if we are doing expanded ensemble */
1207 if (expand->elmcmove > elmcmoveNO)
1209 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1211 /* total energy (for if the temperature changes */
1213 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1215 switch (fep->edHdLPrintEnergy)
1217 case edHdLPrintEnergyPOTENTIAL:
1218 store_energy = enerd->term[F_EPOT];
1220 case edHdLPrintEnergyTOTAL:
1221 case edHdLPrintEnergyYES:
1223 store_energy = enerd->term[F_ETOT];
1225 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1228 if (fep->dhdl_derivatives == edhdlderivativesYES)
1230 for (i = 0; i < efptNR; i++)
1232 if (fep->separate_dvdl[i])
1234 /* assumes F_DVDL is first */
1235 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1239 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1241 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1245 (md->epc != epcNO) &&
1246 (enerd->n_lambda > 0) &&
1247 (fep->init_lambda < 0))
1249 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1250 there are alternate state
1251 lambda and we're not in
1252 compatibility mode */
1254 fprintf(md->fp_dhdl, "\n");
1255 /* and the binary free energy output */
1257 if (md->dhc && bDoDHDL)
1260 for (i = 0; i < efptNR; i++)
1262 if (fep->separate_dvdl[i])
1264 /* assumes F_DVDL is first */
1265 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1269 store_energy = enerd->term[F_ETOT];
1270 /* store_dh is dE */
1271 mde_delta_h_coll_add_dh(md->dhc,
1272 static_cast<double>(state->fep_state),
1276 md->dE + fep->lambda_start_n,
1284 void EnergyOutput::recordNonEnergyStep()
1286 ebin_increase_count(mdebin->ebin, false);
1292 //! Legacy output function
1293 void npr(FILE *log, int n, char c)
1295 for (; (n > 0); n--)
1297 fprintf(log, "%c", c);
1301 //! Legacy output function
1302 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));
1325 void print_ebin_header(FILE *log, int64_t steps, double time)
1329 fprintf(log, " %12s %12s\n"
1331 "Step", "Time", gmx_step_str(steps, buf), time);
1337 // TODO It is too many responsibilities for this function to handle
1338 // both .edr and .log output for both per-time and time-average data.
1339 //! Legacy ebin output function
1340 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1342 int64_t step, double time,
1344 t_mdebin *md, t_fcdata *fcd,
1345 SimulationGroups *groups, t_grpopts *opts,
1348 /*static char **grpnms=NULL;*/
1350 int i, j, n, ni, nj, b;
1353 /* these are for the old-style blocks (1 subblock, only reals), because
1354 there can be only one per ID for these */
1361 if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1365 fprintf(log, "Not enough data recorded to report energy averages\n");
1376 fr.nsteps = md->ebin->nsteps;
1377 fr.dt = md->delta_t;
1378 fr.nsum = md->ebin->nsum;
1379 fr.nre = (bEne) ? md->ebin->nener : 0;
1380 fr.ener = md->ebin->e;
1381 ndisre = bDR ? fcd->disres.npair : 0;
1382 /* Optional additional old-style (real-only) blocks. */
1383 for (i = 0; i < enxNR; i++)
1387 if (bOR && fcd->orires.nr > 0)
1389 diagonalize_orires_tensors(&(fcd->orires));
1390 nr[enxOR] = fcd->orires.nr;
1391 block[enxOR] = fcd->orires.otav;
1393 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1395 block[enxORI] = fcd->orires.oinsl;
1396 id[enxORI] = enxORI;
1397 nr[enxORT] = fcd->orires.nex*12;
1398 block[enxORT] = fcd->orires.eig;
1399 id[enxORT] = enxORT;
1402 /* whether we are going to wrte anything out: */
1403 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1406 /* the old-style blocks go first */
1408 for (i = 0; i < enxNR; i++)
1415 add_blocks_enxframe(&fr, fr.nblock);
1416 for (b = 0; b < fr.nblock; b++)
1418 add_subblocks_enxblock(&(fr.block[b]), 1);
1419 fr.block[b].id = id[b];
1420 fr.block[b].sub[0].nr = nr[b];
1422 fr.block[b].sub[0].type = xdr_datatype_float;
1423 fr.block[b].sub[0].fval = block[b];
1425 fr.block[b].sub[0].type = xdr_datatype_double;
1426 fr.block[b].sub[0].dval = block[b];
1430 /* check for disre block & fill it. */
1435 add_blocks_enxframe(&fr, fr.nblock);
1437 add_subblocks_enxblock(&(fr.block[db]), 2);
1438 fr.block[db].id = enxDISRE;
1439 fr.block[db].sub[0].nr = ndisre;
1440 fr.block[db].sub[1].nr = ndisre;
1442 fr.block[db].sub[0].type = xdr_datatype_float;
1443 fr.block[db].sub[1].type = xdr_datatype_float;
1444 fr.block[db].sub[0].fval = fcd->disres.rt;
1445 fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1447 fr.block[db].sub[0].type = xdr_datatype_double;
1448 fr.block[db].sub[1].type = xdr_datatype_double;
1449 fr.block[db].sub[0].dval = fcd->disres.rt;
1450 fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1453 /* here we can put new-style blocks */
1455 /* Free energy perturbation blocks */
1458 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1461 /* we can now free & reset the data in the blocks */
1464 mde_delta_h_coll_reset(md->dhc);
1467 /* AWH bias blocks. */
1468 if (awh != nullptr) // TODO: add boolean in t_mdebin. See in mdebin.h.
1470 awh->writeToEnergyFrame(step, &fr);
1473 /* do the actual I/O */
1474 do_enx(fp_ene, &fr);
1477 /* We have stored the sums, so reset the sum history */
1478 reset_ebin_sums(md->ebin);
1486 pprint(log, "A V E R A G E S", md);
1492 pprint(log, "R M S - F L U C T U A T I O N S", md);
1496 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1503 for (i = 0; i < opts->ngtc; i++)
1505 if (opts->annealing[i] != eannNO)
1507 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1508 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i]]),
1513 if (mode == eprNORMAL && bOR && fcd->orires.nr > 0)
1515 print_orires_log(log, &(fcd->orires));
1517 fprintf(log, " Energies (%s)\n", unit_energy);
1518 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1521 if (mode == eprAVER)
1525 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1531 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1532 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1534 fprintf(log, " Force Virial (%s)\n", unit_energy);
1535 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1540 fprintf(log, " Total Virial (%s)\n", unit_energy);
1541 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1543 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1544 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1549 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1550 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1556 if (md->print_grpnms == nullptr)
1558 snew(md->print_grpnms, md->nE);
1560 for (i = 0; (i < md->nEg); i++)
1562 ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
1563 for (j = i; (j < md->nEg); j++)
1565 nj = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[j];
1566 sprintf(buf, "%s-%s", *(groups->groupNames[ni]),
1567 *(groups->groupNames[nj]));
1568 md->print_grpnms[n++] = gmx_strdup(buf);
1572 sprintf(buf, "Epot (%s)", unit_energy);
1573 fprintf(log, "%15s ", buf);
1574 for (i = 0; (i < egNR); i++)
1578 fprintf(log, "%12s ", egrp_nm[i]);
1582 for (i = 0; (i < md->nE); i++)
1584 fprintf(log, "%15s", md->print_grpnms[i]);
1585 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1592 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1597 fprintf(log, "%15s %12s %12s %12s\n",
1598 "Group", "Ux", "Uy", "Uz");
1599 for (i = 0; (i < md->nU); i++)
1601 ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
1602 fprintf(log, "%15s", *groups->groupNames[ni]);
1603 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1612 //! Legacy update function
1613 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1615 const t_ebin * const ebin = mdebin->ebin;
1617 enerhist->nsteps = ebin->nsteps;
1618 enerhist->nsum = ebin->nsum;
1619 enerhist->nsteps_sim = ebin->nsteps_sim;
1620 enerhist->nsum_sim = ebin->nsum_sim;
1624 /* This will only actually resize the first time */
1625 enerhist->ener_ave.resize(ebin->nener);
1626 enerhist->ener_sum.resize(ebin->nener);
1628 for (int i = 0; i < ebin->nener; i++)
1630 enerhist->ener_ave[i] = ebin->e[i].eav;
1631 enerhist->ener_sum[i] = ebin->e[i].esum;
1635 if (ebin->nsum_sim > 0)
1637 /* This will only actually resize the first time */
1638 enerhist->ener_sum_sim.resize(ebin->nener);
1640 for (int i = 0; i < ebin->nener; i++)
1642 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1647 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1651 //! Legacy restore function
1652 void restore_energyhistory_from_state(t_mdebin * mdebin,
1653 const energyhistory_t * enerhist)
1655 unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1657 GMX_RELEASE_ASSERT(enerhist, "Need valid history to restore");
1659 if ((enerhist->nsum > 0 && nener != enerhist->ener_sum.size()) ||
1660 (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1662 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1663 nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1666 mdebin->ebin->nsteps = enerhist->nsteps;
1667 mdebin->ebin->nsum = enerhist->nsum;
1668 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1669 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1671 for (int i = 0; i < mdebin->ebin->nener; i++)
1673 mdebin->ebin->e[i].eav =
1674 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1675 mdebin->ebin->e[i].esum =
1676 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1677 mdebin->ebin->e_sim[i].esum =
1678 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1682 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());
1688 EnergyOutput::EnergyOutput()
1693 void EnergyOutput::prepare(ener_file *fp_ene,
1694 const gmx_mtop_t *mtop,
1695 const t_inputrec *ir,
1696 const pull_t *pull_work,
1700 mdebin = init_mdebin(fp_ene, mtop, ir, pull_work, fp_dhdl, isRerun);
1703 EnergyOutput::~EnergyOutput()
1705 done_mdebin(mdebin);
1708 t_ebin *EnergyOutput::getEbin()
1710 return mdebin->ebin;
1713 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
1717 gmx_enerdata_t *enerd,
1726 gmx_ekindata_t *ekind,
1728 const gmx::Constraints *constr)
1730 upd_mdebin(mdebin, bDoDHDL, bSum, time, tmass, enerd, state, fep,
1731 expand, box, svir, fvir, vir, pres, ekind, mu_tot, constr);
1734 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1736 int64_t step, double time,
1739 SimulationGroups *groups, t_grpopts *opts,
1742 print_ebin(fp_ene, bEne, bDR, bOR, log, step, time, mode,
1743 mdebin, fcd, groups, opts, awh);
1746 int EnergyOutput::numEnergyTerms() const
1748 return mdebin->ebin->nener;
1751 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1753 update_energyhistory(enerhist, mdebin);
1756 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1758 restore_energyhistory_from_state(mdebin, &enerhist);