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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Defines code that writes energy-like quantities.
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \author Paul Bauer <paul.bauer.q@gmail.com>
43 * \author Artem Zhmurov <zhmurov@gmail.com>
45 * \ingroup module_mdlib
49 #include "energyoutput.h"
58 #include "gromacs/applied_forces/awh/awh.h"
59 #include "gromacs/fileio/enxio.h"
60 #include "gromacs/fileio/gmxfio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/gmxlib/network.h"
63 #include "gromacs/listed_forces/disre.h"
64 #include "gromacs/listed_forces/orires.h"
65 #include "gromacs/math/functions.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/constr.h"
69 #include "gromacs/mdlib/ebin.h"
70 #include "gromacs/mdlib/mdebin_bar.h"
71 #include "gromacs/mdrunutility/handlerestart.h"
72 #include "gromacs/mdtypes/energyhistory.h"
73 #include "gromacs/mdtypes/fcdata.h"
74 #include "gromacs/mdtypes/group.h"
75 #include "gromacs/mdtypes/inputrec.h"
76 #include "gromacs/mdtypes/md_enums.h"
77 #include "gromacs/mdtypes/state.h"
78 #include "gromacs/pbcutil/pbc.h"
79 #include "gromacs/pulling/pull.h"
80 #include "gromacs/topology/mtop_util.h"
81 #include "gromacs/trajectory/energyframe.h"
82 #include "gromacs/utility/arraysize.h"
83 #include "gromacs/utility/enumerationhelpers.h"
84 #include "gromacs/utility/fatalerror.h"
85 #include "gromacs/utility/mdmodulesnotifiers.h"
86 #include "gromacs/utility/smalloc.h"
87 #include "gromacs/utility/stringutil.h"
89 #include "energydrifttracker.h"
91 //! Labels for energy file quantities
93 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
94 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
96 static constexpr std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
98 static constexpr std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
99 "Box-YX", "Box-ZX", "Box-ZY" };
101 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
102 static const char* vol_nm[] = { "Volume" };
104 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
105 static const char* dens_nm[] = { "Density" };
107 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
108 static const char* pv_nm[] = { "pV" };
110 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
111 static const char* enthalpy_nm[] = { "Enthalpy" };
113 static constexpr std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
114 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
116 const char* enumValueToString(NonBondedEnergyTerms enumValue)
118 static constexpr gmx::EnumerationArray<NonBondedEnergyTerms, const char*> nonBondedEnergyTermTypeNames = {
119 "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14"
121 return nonBondedEnergyTermTypeNames[enumValue];
129 /*! \brief Energy output class
131 * This is the collection of energy averages collected during mdrun, and to
132 * be written out to the .edr file.
134 * \todo Use more std containers.
135 * \todo Write free-energy output also to energy file (after adding more tests)
137 EnergyOutput::EnergyOutput(ener_file* fp_ene,
138 const gmx_mtop_t& mtop,
139 const t_inputrec& inputrec,
140 const pull_t* pull_work,
143 const StartingBehavior startingBehavior,
144 const bool simulationsShareState,
145 const MDModulesNotifiers& mdModulesNotifiers)
147 const char* ener_nm[F_NRE];
148 static const char* vir_nm[] = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
149 "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
150 static const char* pres_nm[] = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
151 "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
152 static const char* surft_nm[] = { "#Surf*SurfTen" };
153 static const char* mu_nm[] = { "Mu-X", "Mu-Y", "Mu-Z" };
154 static const char* vcos_nm[] = { "2CosZ*Vel-X" };
155 static const char* visc_nm[] = { "1/Viscosity" };
156 static const char* baro_nm[] = { "Barostat" };
158 const SimulationGroups* groups;
162 int i, j, ni, nj, n, ncon, nset;
165 if (EI_DYNAMICS(inputrec.eI))
167 delta_t_ = inputrec.delta_t;
174 groups = &mtop.groups;
176 bBHAM = (mtop.ffparams.numTypes() > 0) && (mtop.ffparams.functype[0] == F_BHAM);
177 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
179 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
180 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
181 bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
185 if (ncon > 0 && inputrec.eConstrAlg == ConstraintAlgorithm::Lincs)
195 /* Energy monitoring */
196 for (auto& term : bEInd_)
201 // Setting true only to those energy terms, that have active interactions and
202 // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
203 for (i = 0; i < F_NRE; i++)
205 bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
206 && ((interaction_function[i].flags & IF_VSITE) == 0);
211 bEner_[F_EKIN] = EI_DYNAMICS(inputrec.eI);
212 bEner_[F_ETOT] = EI_DYNAMICS(inputrec.eI);
213 bEner_[F_TEMP] = EI_DYNAMICS(inputrec.eI);
215 bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(&inputrec);
216 bEner_[F_PDISPCORR] = (inputrec.eDispCorr != DispersionCorrectionType::No);
217 bEner_[F_PRES] = true;
220 bEner_[F_LJ] = !bBHAM;
221 bEner_[F_BHAM] = bBHAM;
222 bEner_[F_EQM] = inputrec.bQMMM;
223 bEner_[F_RF_EXCL] = (EEL_RF(inputrec.coulombtype) && inputrec.cutoff_scheme == CutoffScheme::Group);
224 bEner_[F_COUL_RECIP] = EEL_FULL(inputrec.coulombtype);
225 bEner_[F_LJ_RECIP] = EVDW_PME(inputrec.vdwtype);
226 bEner_[F_LJ14] = b14;
227 bEner_[F_COUL14] = b14;
228 bEner_[F_LJC14_Q] = false;
229 bEner_[F_LJC_PAIRS_NB] = false;
232 bEner_[F_DVDL_COUL] = (inputrec.efep != FreeEnergyPerturbationType::No)
233 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Coul];
234 bEner_[F_DVDL_VDW] = (inputrec.efep != FreeEnergyPerturbationType::No)
235 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Vdw];
236 bEner_[F_DVDL_BONDED] = (inputrec.efep != FreeEnergyPerturbationType::No)
237 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Bonded];
238 bEner_[F_DVDL_RESTRAINT] =
239 (inputrec.efep != FreeEnergyPerturbationType::No)
240 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Restraint];
241 bEner_[F_DKDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
242 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Mass];
243 bEner_[F_DVDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
244 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep];
246 bEner_[F_CONSTR] = false;
247 bEner_[F_CONSTRNC] = false;
248 bEner_[F_SETTLE] = false;
250 bEner_[F_COUL_SR] = true;
251 bEner_[F_EPOT] = true;
253 bEner_[F_DISPCORR] = (inputrec.eDispCorr != DispersionCorrectionType::No);
254 bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
255 bEner_[F_ORIRESDEV] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
256 bEner_[F_COM_PULL] = ((inputrec.bPull && pull_have_potential(*pull_work)) || inputrec.bRot);
258 MDModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
259 mdModulesNotifiers.simulationSetupNotifier_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
261 bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
264 // Counting the energy terms that will be printed and saving their names
266 for (i = 0; i < F_NRE; i++)
270 ener_nm[f_nre_] = interaction_function[i].longname;
275 epc_ = isRerun ? PressureCoupling::No : inputrec.epc;
276 bDiagPres_ = !TRICLINIC(inputrec.ref_p) && !isRerun;
277 ref_p_ = (inputrec.ref_p[XX][XX] + inputrec.ref_p[YY][YY] + inputrec.ref_p[ZZ][ZZ]) / DIM;
278 bTricl_ = TRICLINIC(inputrec.compress) || TRICLINIC(inputrec.deform);
279 bDynBox_ = inputrecDynamicBox(&inputrec);
280 etc_ = isRerun ? TemperatureCoupling::No : inputrec.etc;
281 bNHC_trotter_ = inputrecNvtTrotter(&inputrec) && !isRerun;
282 bPrintNHChains_ = inputrec.bPrintNHChains && !isRerun;
283 bMTTK_ = (inputrecNptTrotter(&inputrec) || inputrecNphTrotter(&inputrec)) && !isRerun;
284 bMu_ = inputrecNeedMutot(&inputrec);
288 /* Pass NULL for unit to let get_ebin_space determine the units
289 * for interaction_function[i].longname
291 ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
294 /* This should be called directly after the call for ie_,
295 * such that iconrmsd_ follows directly in the list.
297 iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
301 ib_ = get_ebin_space(ebin_,
302 bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
303 bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
305 ivol_ = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
306 idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
309 ipv_ = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
310 ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
315 ivir_ = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
316 ipres_ = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
317 isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
319 if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
321 ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
325 imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
327 if (inputrec.cos_accel != 0)
329 ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
330 ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
333 /* Energy monitoring */
334 for (auto& term : bEInd_)
338 bEInd_[NonBondedEnergyTerms::CoulombSR] = true;
339 bEInd_[NonBondedEnergyTerms::LJSR] = true;
343 bEInd_[NonBondedEnergyTerms::LJSR] = false;
344 bEInd_[NonBondedEnergyTerms::BuckinghamSR] = true;
348 bEInd_[NonBondedEnergyTerms::LJ14] = true;
349 bEInd_[NonBondedEnergyTerms::Coulomb14] = true;
352 for (auto term : bEInd_)
359 n = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
361 nE_ = (n * (n + 1)) / 2;
368 for (int k = 0; (k < nEc_); k++)
370 snew(gnm[k], STRLEN);
372 for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
374 ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
375 for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
377 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
379 for (auto key : keysOf(bEInd_))
385 enumValueToString(key),
386 *(groups->groupNames[ni]),
387 *(groups->groupNames[nj]));
391 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
395 for (int k = 0; (k < nEc_); k++)
403 gmx_incons("Number of energy terms wrong");
407 nTC_ = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
408 nNHC_ = inputrec.opts.nhchainlength; /* shorthand for number of NH chains */
411 nTCP_ = 1; /* assume only one possible coupling system for barostat
418 if (etc_ == TemperatureCoupling::NoseHoover)
422 mde_n_ = 2 * nNHC_ * nTC_;
428 if (epc_ == PressureCoupling::Mttk)
430 mdeb_n_ = 2 * nNHC_ * nTCP_;
439 tmp_r_.resize(mde_n_);
440 // TODO redo the group name memory management to make it more clear
442 snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
444 for (i = 0; (i < nTC_); i++)
446 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
447 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
448 grpnms[i] = gmx_strdup(buf);
450 itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
451 for (i = 0; i < nTC_; i++)
457 if (etc_ == TemperatureCoupling::NoseHoover)
463 for (i = 0; (i < nTC_); i++)
465 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
466 bufi = *(groups->groupNames[ni]);
467 for (j = 0; (j < nNHC_); j++)
469 sprintf(buf, "Xi-%d-%s", j, bufi);
470 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
471 sprintf(buf, "vXi-%d-%s", j, bufi);
472 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
475 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
479 for (i = 0; (i < nTCP_); i++)
481 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
482 for (j = 0; (j < nNHC_); j++)
484 sprintf(buf, "Xi-%d-%s", j, bufi);
485 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
486 sprintf(buf, "vXi-%d-%s", j, bufi);
487 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
490 itcb_ = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
496 for (i = 0; (i < nTC_); i++)
498 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
499 bufi = *(groups->groupNames[ni]);
500 sprintf(buf, "Xi-%s", bufi);
501 grpnms[2 * i] = gmx_strdup(buf);
502 sprintf(buf, "vXi-%s", bufi);
503 grpnms[2 * i + 1] = gmx_strdup(buf);
505 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
510 else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
511 || etc_ == TemperatureCoupling::VRescale)
513 for (i = 0; (i < nTC_); i++)
515 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
516 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
517 grpnms[i] = gmx_strdup(buf);
519 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, "");
523 for (i = 0; i < allocated; i++)
529 /* Note that fp_ene should be valid on the master rank and null otherwise */
530 if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
532 do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
535 /* check whether we're going to write dh histograms */
537 if (inputrec.fepvals->separate_dhdl_file == SeparateDhdlFile::No)
539 /* Currently dh histograms are only written with dynamics */
540 if (EI_DYNAMICS(inputrec.eI))
542 dhc_ = std::make_unique<t_mde_delta_h_coll>(inputrec);
545 dE_.resize(inputrec.fepvals->n_lambda);
550 dE_.resize(inputrec.fepvals->n_lambda);
552 if (inputrec.bSimTemp)
554 temperatures_ = inputrec.simtempvals->temperatures;
557 if (EI_MD(inputrec.eI) && !simulationsShareState)
559 conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop.natoms);
563 EnergyOutput::~EnergyOutput()
570 /*! \brief Print a lambda vector to a string
572 * \param[in] fep The inputrec's FEP input data
573 * \param[in] i The index of the lambda vector
574 * \param[in] get_native_lambda Whether to print the native lambda
575 * \param[in] get_names Whether to print the names rather than the values
576 * \param[in,out] str The pre-allocated string buffer to print to.
578 static void print_lambda_vector(t_lambda* fep, int i, bool get_native_lambda, bool get_names, char* str)
583 for (auto j : keysOf(fep->separate_dvdl))
585 if (fep->separate_dvdl[j])
590 str[0] = 0; /* reset the string */
593 str += sprintf(str, "("); /* set the opening parenthesis*/
595 for (auto j : keysOf(fep->separate_dvdl))
597 if (fep->separate_dvdl[j])
601 if (get_native_lambda && fep->init_lambda >= 0)
603 str += sprintf(str, "%.4f", fep->init_lambda);
607 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
612 str += sprintf(str, "%s", enumValueToStringSingular(j));
614 /* print comma for the next item */
617 str += sprintf(str, ", ");
624 /* and add the closing parenthesis */
629 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv)
632 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
633 *lambdastate = "\\lambda state";
634 int i, nsets, nsets_de, nsetsbegin;
635 int n_lambda_terms = 0;
636 t_lambda* fep = ir->fepvals.get(); /* for simplicity */
637 t_expanded* expand = ir->expandedvals.get();
638 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
643 bool write_pV = false;
645 /* count the number of different lambda terms */
646 for (auto i : keysOf(fep->separate_dvdl))
648 if (fep->separate_dvdl[i])
654 std::string title, label_x, label_y;
655 if (fep->n_lambda == 0)
657 title = gmx::formatString("%s", dhdl);
658 label_x = gmx::formatString("Time (ps)");
659 label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
663 title = gmx::formatString("%s and %s", dhdl, deltag);
664 label_x = gmx::formatString("Time (ps)");
665 label_y = gmx::formatString(
666 "%s and %s (%s %s)", dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
668 fp = gmx_fio_fopen(filename, "w+");
669 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
674 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
676 if ((ir->efep != FreeEnergyPerturbationType::SlowGrowth)
677 && (ir->efep != FreeEnergyPerturbationType::Expanded))
679 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
681 /* compatibility output */
682 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
686 print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
687 print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
688 buf += gmx::formatString(
689 "%s %d: %s = %s", lambdastate, fep->init_fep_state, lambda_name_str, lambda_vec_str);
692 xvgr_subtitle(fp, buf.c_str(), oenv);
696 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
698 nsets_dhdl = n_lambda_terms;
700 /* count the number of delta_g states */
701 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
703 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
705 if (fep->n_lambda > 0 && (expand->elmcmove > LambdaMoveCalculation::No))
707 nsets += 1; /*add fep state for expanded ensemble */
710 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
712 nsets += 1; /* add energy to the dhdl as well */
716 if ((ir->epc != PressureCoupling::No) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
718 nsetsextend += 1; /* for PV term, other terms possible if required for
719 the reduced potential (only needed with foreign
720 lambda, and only output when init_lambda is not
721 set in order to maintain compatibility of the
725 std::vector<std::string> setname(nsetsextend);
727 if (expand->elmcmove > LambdaMoveCalculation::No)
729 /* state for the fep_vals, if we have alchemical sampling */
730 setname[s++] = "Thermodynamic state";
733 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
736 switch (fep->edHdLPrintEnergy)
738 case FreeEnergyPrintEnergy::Potential:
739 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
741 case FreeEnergyPrintEnergy::Total:
742 case FreeEnergyPrintEnergy::Yes:
743 default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
745 setname[s++] = energy;
748 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
750 for (auto i : keysOf(fep->separate_dvdl))
752 if (fep->separate_dvdl[i])
754 std::string derivative;
755 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
757 /* compatibility output */
758 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
762 double lam = fep->init_lambda;
763 if (fep->init_lambda < 0)
765 lam = fep->all_lambda[i][fep->init_fep_state];
767 derivative = gmx::formatString("%s %s = %.4f", dhdl, enumValueToStringSingular(i), lam);
769 setname[s++] = derivative;
774 if (fep->n_lambda > 0)
776 /* g_bar has to determine the lambda values used in this simulation
777 * from this xvg legend.
780 if (expand->elmcmove > LambdaMoveCalculation::No)
782 nsetsbegin = 1; /* for including the expanded ensemble */
789 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
793 nsetsbegin += nsets_dhdl;
795 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
797 print_lambda_vector(fep, i, false, false, lambda_vec_str);
799 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
801 /* for compatible dhdl.xvg files */
802 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
806 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
811 /* print the temperature for this state if doing simulated annealing */
812 buf += gmx::formatString(
813 "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
819 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
822 xvgrLegend(fp, setname, oenv);
831 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
835 const gmx_enerdata_t* enerd,
837 const t_expanded* expand,
839 PTCouplingArrays ptCouplingArrays,
843 const gmx_ekindata_t* ekind,
845 const gmx::Constraints* constr)
847 int j, k, kk, n, gid;
848 real crmsd[2], tmp6[6];
849 real bs[tricl_boxs_nm.size()], vol, dens, enthalpy;
850 real eee[static_cast<int>(NonBondedEnergyTerms::Count)];
851 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, double> store_dhdl;
852 real store_energy = 0;
854 real pv = 0.0; // static analyzer warns about uninitialized variable warnings here.
856 /* Do NOT use the box in the state variable, but the separate box provided
857 * as an argument. This is because we sometimes need to write the box from
858 * the last timestep to match the trajectory frames.
860 add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
863 crmsd[0] = constr->rmsd();
864 add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
877 nboxs = tricl_boxs_nm.size();
884 nboxs = boxs_nm.size();
886 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
887 dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
888 add_ebin(ebin_, ib_, nboxs, bs, bSum);
889 add_ebin(ebin_, ivol_, 1, &vol, bSum);
890 add_ebin(ebin_, idens_, 1, &dens, bSum);
894 /* This is pV (in kJ/mol). The pressure is the reference pressure,
895 not the instantaneous pressure */
896 pv = vol * ref_p_ / gmx::c_presfac;
898 add_ebin(ebin_, ipv_, 1, &pv, bSum);
899 enthalpy = pv + enerd->term[F_ETOT];
900 add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
905 add_ebin(ebin_, ivir_, 9, vir[0], bSum);
906 add_ebin(ebin_, ipres_, 9, pres[0], bSum);
907 tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
908 add_ebin(ebin_, isurft_, 1, &tmp, bSum);
910 if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
912 tmp6[0] = ptCouplingArrays.boxv[XX][XX];
913 tmp6[1] = ptCouplingArrays.boxv[YY][YY];
914 tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
915 tmp6[3] = ptCouplingArrays.boxv[YY][XX];
916 tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
917 tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
918 add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
922 add_ebin(ebin_, imu_, 3, mu_tot, bSum);
924 if (ekind && ekind->cosacc.cos_accel != 0)
926 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
927 dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
928 add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
929 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
931 / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * gmx::c_pico) * dens
932 * gmx::square(box[ZZ][ZZ] * gmx::c_nano / (2 * M_PI)));
933 add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
938 for (int i = 0; (i < nEg_); i++)
940 for (j = i; (j < nEg_); j++)
942 gid = GID(i, j, nEg_);
943 for (k = kk = 0; (k < static_cast<int>(NonBondedEnergyTerms::Count)); k++)
947 eee[kk++] = enerd->grpp.energyGroupPairTerms[k][gid];
950 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
958 for (int i = 0; (i < nTC_); i++)
960 tmp_r_[i] = ekind->tcstat[i].T;
962 add_ebin(ebin_, itemp_, nTC_, tmp_r_.data(), bSum);
964 if (etc_ == TemperatureCoupling::NoseHoover)
966 /* whether to print Nose-Hoover chains: */
971 for (int i = 0; (i < nTC_); i++)
973 for (j = 0; j < nNHC_; j++)
976 tmp_r_[2 * k] = ptCouplingArrays.nosehoover_xi[k];
977 tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
980 add_ebin(ebin_, itc_, mde_n_, tmp_r_.data(), bSum);
984 for (int i = 0; (i < nTCP_); i++)
986 for (j = 0; j < nNHC_; j++)
989 tmp_r_[2 * k] = ptCouplingArrays.nhpres_xi[k];
990 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
993 add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_.data(), bSum);
998 for (int i = 0; (i < nTC_); i++)
1000 tmp_r_[2 * i] = ptCouplingArrays.nosehoover_xi[i];
1001 tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1003 add_ebin(ebin_, itc_, mde_n_, tmp_r_.data(), bSum);
1007 else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
1008 || etc_ == TemperatureCoupling::VRescale)
1010 for (int i = 0; (i < nTC_); i++)
1012 tmp_r_[i] = ekind->tcstat[i].lambda;
1014 add_ebin(ebin_, itc_, nTC_, tmp_r_.data(), bSum);
1018 ebin_increase_count(1, ebin_, bSum);
1020 // BAR + thermodynamic integration values
1021 if ((fp_dhdl_ || dhc_) && bDoDHDL)
1023 const auto& foreignTerms = enerd->foreignLambdaTerms;
1024 for (int i = 0; i < foreignTerms.numLambdas(); i++)
1026 /* zero for simulated tempering */
1027 dE_[i] = foreignTerms.deltaH(i);
1028 if (!temperatures_.empty())
1030 GMX_RELEASE_ASSERT(gmx::ssize(temperatures_) > fep_state,
1031 "Number of lambdas in state is bigger then in input record");
1033 gmx::ssize(temperatures_) >= foreignTerms.numLambdas(),
1034 "Number of lambdas in energy data is bigger then in input record");
1035 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1036 /* is this even useful to have at all? */
1037 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1043 fprintf(fp_dhdl_, "%.4f", time);
1044 /* the current free energy state */
1046 /* print the current state if we are doing expanded ensemble */
1047 if (expand->elmcmove > LambdaMoveCalculation::No)
1049 fprintf(fp_dhdl_, " %4d", fep_state);
1051 /* total energy (for if the temperature changes */
1053 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
1055 switch (fep->edHdLPrintEnergy)
1057 case FreeEnergyPrintEnergy::Potential:
1058 store_energy = enerd->term[F_EPOT];
1060 case FreeEnergyPrintEnergy::Total:
1061 case FreeEnergyPrintEnergy::Yes:
1062 default: store_energy = enerd->term[F_ETOT];
1064 fprintf(fp_dhdl_, " %#.8g", store_energy);
1067 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
1069 for (auto i : keysOf(fep->separate_dvdl))
1071 if (fep->separate_dvdl[i])
1073 /* assumes F_DVDL is first */
1074 fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + static_cast<int>(i)]);
1078 for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1080 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1082 if (bDynBox_ && bDiagPres_ && (epc_ != PressureCoupling::No)
1083 && foreignTerms.numLambdas() > 0 && (fep->init_lambda < 0))
1085 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1086 there are alternate state
1087 lambda and we're not in
1088 compatibility mode */
1090 fprintf(fp_dhdl_, "\n");
1091 /* and the binary free energy output */
1093 if (dhc_ && bDoDHDL)
1096 for (auto i : keysOf(fep->separate_dvdl))
1098 if (fep->separate_dvdl[i])
1100 /* assumes F_DVDL is first */
1101 store_dhdl[idhdl] = enerd->term[F_DVDL + static_cast<int>(i)];
1105 store_energy = enerd->term[F_ETOT];
1106 /* store_dh is dE */
1107 mde_delta_h_coll_add_dh(dhc_.get(),
1108 static_cast<double>(fep_state),
1112 dE_.data() + fep->lambda_start_n,
1117 if (conservedEnergyTracker_)
1119 conservedEnergyTracker_->addPoint(
1120 time, bEner_[F_ECONSERVED] ? enerd->term[F_ECONSERVED] : enerd->term[F_ETOT]);
1124 void EnergyOutput::recordNonEnergyStep()
1126 ebin_increase_count(1, ebin_, false);
1129 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1138 gmx_step_str(steps, buf),
1142 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1156 fr.nsteps = ebin_->nsteps;
1158 fr.nsum = ebin_->nsum;
1159 fr.nre = (bEne) ? ebin_->nener : 0;
1161 int ndisre = bDR ? fcd->disres->npair : 0;
1162 /* these are for the old-style blocks (1 subblock, only reals), because
1163 there can be only one per ID for these */
1167 /* Optional additional old-style (real-only) blocks. */
1168 for (int i = 0; i < enxNR; i++)
1173 if (bOR && fcd->orires)
1175 t_oriresdata& orires = *fcd->orires;
1176 diagonalize_orires_tensors(&orires);
1177 nr[enxOR] = orires.numRestraints;
1178 block[enxOR] = orires.orientationsTimeAndEnsembleAv.data();
1180 nr[enxORI] = (orires.orientations.data() != orires.orientationsTimeAndEnsembleAv.data())
1181 ? orires.numRestraints
1183 block[enxORI] = orires.orientations.data();
1184 id[enxORI] = enxORI;
1185 nr[enxORT] = ssize(orires.eigenOutput);
1186 block[enxORT] = orires.eigenOutput.data();
1187 id[enxORT] = enxORT;
1190 /* whether we are going to write anything out: */
1191 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1193 /* the old-style blocks go first */
1195 for (int i = 0; i < enxNR; i++)
1202 add_blocks_enxframe(&fr, fr.nblock);
1203 for (int b = 0; b < fr.nblock; b++)
1205 add_subblocks_enxblock(&(fr.block[b]), 1);
1206 fr.block[b].id = id[b];
1207 fr.block[b].sub[0].nr = nr[b];
1209 fr.block[b].sub[0].type = XdrDataType::Float;
1210 fr.block[b].sub[0].fval = block[b];
1212 fr.block[b].sub[0].type = XdrDataType::Double;
1213 fr.block[b].sub[0].dval = block[b];
1217 /* check for disre block & fill it. */
1222 add_blocks_enxframe(&fr, fr.nblock);
1224 add_subblocks_enxblock(&(fr.block[db]), 2);
1225 const t_disresdata& disres = *fcd->disres;
1226 fr.block[db].id = enxDISRE;
1227 fr.block[db].sub[0].nr = ndisre;
1228 fr.block[db].sub[1].nr = ndisre;
1230 fr.block[db].sub[0].type = XdrDataType::Float;
1231 fr.block[db].sub[1].type = XdrDataType::Float;
1232 fr.block[db].sub[0].fval = disres.rt;
1233 fr.block[db].sub[1].fval = disres.rm3tav;
1235 fr.block[db].sub[0].type = XdrDataType::Double;
1236 fr.block[db].sub[1].type = XdrDataType::Double;
1237 fr.block[db].sub[0].dval = disres.rt;
1238 fr.block[db].sub[1].dval = disres.rm3tav;
1241 /* here we can put new-style blocks */
1243 /* Free energy perturbation blocks */
1246 mde_delta_h_coll_handle_block(dhc_.get(), &fr, fr.nblock);
1249 /* we can now free & reset the data in the blocks */
1252 mde_delta_h_coll_reset(dhc_.get());
1255 /* AWH bias blocks. */
1256 if (awh != nullptr) // TODO: add boolean flag.
1258 awh->writeToEnergyFrame(step, &fr);
1261 /* do the actual I/O */
1262 do_enx(fp_ene, &fr);
1265 /* We have stored the sums, so reset the sum history */
1266 reset_ebin_sums(ebin_);
1272 if (bOR && fcd->orires)
1274 print_orires_log(log, fcd->orires.get());
1277 fprintf(log, " Energies (%s)\n", unit_energy);
1278 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1283 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts)
1289 for (int i = 0; i < opts->ngtc; i++)
1291 if (opts->annealing[i] != SimulatedAnnealing::No)
1294 "Current ref_t for group %s: %8.1f\n",
1295 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1304 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1306 if (ebin_->nsum_sim <= 0)
1310 fprintf(log, "Not enough data recorded to report energy averages\n");
1317 char buf1[22], buf2[22];
1319 fprintf(log, "\t<====== ############### ==>\n");
1320 fprintf(log, "\t<==== A V E R A G E S ====>\n");
1321 fprintf(log, "\t<== ############### ======>\n\n");
1324 "\tStatistics over %s steps using %s frames\n",
1325 gmx_step_str(ebin_->nsteps_sim, buf1),
1326 gmx_step_str(ebin_->nsum_sim, buf2));
1329 fprintf(log, " Energies (%s)\n", unit_energy);
1330 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1335 pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1340 fprintf(log, " Total Virial (%s)\n", unit_energy);
1341 pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1343 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1344 pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1349 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1350 pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1356 int padding = 8 - strlen(unit_energy);
1357 fprintf(log, "%*sEpot (%s) ", padding, "", unit_energy);
1358 for (auto key : keysOf(bEInd_))
1362 fprintf(log, "%12s ", enumValueToString(key));
1368 for (int i = 0; (i < nEg_); i++)
1370 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1371 for (int j = i; (j < nEg_); j++)
1373 int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1375 14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1376 fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1377 pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1385 pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1391 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1393 const t_ebin* const ebin = ebin_;
1395 enerhist->nsteps = ebin->nsteps;
1396 enerhist->nsum = ebin->nsum;
1397 enerhist->nsteps_sim = ebin->nsteps_sim;
1398 enerhist->nsum_sim = ebin->nsum_sim;
1402 /* This will only actually resize the first time */
1403 enerhist->ener_ave.resize(ebin->nener);
1404 enerhist->ener_sum.resize(ebin->nener);
1406 for (int i = 0; i < ebin->nener; i++)
1408 enerhist->ener_ave[i] = ebin->e[i].eav;
1409 enerhist->ener_sum[i] = ebin->e[i].esum;
1413 if (ebin->nsum_sim > 0)
1415 /* This will only actually resize the first time */
1416 enerhist->ener_sum_sim.resize(ebin->nener);
1418 for (int i = 0; i < ebin->nener; i++)
1420 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1425 mde_delta_h_coll_update_energyhistory(dhc_.get(), enerhist);
1429 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1431 unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1433 if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1434 || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1437 "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1440 enerhist.ener_sum.size(),
1441 enerhist.ener_sum_sim.size());
1444 ebin_->nsteps = enerhist.nsteps;
1445 ebin_->nsum = enerhist.nsum;
1446 ebin_->nsteps_sim = enerhist.nsteps_sim;
1447 ebin_->nsum_sim = enerhist.nsum_sim;
1449 for (int i = 0; i < ebin_->nener; i++)
1451 ebin_->e[i].eav = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1452 ebin_->e[i].esum = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1453 ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1457 mde_delta_h_coll_restore_energyhistory(dhc_.get(), enerhist.deltaHForeignLambdas.get());
1461 int EnergyOutput::numEnergyTerms() const
1463 return ebin_->nener;
1466 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1468 if (fplog == nullptr)
1473 if (conservedEnergyTracker_)
1475 std::string partName = formatString("simulation part #%d", simulationPart);
1476 fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1478 else if (usingMdIntegrator)
1481 "\nCannot report drift of the conserved energy quantity because simulations share "