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/applied_forces/awh/read_params.h"
60 #include "gromacs/fileio/enxio.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/fileio/xvgr.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/listed_forces/disre.h"
65 #include "gromacs/listed_forces/orires.h"
66 #include "gromacs/math/functions.h"
67 #include "gromacs/math/units.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/ebin.h"
71 #include "gromacs/mdlib/mdebin_bar.h"
72 #include "gromacs/mdrunutility/handlerestart.h"
73 #include "gromacs/mdtypes/energyhistory.h"
74 #include "gromacs/mdtypes/fcdata.h"
75 #include "gromacs/mdtypes/group.h"
76 #include "gromacs/mdtypes/inputrec.h"
77 #include "gromacs/mdtypes/md_enums.h"
78 #include "gromacs/mdtypes/state.h"
79 #include "gromacs/pbcutil/pbc.h"
80 #include "gromacs/pulling/pull.h"
81 #include "gromacs/topology/mtop_util.h"
82 #include "gromacs/trajectory/energyframe.h"
83 #include "gromacs/utility/arraysize.h"
84 #include "gromacs/utility/enumerationhelpers.h"
85 #include "gromacs/utility/fatalerror.h"
86 #include "gromacs/utility/mdmodulesnotifiers.h"
87 #include "gromacs/utility/smalloc.h"
88 #include "gromacs/utility/stringutil.h"
90 #include "energydrifttracker.h"
92 //! Labels for energy file quantities
94 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
95 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
97 static constexpr std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
99 static constexpr std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
100 "Box-YX", "Box-ZX", "Box-ZY" };
102 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
103 static const char* vol_nm[] = { "Volume" };
105 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
106 static const char* dens_nm[] = { "Density" };
108 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
109 static const char* pv_nm[] = { "pV" };
111 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
112 static const char* enthalpy_nm[] = { "Enthalpy" };
114 static constexpr std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
115 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
117 const char* enumValueToString(NonBondedEnergyTerms enumValue)
119 static constexpr gmx::EnumerationArray<NonBondedEnergyTerms, const char*> nonBondedEnergyTermTypeNames = {
120 "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14"
122 return nonBondedEnergyTermTypeNames[enumValue];
127 static bool haveFepLambdaMoves(const t_inputrec& inputrec)
129 return (inputrec.bExpanded && inputrec.expandedvals->elmcmove > LambdaMoveCalculation::No)
130 || (inputrec.efep != FreeEnergyPerturbationType::No && inputrec.bDoAwh
131 && awhHasFepLambdaDimension(*inputrec.awhParams));
137 /*! \brief Energy output class
139 * This is the collection of energy averages collected during mdrun, and to
140 * be written out to the .edr file.
142 * \todo Use more std containers.
143 * \todo Write free-energy output also to energy file (after adding more tests)
145 EnergyOutput::EnergyOutput(ener_file* fp_ene,
146 const gmx_mtop_t& mtop,
147 const t_inputrec& inputrec,
148 const pull_t* pull_work,
151 const StartingBehavior startingBehavior,
152 const bool simulationsShareState,
153 const MDModulesNotifiers& mdModulesNotifiers) :
154 haveFepLambdaMoves_(haveFepLambdaMoves(inputrec))
156 const char* ener_nm[F_NRE];
157 static const char* vir_nm[] = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
158 "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
159 static const char* pres_nm[] = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
160 "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
161 static const char* surft_nm[] = { "#Surf*SurfTen" };
162 static const char* mu_nm[] = { "Mu-X", "Mu-Y", "Mu-Z" };
163 static const char* vcos_nm[] = { "2CosZ*Vel-X" };
164 static const char* visc_nm[] = { "1/Viscosity" };
165 static const char* baro_nm[] = { "Barostat" };
167 const SimulationGroups* groups;
171 int i, j, ni, nj, n, ncon, nset;
174 if (EI_DYNAMICS(inputrec.eI))
176 delta_t_ = inputrec.delta_t;
183 groups = &mtop.groups;
185 bBHAM = (mtop.ffparams.numTypes() > 0) && (mtop.ffparams.functype[0] == F_BHAM);
186 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
188 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
189 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
190 bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
194 if (ncon > 0 && inputrec.eConstrAlg == ConstraintAlgorithm::Lincs)
204 /* Energy monitoring */
205 for (auto& term : bEInd_)
210 // Setting true only to those energy terms, that have active interactions and
211 // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
212 for (i = 0; i < F_NRE; i++)
214 bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
215 && ((interaction_function[i].flags & IF_VSITE) == 0);
220 bEner_[F_EKIN] = EI_DYNAMICS(inputrec.eI);
221 bEner_[F_ETOT] = EI_DYNAMICS(inputrec.eI);
222 bEner_[F_TEMP] = EI_DYNAMICS(inputrec.eI);
224 bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(&inputrec);
225 bEner_[F_PDISPCORR] = (inputrec.eDispCorr != DispersionCorrectionType::No);
226 bEner_[F_PRES] = true;
229 bEner_[F_LJ] = !bBHAM;
230 bEner_[F_BHAM] = bBHAM;
231 bEner_[F_RF_EXCL] = (EEL_RF(inputrec.coulombtype) && inputrec.cutoff_scheme == CutoffScheme::Group);
232 bEner_[F_COUL_RECIP] = EEL_FULL(inputrec.coulombtype);
233 bEner_[F_LJ_RECIP] = EVDW_PME(inputrec.vdwtype);
234 bEner_[F_LJ14] = b14;
235 bEner_[F_COUL14] = b14;
236 bEner_[F_LJC14_Q] = false;
237 bEner_[F_LJC_PAIRS_NB] = false;
240 bEner_[F_DVDL_COUL] = (inputrec.efep != FreeEnergyPerturbationType::No)
241 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Coul];
242 bEner_[F_DVDL_VDW] = (inputrec.efep != FreeEnergyPerturbationType::No)
243 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Vdw];
244 bEner_[F_DVDL_BONDED] = (inputrec.efep != FreeEnergyPerturbationType::No)
245 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Bonded];
246 bEner_[F_DVDL_RESTRAINT] =
247 (inputrec.efep != FreeEnergyPerturbationType::No)
248 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Restraint];
249 bEner_[F_DKDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
250 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Mass];
251 bEner_[F_DVDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
252 && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep];
254 bEner_[F_CONSTR] = false;
255 bEner_[F_CONSTRNC] = false;
256 bEner_[F_SETTLE] = false;
258 bEner_[F_COUL_SR] = true;
259 bEner_[F_EPOT] = true;
261 bEner_[F_DISPCORR] = (inputrec.eDispCorr != DispersionCorrectionType::No);
262 bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
263 bEner_[F_ORIRESDEV] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
264 bEner_[F_COM_PULL] = ((inputrec.bPull && pull_have_potential(*pull_work)) || inputrec.bRot);
266 // Check MDModules for any energy output
267 MDModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
268 mdModulesNotifiers.simulationSetupNotifier_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
270 bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
272 MDModulesEnergyOutputToQMMMRequestChecker mdModulesAddOutputToQMMMFieldRequest;
273 mdModulesNotifiers.simulationSetupNotifier_.notify(&mdModulesAddOutputToQMMMFieldRequest);
275 bEner_[F_EQM] = mdModulesAddOutputToQMMMFieldRequest.energyOutputToQMMM_;
277 // Counting the energy terms that will be printed and saving their names
279 for (i = 0; i < F_NRE; i++)
283 ener_nm[f_nre_] = interaction_function[i].longname;
288 epc_ = isRerun ? PressureCoupling::No : inputrec.epc;
289 bDiagPres_ = !TRICLINIC(inputrec.ref_p) && !isRerun;
290 ref_p_ = (inputrec.ref_p[XX][XX] + inputrec.ref_p[YY][YY] + inputrec.ref_p[ZZ][ZZ]) / DIM;
291 bTricl_ = TRICLINIC(inputrec.compress) || TRICLINIC(inputrec.deform);
292 bDynBox_ = inputrecDynamicBox(&inputrec);
293 etc_ = isRerun ? TemperatureCoupling::No : inputrec.etc;
294 bNHC_trotter_ = inputrecNvtTrotter(&inputrec) && !isRerun;
295 bPrintNHChains_ = inputrec.bPrintNHChains && !isRerun;
296 bMTTK_ = (inputrecNptTrotter(&inputrec) || inputrecNphTrotter(&inputrec)) && !isRerun;
297 bMu_ = inputrecNeedMutot(&inputrec);
301 /* Pass NULL for unit to let get_ebin_space determine the units
302 * for interaction_function[i].longname
304 ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
307 /* This should be called directly after the call for ie_,
308 * such that iconrmsd_ follows directly in the list.
310 iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
314 ib_ = get_ebin_space(ebin_,
315 bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
316 bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
318 ivol_ = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
319 idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
322 ipv_ = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
323 ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
328 ivir_ = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
329 ipres_ = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
330 isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
332 if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
334 ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
338 imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
340 if (inputrec.cos_accel != 0)
342 ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
343 ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
346 /* Energy monitoring */
347 for (auto& term : bEInd_)
351 bEInd_[NonBondedEnergyTerms::CoulombSR] = true;
352 bEInd_[NonBondedEnergyTerms::LJSR] = true;
356 bEInd_[NonBondedEnergyTerms::LJSR] = false;
357 bEInd_[NonBondedEnergyTerms::BuckinghamSR] = true;
361 bEInd_[NonBondedEnergyTerms::LJ14] = true;
362 bEInd_[NonBondedEnergyTerms::Coulomb14] = true;
365 for (auto term : bEInd_)
372 n = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
374 nE_ = (n * (n + 1)) / 2;
381 for (int k = 0; (k < nEc_); k++)
383 snew(gnm[k], STRLEN);
385 for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
387 ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
388 for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
390 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
392 for (auto key : keysOf(bEInd_))
398 enumValueToString(key),
399 *(groups->groupNames[ni]),
400 *(groups->groupNames[nj]));
404 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
408 for (int k = 0; (k < nEc_); k++)
416 gmx_incons("Number of energy terms wrong");
420 nTC_ = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
421 nNHC_ = inputrec.opts.nhchainlength; /* shorthand for number of NH chains */
424 nTCP_ = 1; /* assume only one possible coupling system for barostat
431 if (etc_ == TemperatureCoupling::NoseHoover)
435 mde_n_ = 2 * nNHC_ * nTC_;
441 if (epc_ == PressureCoupling::Mttk)
443 mdeb_n_ = 2 * nNHC_ * nTCP_;
452 tmp_r_.resize(mde_n_);
453 // TODO redo the group name memory management to make it more clear
455 snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
457 for (i = 0; (i < nTC_); i++)
459 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
460 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
461 grpnms[i] = gmx_strdup(buf);
463 itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
464 for (i = 0; i < nTC_; i++)
470 if (etc_ == TemperatureCoupling::NoseHoover)
476 for (i = 0; (i < nTC_); i++)
478 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
479 bufi = *(groups->groupNames[ni]);
480 for (j = 0; (j < nNHC_); j++)
482 sprintf(buf, "Xi-%d-%s", j, bufi);
483 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
484 sprintf(buf, "vXi-%d-%s", j, bufi);
485 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
488 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
492 for (i = 0; (i < nTCP_); i++)
494 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
495 for (j = 0; (j < nNHC_); j++)
497 sprintf(buf, "Xi-%d-%s", j, bufi);
498 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
499 sprintf(buf, "vXi-%d-%s", j, bufi);
500 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
503 itcb_ = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
509 for (i = 0; (i < nTC_); i++)
511 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
512 bufi = *(groups->groupNames[ni]);
513 sprintf(buf, "Xi-%s", bufi);
514 grpnms[2 * i] = gmx_strdup(buf);
515 sprintf(buf, "vXi-%s", bufi);
516 grpnms[2 * i + 1] = gmx_strdup(buf);
518 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
523 else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
524 || etc_ == TemperatureCoupling::VRescale)
526 for (i = 0; (i < nTC_); i++)
528 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
529 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
530 grpnms[i] = gmx_strdup(buf);
532 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, "");
536 for (i = 0; i < allocated; i++)
542 /* Note that fp_ene should be valid on the master rank and null otherwise */
543 if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
545 do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
548 /* check whether we're going to write dh histograms */
550 if (inputrec.fepvals->separate_dhdl_file == SeparateDhdlFile::No)
552 /* Currently dh histograms are only written with dynamics */
553 if (EI_DYNAMICS(inputrec.eI))
555 dhc_ = std::make_unique<t_mde_delta_h_coll>(inputrec);
558 dE_.resize(inputrec.fepvals->n_lambda);
563 dE_.resize(inputrec.fepvals->n_lambda);
565 if (inputrec.bSimTemp)
567 temperatures_ = inputrec.simtempvals->temperatures;
570 if (EI_MD(inputrec.eI) && !simulationsShareState)
572 conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop.natoms);
576 EnergyOutput::~EnergyOutput()
583 /*! \brief Print a lambda vector to a string
585 * \param[in] fep The inputrec's FEP input data
586 * \param[in] i The index of the lambda vector
587 * \param[in] get_native_lambda Whether to print the native lambda
588 * \param[in] get_names Whether to print the names rather than the values
589 * \param[in,out] str The pre-allocated string buffer to print to.
591 static void print_lambda_vector(t_lambda* fep, int i, bool get_native_lambda, bool get_names, char* str)
596 for (auto j : keysOf(fep->separate_dvdl))
598 if (fep->separate_dvdl[j])
603 str[0] = 0; /* reset the string */
606 str += sprintf(str, "("); /* set the opening parenthesis*/
608 for (auto j : keysOf(fep->separate_dvdl))
610 if (fep->separate_dvdl[j])
614 if (get_native_lambda && fep->init_lambda >= 0)
616 str += sprintf(str, "%.4f", fep->init_lambda);
620 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
625 str += sprintf(str, "%s", enumValueToStringSingular(j));
627 /* print comma for the next item */
630 str += sprintf(str, ", ");
637 /* and add the closing parenthesis */
642 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv)
645 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
646 *lambdastate = "\\lambda state";
647 int i, nsets, nsets_de, nsetsbegin;
648 int n_lambda_terms = 0;
649 t_lambda* fep = ir->fepvals.get(); /* for simplicity */
650 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
655 bool write_pV = false;
657 /* count the number of different lambda terms */
658 for (auto i : keysOf(fep->separate_dvdl))
660 if (fep->separate_dvdl[i])
666 std::string title, label_x, label_y;
667 if (fep->n_lambda == 0)
669 title = gmx::formatString("%s", dhdl);
670 label_x = gmx::formatString("Time (ps)");
671 label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
675 title = gmx::formatString("%s and %s", dhdl, deltag);
676 label_x = gmx::formatString("Time (ps)");
677 label_y = gmx::formatString(
678 "%s and %s (%s %s)", dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
680 fp = gmx_fio_fopen(filename, "w+");
681 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
686 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
688 if ((ir->efep != FreeEnergyPerturbationType::SlowGrowth)
689 && (ir->efep != FreeEnergyPerturbationType::Expanded)
690 && !(ir->bDoAwh && awhHasFepLambdaDimension(*ir->awhParams)))
692 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
694 /* compatibility output */
695 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
699 print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
700 print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
701 buf += gmx::formatString(
702 "%s %d: %s = %s", lambdastate, fep->init_fep_state, lambda_name_str, lambda_vec_str);
705 xvgr_subtitle(fp, buf.c_str(), oenv);
709 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
711 nsets_dhdl = n_lambda_terms;
713 /* count the number of delta_g states */
714 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
716 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
718 if (haveFepLambdaMoves(*ir))
720 nsets += 1; /*add fep state for expanded ensemble */
723 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
725 nsets += 1; /* add energy to the dhdl as well */
729 if ((ir->epc != PressureCoupling::No) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
731 nsetsextend += 1; /* for PV term, other terms possible if required for
732 the reduced potential (only needed with foreign
733 lambda, and only output when init_lambda is not
734 set in order to maintain compatibility of the
738 std::vector<std::string> setname(nsetsextend);
740 if (haveFepLambdaMoves(*ir))
742 /* state for the fep_vals, if we have alchemical sampling */
743 setname[s++] = "Thermodynamic state";
746 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
749 switch (fep->edHdLPrintEnergy)
751 case FreeEnergyPrintEnergy::Potential:
752 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
754 case FreeEnergyPrintEnergy::Total:
755 case FreeEnergyPrintEnergy::Yes:
756 default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
758 setname[s++] = energy;
761 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
763 for (auto i : keysOf(fep->separate_dvdl))
765 if (fep->separate_dvdl[i])
767 std::string derivative;
768 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
770 /* compatibility output */
771 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
775 double lam = fep->init_lambda;
776 if (fep->init_lambda < 0)
778 lam = fep->all_lambda[i][fep->init_fep_state];
780 derivative = gmx::formatString("%s %s = %.4f", dhdl, enumValueToStringSingular(i), lam);
782 setname[s++] = derivative;
787 if (fep->n_lambda > 0)
789 /* g_bar has to determine the lambda values used in this simulation
790 * from this xvg legend.
793 if (haveFepLambdaMoves(*ir))
795 nsetsbegin = 1; /* for including the expanded ensemble */
802 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
806 nsetsbegin += nsets_dhdl;
808 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
810 print_lambda_vector(fep, i, false, false, lambda_vec_str);
812 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
814 /* for compatible dhdl.xvg files */
815 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
819 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
824 /* print the temperature for this state if doing simulated annealing */
825 buf += gmx::formatString(
826 "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
832 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
835 xvgrLegend(fp, setname, oenv);
844 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
848 const gmx_enerdata_t* enerd,
851 PTCouplingArrays ptCouplingArrays,
855 const gmx_ekindata_t* ekind,
857 const gmx::Constraints* constr)
859 int j, k, kk, n, gid;
860 real crmsd[2], tmp6[6];
861 real bs[tricl_boxs_nm.size()], vol, dens, enthalpy;
862 real eee[static_cast<int>(NonBondedEnergyTerms::Count)];
863 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, double> store_dhdl;
864 real store_energy = 0;
866 real pv = 0.0; // static analyzer warns about uninitialized variable warnings here.
868 /* Do NOT use the box in the state variable, but the separate box provided
869 * as an argument. This is because we sometimes need to write the box from
870 * the last timestep to match the trajectory frames.
872 add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
875 crmsd[0] = constr->rmsd();
876 add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
889 nboxs = tricl_boxs_nm.size();
896 nboxs = boxs_nm.size();
898 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
899 dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
900 add_ebin(ebin_, ib_, nboxs, bs, bSum);
901 add_ebin(ebin_, ivol_, 1, &vol, bSum);
902 add_ebin(ebin_, idens_, 1, &dens, bSum);
906 /* This is pV (in kJ/mol). The pressure is the reference pressure,
907 not the instantaneous pressure */
908 pv = vol * ref_p_ / gmx::c_presfac;
910 add_ebin(ebin_, ipv_, 1, &pv, bSum);
911 enthalpy = pv + enerd->term[F_ETOT];
912 add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
917 add_ebin(ebin_, ivir_, 9, vir[0], bSum);
918 add_ebin(ebin_, ipres_, 9, pres[0], bSum);
919 tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
920 add_ebin(ebin_, isurft_, 1, &tmp, bSum);
922 if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
924 tmp6[0] = ptCouplingArrays.boxv[XX][XX];
925 tmp6[1] = ptCouplingArrays.boxv[YY][YY];
926 tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
927 tmp6[3] = ptCouplingArrays.boxv[YY][XX];
928 tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
929 tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
930 add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
934 add_ebin(ebin_, imu_, 3, mu_tot, bSum);
936 if (ekind && ekind->cosacc.cos_accel != 0)
938 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
939 dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
940 add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
941 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
943 / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * gmx::c_pico) * dens
944 * gmx::square(box[ZZ][ZZ] * gmx::c_nano / (2 * M_PI)));
945 add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
950 for (int i = 0; (i < nEg_); i++)
952 for (j = i; (j < nEg_); j++)
954 gid = GID(i, j, nEg_);
955 for (k = kk = 0; (k < static_cast<int>(NonBondedEnergyTerms::Count)); k++)
959 eee[kk++] = enerd->grpp.energyGroupPairTerms[k][gid];
962 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
970 for (int i = 0; (i < nTC_); i++)
972 tmp_r_[i] = ekind->tcstat[i].T;
974 add_ebin(ebin_, itemp_, nTC_, tmp_r_.data(), bSum);
976 if (etc_ == TemperatureCoupling::NoseHoover)
978 /* whether to print Nose-Hoover chains: */
983 for (int i = 0; (i < nTC_); i++)
985 for (j = 0; j < nNHC_; j++)
988 tmp_r_[2 * k] = ptCouplingArrays.nosehoover_xi[k];
989 tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
992 add_ebin(ebin_, itc_, mde_n_, tmp_r_.data(), bSum);
996 for (int i = 0; (i < nTCP_); i++)
998 for (j = 0; j < nNHC_; j++)
1001 tmp_r_[2 * k] = ptCouplingArrays.nhpres_xi[k];
1002 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
1005 add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_.data(), bSum);
1010 for (int i = 0; (i < nTC_); i++)
1012 tmp_r_[2 * i] = ptCouplingArrays.nosehoover_xi[i];
1013 tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1015 add_ebin(ebin_, itc_, mde_n_, tmp_r_.data(), bSum);
1019 else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
1020 || etc_ == TemperatureCoupling::VRescale)
1022 for (int i = 0; (i < nTC_); i++)
1024 tmp_r_[i] = ekind->tcstat[i].lambda;
1026 add_ebin(ebin_, itc_, nTC_, tmp_r_.data(), bSum);
1030 ebin_increase_count(1, ebin_, bSum);
1032 // BAR + thermodynamic integration values
1033 if ((fp_dhdl_ || dhc_) && bDoDHDL)
1035 const auto& foreignTerms = enerd->foreignLambdaTerms;
1036 for (int i = 0; i < foreignTerms.numLambdas(); i++)
1038 /* zero for simulated tempering */
1039 dE_[i] = foreignTerms.deltaH(i);
1040 if (!temperatures_.empty())
1042 GMX_RELEASE_ASSERT(gmx::ssize(temperatures_) > fep_state,
1043 "Number of lambdas in state is bigger then in input record");
1045 gmx::ssize(temperatures_) >= foreignTerms.numLambdas(),
1046 "Number of lambdas in energy data is bigger then in input record");
1047 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1048 /* is this even useful to have at all? */
1049 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1055 fprintf(fp_dhdl_, "%.4f", time);
1056 /* the current free energy state */
1058 /* print the current state if we are doing expanded ensemble */
1059 if (haveFepLambdaMoves_)
1061 fprintf(fp_dhdl_, " %4d", fep_state);
1063 /* total energy (for if the temperature changes */
1065 if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
1067 switch (fep->edHdLPrintEnergy)
1069 case FreeEnergyPrintEnergy::Potential:
1070 store_energy = enerd->term[F_EPOT];
1072 case FreeEnergyPrintEnergy::Total:
1073 case FreeEnergyPrintEnergy::Yes:
1074 default: store_energy = enerd->term[F_ETOT];
1076 fprintf(fp_dhdl_, " %#.8g", store_energy);
1079 if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
1081 for (auto i : keysOf(fep->separate_dvdl))
1083 if (fep->separate_dvdl[i])
1085 /* assumes F_DVDL is first */
1086 fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + static_cast<int>(i)]);
1090 for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1092 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1094 if (bDynBox_ && bDiagPres_ && (epc_ != PressureCoupling::No)
1095 && foreignTerms.numLambdas() > 0 && (fep->init_lambda < 0))
1097 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1098 there are alternate state
1099 lambda and we're not in
1100 compatibility mode */
1102 fprintf(fp_dhdl_, "\n");
1103 /* and the binary free energy output */
1105 if (dhc_ && bDoDHDL)
1108 for (auto i : keysOf(fep->separate_dvdl))
1110 if (fep->separate_dvdl[i])
1112 /* assumes F_DVDL is first */
1113 store_dhdl[idhdl] = enerd->term[F_DVDL + static_cast<int>(i)];
1117 store_energy = enerd->term[F_ETOT];
1118 /* store_dh is dE */
1119 mde_delta_h_coll_add_dh(dhc_.get(),
1120 static_cast<double>(fep_state),
1124 dE_.data() + fep->lambda_start_n,
1129 if (conservedEnergyTracker_)
1131 conservedEnergyTracker_->addPoint(
1132 time, bEner_[F_ECONSERVED] ? enerd->term[F_ECONSERVED] : enerd->term[F_ETOT]);
1136 void EnergyOutput::recordNonEnergyStep()
1138 ebin_increase_count(1, ebin_, false);
1141 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1150 gmx_step_str(steps, buf),
1154 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1168 fr.nsteps = ebin_->nsteps;
1170 fr.nsum = ebin_->nsum;
1171 fr.nre = (bEne) ? ebin_->nener : 0;
1173 int ndisre = bDR ? fcd->disres->npair : 0;
1174 /* these are for the old-style blocks (1 subblock, only reals), because
1175 there can be only one per ID for these */
1179 /* Optional additional old-style (real-only) blocks. */
1180 for (int i = 0; i < enxNR; i++)
1185 if (bOR && fcd->orires)
1187 t_oriresdata& orires = *fcd->orires;
1188 diagonalize_orires_tensors(&orires);
1189 nr[enxOR] = orires.numRestraints;
1190 block[enxOR] = orires.orientationsTimeAndEnsembleAv.data();
1192 nr[enxORI] = (orires.orientations.data() != orires.orientationsTimeAndEnsembleAv.data())
1193 ? orires.numRestraints
1195 block[enxORI] = orires.orientations.data();
1196 id[enxORI] = enxORI;
1197 nr[enxORT] = ssize(orires.eigenOutput);
1198 block[enxORT] = orires.eigenOutput.data();
1199 id[enxORT] = enxORT;
1202 /* whether we are going to write anything out: */
1203 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1205 /* the old-style blocks go first */
1207 for (int i = 0; i < enxNR; i++)
1214 add_blocks_enxframe(&fr, fr.nblock);
1215 for (int b = 0; b < fr.nblock; b++)
1217 add_subblocks_enxblock(&(fr.block[b]), 1);
1218 fr.block[b].id = id[b];
1219 fr.block[b].sub[0].nr = nr[b];
1221 fr.block[b].sub[0].type = XdrDataType::Float;
1222 fr.block[b].sub[0].fval = block[b];
1224 fr.block[b].sub[0].type = XdrDataType::Double;
1225 fr.block[b].sub[0].dval = block[b];
1229 /* check for disre block & fill it. */
1234 add_blocks_enxframe(&fr, fr.nblock);
1236 add_subblocks_enxblock(&(fr.block[db]), 2);
1237 const t_disresdata& disres = *fcd->disres;
1238 fr.block[db].id = enxDISRE;
1239 fr.block[db].sub[0].nr = ndisre;
1240 fr.block[db].sub[1].nr = ndisre;
1242 fr.block[db].sub[0].type = XdrDataType::Float;
1243 fr.block[db].sub[1].type = XdrDataType::Float;
1244 fr.block[db].sub[0].fval = disres.rt;
1245 fr.block[db].sub[1].fval = disres.rm3tav;
1247 fr.block[db].sub[0].type = XdrDataType::Double;
1248 fr.block[db].sub[1].type = XdrDataType::Double;
1249 fr.block[db].sub[0].dval = disres.rt;
1250 fr.block[db].sub[1].dval = disres.rm3tav;
1253 /* here we can put new-style blocks */
1255 /* Free energy perturbation blocks */
1258 mde_delta_h_coll_handle_block(dhc_.get(), &fr, fr.nblock);
1261 /* we can now free & reset the data in the blocks */
1264 mde_delta_h_coll_reset(dhc_.get());
1267 /* AWH bias blocks. */
1268 if (awh != nullptr) // TODO: add boolean flag.
1270 awh->writeToEnergyFrame(step, &fr);
1273 /* do the actual I/O */
1274 do_enx(fp_ene, &fr);
1277 /* We have stored the sums, so reset the sum history */
1278 reset_ebin_sums(ebin_);
1284 if (bOR && fcd->orires)
1286 print_orires_log(log, fcd->orires.get());
1289 fprintf(log, " Energies (%s)\n", unit_energy);
1290 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1295 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts)
1301 for (int i = 0; i < opts->ngtc; i++)
1303 if (opts->annealing[i] != SimulatedAnnealing::No)
1306 "Current ref_t for group %s: %8.1f\n",
1307 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1316 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1318 if (ebin_->nsum_sim <= 0)
1322 fprintf(log, "Not enough data recorded to report energy averages\n");
1329 char buf1[22], buf2[22];
1331 fprintf(log, "\t<====== ############### ==>\n");
1332 fprintf(log, "\t<==== A V E R A G E S ====>\n");
1333 fprintf(log, "\t<== ############### ======>\n\n");
1336 "\tStatistics over %s steps using %s frames\n",
1337 gmx_step_str(ebin_->nsteps_sim, buf1),
1338 gmx_step_str(ebin_->nsum_sim, buf2));
1341 fprintf(log, " Energies (%s)\n", unit_energy);
1342 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1347 pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1352 fprintf(log, " Total Virial (%s)\n", unit_energy);
1353 pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1355 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1356 pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1361 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1362 pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1368 int padding = 8 - strlen(unit_energy);
1369 fprintf(log, "%*sEpot (%s) ", padding, "", unit_energy);
1370 for (auto key : keysOf(bEInd_))
1374 fprintf(log, "%12s ", enumValueToString(key));
1380 for (int i = 0; (i < nEg_); i++)
1382 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1383 for (int j = i; (j < nEg_); j++)
1385 int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1387 14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1388 fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1389 pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1397 pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1403 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1405 const t_ebin* const ebin = ebin_;
1407 enerhist->nsteps = ebin->nsteps;
1408 enerhist->nsum = ebin->nsum;
1409 enerhist->nsteps_sim = ebin->nsteps_sim;
1410 enerhist->nsum_sim = ebin->nsum_sim;
1414 /* This will only actually resize the first time */
1415 enerhist->ener_ave.resize(ebin->nener);
1416 enerhist->ener_sum.resize(ebin->nener);
1418 for (int i = 0; i < ebin->nener; i++)
1420 enerhist->ener_ave[i] = ebin->e[i].eav;
1421 enerhist->ener_sum[i] = ebin->e[i].esum;
1425 if (ebin->nsum_sim > 0)
1427 /* This will only actually resize the first time */
1428 enerhist->ener_sum_sim.resize(ebin->nener);
1430 for (int i = 0; i < ebin->nener; i++)
1432 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1437 mde_delta_h_coll_update_energyhistory(dhc_.get(), enerhist);
1441 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1443 unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1445 if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1446 || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1449 "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1452 enerhist.ener_sum.size(),
1453 enerhist.ener_sum_sim.size());
1456 ebin_->nsteps = enerhist.nsteps;
1457 ebin_->nsum = enerhist.nsum;
1458 ebin_->nsteps_sim = enerhist.nsteps_sim;
1459 ebin_->nsum_sim = enerhist.nsum_sim;
1461 for (int i = 0; i < ebin_->nener; i++)
1463 ebin_->e[i].eav = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1464 ebin_->e[i].esum = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1465 ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1469 mde_delta_h_coll_restore_energyhistory(dhc_.get(), enerhist.deltaHForeignLambdas.get());
1473 int EnergyOutput::numEnergyTerms() const
1475 return ebin_->nener;
1478 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1480 if (fplog == nullptr)
1485 if (conservedEnergyTracker_)
1487 std::string partName = formatString("simulation part #%d", simulationPart);
1488 fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1490 else if (usingMdIntegrator)
1493 "\nCannot report drift of the conserved energy quantity because simulations share "