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/fatalerror.h"
84 #include "gromacs/utility/mdmodulenotification.h"
85 #include "gromacs/utility/smalloc.h"
86 #include "gromacs/utility/stringutil.h"
88 #include "energydrifttracker.h"
90 //! Labels for energy file quantities
92 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
94 static std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
96 static std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
97 "Box-YX", "Box-ZX", "Box-ZY" };
99 static const char* vol_nm[] = { "Volume" };
101 static const char* dens_nm[] = { "Density" };
103 static const char* pv_nm[] = { "pV" };
105 static const char* enthalpy_nm[] = { "Enthalpy" };
107 static std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
108 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
110 const char* egrp_nm[egNR + 1] = { "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14", nullptr };
116 /*! \brief Energy output class
118 * This is the collection of energy averages collected during mdrun, and to
119 * be written out to the .edr file.
121 * \todo Use more std containers.
122 * \todo Remove GMX_CONSTRAINTVIR
123 * \todo Write free-energy output also to energy file (after adding more tests)
125 EnergyOutput::EnergyOutput(ener_file* fp_ene,
126 const gmx_mtop_t* mtop,
127 const t_inputrec* ir,
128 const pull_t* pull_work,
131 const StartingBehavior startingBehavior,
132 const bool simulationsShareState,
133 const MdModulesNotifier& mdModulesNotifier)
135 const char* ener_nm[F_NRE];
136 static const char* vir_nm[] = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
137 "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
138 static const char* sv_nm[] = { "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
139 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
140 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ" };
141 static const char* fv_nm[] = { "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
142 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
143 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ" };
144 static const char* pres_nm[] = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
145 "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
146 static const char* surft_nm[] = { "#Surf*SurfTen" };
147 static const char* mu_nm[] = { "Mu-X", "Mu-Y", "Mu-Z" };
148 static const char* vcos_nm[] = { "2CosZ*Vel-X" };
149 static const char* visc_nm[] = { "1/Viscosity" };
150 static const char* baro_nm[] = { "Barostat" };
152 const SimulationGroups* groups;
156 int i, j, ni, nj, n, k, kk, ncon, nset;
159 if (EI_DYNAMICS(ir->eI))
161 delta_t_ = ir->delta_t;
168 groups = &mtop->groups;
170 bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
171 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
173 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
174 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
175 bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
180 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
184 bConstrVir_ = (getenv("GMX_CONSTRAINTVIR") != nullptr);
191 /* Energy monitoring */
192 for (i = 0; i < egNR; i++)
197 // Setting true only to those energy terms, that have active interactions and
198 // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
199 for (i = 0; i < F_NRE; i++)
201 bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
202 && ((interaction_function[i].flags & IF_VSITE) == 0);
207 bEner_[F_EKIN] = EI_DYNAMICS(ir->eI);
208 bEner_[F_ETOT] = EI_DYNAMICS(ir->eI);
209 bEner_[F_TEMP] = EI_DYNAMICS(ir->eI);
211 bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(ir);
212 bEner_[F_PDISPCORR] = (ir->eDispCorr != edispcNO);
213 bEner_[F_PRES] = true;
216 bEner_[F_LJ] = !bBHAM;
217 bEner_[F_BHAM] = bBHAM;
218 bEner_[F_EQM] = ir->bQMMM;
219 bEner_[F_RF_EXCL] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
220 bEner_[F_COUL_RECIP] = EEL_FULL(ir->coulombtype);
221 bEner_[F_LJ_RECIP] = EVDW_PME(ir->vdwtype);
222 bEner_[F_LJ14] = b14;
223 bEner_[F_COUL14] = b14;
224 bEner_[F_LJC14_Q] = false;
225 bEner_[F_LJC_PAIRS_NB] = false;
228 bEner_[F_DVDL_COUL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptCOUL];
229 bEner_[F_DVDL_VDW] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptVDW];
230 bEner_[F_DVDL_BONDED] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptBONDED];
231 bEner_[F_DVDL_RESTRAINT] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptRESTRAINT];
232 bEner_[F_DKDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptMASS];
233 bEner_[F_DVDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptFEP];
235 bEner_[F_CONSTR] = false;
236 bEner_[F_CONSTRNC] = false;
237 bEner_[F_SETTLE] = false;
239 bEner_[F_COUL_SR] = true;
240 bEner_[F_EPOT] = true;
242 bEner_[F_DISPCORR] = (ir->eDispCorr != edispcNO);
243 bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
244 bEner_[F_ORIRESDEV] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
245 bEner_[F_COM_PULL] = ((ir->bPull && pull_have_potential(*pull_work)) || ir->bRot);
247 MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
248 mdModulesNotifier.simulationSetupNotifications_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
250 bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
253 // Counting the energy terms that will be printed and saving their names
255 for (i = 0; i < F_NRE; i++)
259 ener_nm[f_nre_] = interaction_function[i].longname;
264 epc_ = isRerun ? epcNO : ir->epc;
265 bDiagPres_ = !TRICLINIC(ir->ref_p) && !isRerun;
266 ref_p_ = (ir->ref_p[XX][XX] + ir->ref_p[YY][YY] + ir->ref_p[ZZ][ZZ]) / DIM;
267 bTricl_ = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
268 bDynBox_ = inputrecDynamicBox(ir);
269 etc_ = isRerun ? etcNO : ir->etc;
270 bNHC_trotter_ = inputrecNvtTrotter(ir) && !isRerun;
271 bPrintNHChains_ = ir->bPrintNHChains && !isRerun;
272 bMTTK_ = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
273 bMu_ = inputrecNeedMutot(ir);
277 /* Pass NULL for unit to let get_ebin_space determine the units
278 * for interaction_function[i].longname
280 ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
283 /* This should be called directly after the call for ie_,
284 * such that iconrmsd_ follows directly in the list.
286 iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
290 ib_ = get_ebin_space(ebin_,
291 bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
292 bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
294 ivol_ = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
295 idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
298 ipv_ = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
299 ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
304 isvir_ = get_ebin_space(ebin_, asize(sv_nm), sv_nm, unit_energy);
305 ifvir_ = get_ebin_space(ebin_, asize(fv_nm), fv_nm, unit_energy);
309 ivir_ = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
310 ipres_ = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
311 isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
313 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
315 ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
319 imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
321 if (ir->cos_accel != 0)
323 ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
324 ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
327 /* Energy monitoring */
328 for (i = 0; i < egNR; i++)
332 bEInd_[egCOULSR] = true;
333 bEInd_[egLJSR] = true;
337 bEInd_[egLJSR] = false;
338 bEInd_[egBHAMSR] = true;
342 bEInd_[egLJ14] = true;
343 bEInd_[egCOUL14] = true;
346 for (i = 0; (i < egNR); i++)
353 n = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
355 nE_ = (n * (n + 1)) / 2;
362 for (k = 0; (k < nEc_); k++)
364 snew(gnm[k], STRLEN);
366 for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
368 ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
369 for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
371 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
372 for (k = kk = 0; (k < egNR); k++)
379 *(groups->groupNames[ni]),
380 *(groups->groupNames[nj]));
384 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
388 for (k = 0; (k < nEc_); k++)
396 gmx_incons("Number of energy terms wrong");
400 nTC_ = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
401 nNHC_ = ir->opts.nhchainlength; /* shorthand for number of NH chains */
404 nTCP_ = 1; /* assume only one possible coupling system for barostat
411 if (etc_ == etcNOSEHOOVER)
415 mde_n_ = 2 * nNHC_ * nTC_;
423 mdeb_n_ = 2 * nNHC_ * nTCP_;
432 snew(tmp_r_, mde_n_);
433 // TODO redo the group name memory management to make it more clear
435 snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
437 for (i = 0; (i < nTC_); i++)
439 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
440 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
441 grpnms[i] = gmx_strdup(buf);
443 itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
444 for (i = 0; i < nTC_; i++)
450 if (etc_ == etcNOSEHOOVER)
456 for (i = 0; (i < nTC_); i++)
458 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
459 bufi = *(groups->groupNames[ni]);
460 for (j = 0; (j < nNHC_); j++)
462 sprintf(buf, "Xi-%d-%s", j, bufi);
463 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
464 sprintf(buf, "vXi-%d-%s", j, bufi);
465 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
468 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
472 for (i = 0; (i < nTCP_); i++)
474 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
475 for (j = 0; (j < nNHC_); j++)
477 sprintf(buf, "Xi-%d-%s", j, bufi);
478 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
479 sprintf(buf, "vXi-%d-%s", j, bufi);
480 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
483 itcb_ = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
489 for (i = 0; (i < nTC_); i++)
491 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
492 bufi = *(groups->groupNames[ni]);
493 sprintf(buf, "Xi-%s", bufi);
494 grpnms[2 * i] = gmx_strdup(buf);
495 sprintf(buf, "vXi-%s", bufi);
496 grpnms[2 * i + 1] = gmx_strdup(buf);
498 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
503 else if (etc_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
505 for (i = 0; (i < nTC_); i++)
507 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
508 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
509 grpnms[i] = gmx_strdup(buf);
511 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, "");
515 for (i = 0; i < allocated; i++)
521 /* Note that fp_ene should be valid on the master rank and null otherwise */
522 if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
524 do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
527 /* check whether we're going to write dh histograms */
529 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
531 /* Currently dh histograms are only written with dynamics */
532 if (EI_DYNAMICS(ir->eI))
536 mde_delta_h_coll_init(dhc_, ir);
539 snew(dE_, ir->fepvals->n_lambda);
544 snew(dE_, ir->fepvals->n_lambda);
549 snew(temperatures_, ir->fepvals->n_lambda);
550 numTemperatures_ = ir->fepvals->n_lambda;
551 for (i = 0; i < ir->fepvals->n_lambda; i++)
553 temperatures_[i] = ir->simtempvals->temperatures[i];
558 numTemperatures_ = 0;
561 if (EI_MD(ir->eI) && !simulationsShareState)
563 conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop->natoms);
567 EnergyOutput::~EnergyOutput()
573 done_mde_delta_h_coll(dhc_);
575 if (numTemperatures_ > 0)
577 sfree(temperatures_);
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 (j = 0; j < efptNR; j++)
598 if (fep->separate_dvdl[j])
603 str[0] = 0; /* reset the string */
606 str += sprintf(str, "("); /* set the opening parenthesis*/
608 for (j = 0; j < efptNR; j++)
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", efpt_singular_names[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; /* for simplicity */
650 t_expanded* expand = ir->expandedvals;
651 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
656 bool write_pV = false;
658 /* count the number of different lambda terms */
659 for (i = 0; i < efptNR; i++)
661 if (fep->separate_dvdl[i])
667 std::string title, label_x, label_y;
668 if (fep->n_lambda == 0)
670 title = gmx::formatString("%s", dhdl);
671 label_x = gmx::formatString("Time (ps)");
672 label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
676 title = gmx::formatString("%s and %s", dhdl, deltag);
677 label_x = gmx::formatString("Time (ps)");
678 label_y = gmx::formatString(
679 "%s and %s (%s %s)", dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
681 fp = gmx_fio_fopen(filename, "w+");
682 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
687 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
689 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
691 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
693 /* compatibility output */
694 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
698 print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
699 print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
700 buf += gmx::formatString(
701 "%s %d: %s = %s", lambdastate, fep->init_fep_state, lambda_name_str, lambda_vec_str);
704 xvgr_subtitle(fp, buf.c_str(), oenv);
708 if (fep->dhdl_derivatives == edhdlderivativesYES)
710 nsets_dhdl = n_lambda_terms;
712 /* count the number of delta_g states */
713 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
715 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
717 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
719 nsets += 1; /*add fep state for expanded ensemble */
722 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
724 nsets += 1; /* add energy to the dhdl as well */
728 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
730 nsetsextend += 1; /* for PV term, other terms possible if required for
731 the reduced potential (only needed with foreign
732 lambda, and only output when init_lambda is not
733 set in order to maintain compatibility of the
737 std::vector<std::string> setname(nsetsextend);
739 if (expand->elmcmove > elmcmoveNO)
741 /* state for the fep_vals, if we have alchemical sampling */
742 setname[s++] = "Thermodynamic state";
745 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
748 switch (fep->edHdLPrintEnergy)
750 case edHdLPrintEnergyPOTENTIAL:
751 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
753 case edHdLPrintEnergyTOTAL:
754 case edHdLPrintEnergyYES:
755 default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
757 setname[s++] = energy;
760 if (fep->dhdl_derivatives == edhdlderivativesYES)
762 for (i = 0; i < efptNR; i++)
764 if (fep->separate_dvdl[i])
766 std::string derivative;
767 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
769 /* compatibility output */
770 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
774 double lam = fep->init_lambda;
775 if (fep->init_lambda < 0)
777 lam = fep->all_lambda[i][fep->init_fep_state];
779 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i], lam);
781 setname[s++] = derivative;
786 if (fep->n_lambda > 0)
788 /* g_bar has to determine the lambda values used in this simulation
789 * from this xvg legend.
792 if (expand->elmcmove > elmcmoveNO)
794 nsetsbegin = 1; /* for including the expanded ensemble */
801 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
805 nsetsbegin += nsets_dhdl;
807 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
809 print_lambda_vector(fep, i, false, false, lambda_vec_str);
811 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
813 /* for compatible dhdl.xvg files */
814 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
818 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
823 /* print the temperature for this state if doing simulated annealing */
824 buf += gmx::formatString(
825 "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
831 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
834 xvgrLegend(fp, setname, oenv);
843 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
847 const gmx_enerdata_t* enerd,
849 const t_expanded* expand,
851 PTCouplingArrays ptCouplingArrays,
857 const gmx_ekindata_t* ekind,
859 const gmx::Constraints* constr)
861 int j, k, kk, n, gid;
862 real crmsd[2], tmp6[6];
863 real bs[tricl_boxs_nm.size()], vol, dens, pv, enthalpy;
865 double store_dhdl[efptNR];
866 real store_energy = 0;
869 /* Do NOT use the box in the state variable, but the separate box provided
870 * as an argument. This is because we sometimes need to write the box from
871 * the last timestep to match the trajectory frames.
873 add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
876 crmsd[0] = constr->rmsd();
877 add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
890 nboxs = tricl_boxs_nm.size();
897 nboxs = boxs_nm.size();
899 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
900 dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
901 add_ebin(ebin_, ib_, nboxs, bs, bSum);
902 add_ebin(ebin_, ivol_, 1, &vol, bSum);
903 add_ebin(ebin_, idens_, 1, &dens, bSum);
907 /* This is pV (in kJ/mol). The pressure is the reference pressure,
908 not the instantaneous pressure */
909 pv = vol * ref_p_ / PRESFAC;
911 add_ebin(ebin_, ipv_, 1, &pv, bSum);
912 enthalpy = pv + enerd->term[F_ETOT];
913 add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
918 add_ebin(ebin_, isvir_, 9, svir[0], bSum);
919 add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
923 add_ebin(ebin_, ivir_, 9, vir[0], bSum);
924 add_ebin(ebin_, ipres_, 9, pres[0], bSum);
925 tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
926 add_ebin(ebin_, isurft_, 1, &tmp, bSum);
928 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
930 tmp6[0] = ptCouplingArrays.boxv[XX][XX];
931 tmp6[1] = ptCouplingArrays.boxv[YY][YY];
932 tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
933 tmp6[3] = ptCouplingArrays.boxv[YY][XX];
934 tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
935 tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
936 add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
940 add_ebin(ebin_, imu_, 3, mu_tot, bSum);
942 if (ekind && ekind->cosacc.cos_accel != 0)
944 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
945 dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
946 add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
947 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
949 / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * PICO) * dens
950 * gmx::square(box[ZZ][ZZ] * NANO / (2 * M_PI)));
951 add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
956 for (int i = 0; (i < nEg_); i++)
958 for (j = i; (j < nEg_); j++)
960 gid = GID(i, j, nEg_);
961 for (k = kk = 0; (k < egNR); k++)
965 eee[kk++] = enerd->grpp.ener[k][gid];
968 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
976 for (int i = 0; (i < nTC_); i++)
978 tmp_r_[i] = ekind->tcstat[i].T;
980 add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
982 if (etc_ == etcNOSEHOOVER)
984 /* whether to print Nose-Hoover chains: */
989 for (int i = 0; (i < nTC_); i++)
991 for (j = 0; j < nNHC_; j++)
994 tmp_r_[2 * k] = ptCouplingArrays.nosehoover_xi[k];
995 tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
998 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1002 for (int i = 0; (i < nTCP_); i++)
1004 for (j = 0; j < nNHC_; j++)
1007 tmp_r_[2 * k] = ptCouplingArrays.nhpres_xi[k];
1008 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
1011 add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1016 for (int i = 0; (i < nTC_); i++)
1018 tmp_r_[2 * i] = ptCouplingArrays.nosehoover_xi[i];
1019 tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1021 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1025 else if (etc_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
1027 for (int i = 0; (i < nTC_); i++)
1029 tmp_r_[i] = ekind->tcstat[i].lambda;
1031 add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1035 ebin_increase_count(1, ebin_, bSum);
1037 // BAR + thermodynamic integration values
1038 if ((fp_dhdl_ || dhc_) && bDoDHDL)
1040 const auto& foreignTerms = enerd->foreignLambdaTerms;
1041 for (int i = 0; i < foreignTerms.numLambdas(); i++)
1043 /* zero for simulated tempering */
1044 dE_[i] = foreignTerms.deltaH(i);
1045 if (numTemperatures_ > 0)
1047 GMX_RELEASE_ASSERT(numTemperatures_ > fep_state,
1048 "Number of lambdas in state is bigger then in input record");
1050 numTemperatures_ >= foreignTerms.numLambdas(),
1051 "Number of lambdas in energy data is bigger then in input record");
1052 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1053 /* is this even useful to have at all? */
1054 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1060 fprintf(fp_dhdl_, "%.4f", time);
1061 /* the current free energy state */
1063 /* print the current state if we are doing expanded ensemble */
1064 if (expand->elmcmove > elmcmoveNO)
1066 fprintf(fp_dhdl_, " %4d", fep_state);
1068 /* total energy (for if the temperature changes */
1070 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1072 switch (fep->edHdLPrintEnergy)
1074 case edHdLPrintEnergyPOTENTIAL: store_energy = enerd->term[F_EPOT]; break;
1075 case edHdLPrintEnergyTOTAL:
1076 case edHdLPrintEnergyYES:
1077 default: store_energy = enerd->term[F_ETOT];
1079 fprintf(fp_dhdl_, " %#.8g", store_energy);
1082 if (fep->dhdl_derivatives == edhdlderivativesYES)
1084 for (int i = 0; i < efptNR; i++)
1086 if (fep->separate_dvdl[i])
1088 /* assumes F_DVDL is first */
1089 fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + i]);
1093 for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1095 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1097 if (bDynBox_ && bDiagPres_ && (epc_ != epcNO) && foreignTerms.numLambdas() > 0
1098 && (fep->init_lambda < 0))
1100 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1101 there are alternate state
1102 lambda and we're not in
1103 compatibility mode */
1105 fprintf(fp_dhdl_, "\n");
1106 /* and the binary free energy output */
1108 if (dhc_ && bDoDHDL)
1111 for (int i = 0; i < efptNR; i++)
1113 if (fep->separate_dvdl[i])
1115 /* assumes F_DVDL is first */
1116 store_dhdl[idhdl] = enerd->term[F_DVDL + i];
1120 store_energy = enerd->term[F_ETOT];
1121 /* store_dh is dE */
1122 mde_delta_h_coll_add_dh(
1123 dhc_, static_cast<double>(fep_state), store_energy, pv, store_dhdl, dE_ + fep->lambda_start_n, time);
1127 if (conservedEnergyTracker_)
1129 conservedEnergyTracker_->addPoint(
1130 time, bEner_[F_ECONSERVED] ? enerd->term[F_ECONSERVED] : enerd->term[F_ETOT]);
1134 void EnergyOutput::recordNonEnergyStep()
1136 ebin_increase_count(1, ebin_, false);
1139 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1148 gmx_step_str(steps, buf),
1152 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1166 fr.nsteps = ebin_->nsteps;
1168 fr.nsum = ebin_->nsum;
1169 fr.nre = (bEne) ? ebin_->nener : 0;
1171 int ndisre = bDR ? fcd->disres->npair : 0;
1172 /* these are for the old-style blocks (1 subblock, only reals), because
1173 there can be only one per ID for these */
1177 /* Optional additional old-style (real-only) blocks. */
1178 for (int i = 0; i < enxNR; i++)
1183 if (bOR && fcd->orires->nr > 0)
1185 t_oriresdata& orires = *fcd->orires;
1186 diagonalize_orires_tensors(&orires);
1187 nr[enxOR] = orires.nr;
1188 block[enxOR] = orires.otav;
1190 nr[enxORI] = (orires.oinsl != orires.otav) ? orires.nr : 0;
1191 block[enxORI] = orires.oinsl;
1192 id[enxORI] = enxORI;
1193 nr[enxORT] = orires.nex * 12;
1194 block[enxORT] = orires.eig;
1195 id[enxORT] = enxORT;
1198 /* whether we are going to write anything out: */
1199 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1201 /* the old-style blocks go first */
1203 for (int i = 0; i < enxNR; i++)
1210 add_blocks_enxframe(&fr, fr.nblock);
1211 for (int b = 0; b < fr.nblock; b++)
1213 add_subblocks_enxblock(&(fr.block[b]), 1);
1214 fr.block[b].id = id[b];
1215 fr.block[b].sub[0].nr = nr[b];
1217 fr.block[b].sub[0].type = xdr_datatype_float;
1218 fr.block[b].sub[0].fval = block[b];
1220 fr.block[b].sub[0].type = xdr_datatype_double;
1221 fr.block[b].sub[0].dval = block[b];
1225 /* check for disre block & fill it. */
1230 add_blocks_enxframe(&fr, fr.nblock);
1232 add_subblocks_enxblock(&(fr.block[db]), 2);
1233 const t_disresdata& disres = *fcd->disres;
1234 fr.block[db].id = enxDISRE;
1235 fr.block[db].sub[0].nr = ndisre;
1236 fr.block[db].sub[1].nr = ndisre;
1238 fr.block[db].sub[0].type = xdr_datatype_float;
1239 fr.block[db].sub[1].type = xdr_datatype_float;
1240 fr.block[db].sub[0].fval = disres.rt;
1241 fr.block[db].sub[1].fval = disres.rm3tav;
1243 fr.block[db].sub[0].type = xdr_datatype_double;
1244 fr.block[db].sub[1].type = xdr_datatype_double;
1245 fr.block[db].sub[0].dval = disres.rt;
1246 fr.block[db].sub[1].dval = disres.rm3tav;
1249 /* here we can put new-style blocks */
1251 /* Free energy perturbation blocks */
1254 mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1257 /* we can now free & reset the data in the blocks */
1260 mde_delta_h_coll_reset(dhc_);
1263 /* AWH bias blocks. */
1264 if (awh != nullptr) // TODO: add boolean flag.
1266 awh->writeToEnergyFrame(step, &fr);
1269 /* do the actual I/O */
1270 do_enx(fp_ene, &fr);
1273 /* We have stored the sums, so reset the sum history */
1274 reset_ebin_sums(ebin_);
1280 if (bOR && fcd->orires->nr > 0)
1282 print_orires_log(log, fcd->orires);
1285 fprintf(log, " Energies (%s)\n", unit_energy);
1286 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1291 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts)
1297 for (int i = 0; i < opts->ngtc; i++)
1299 if (opts->annealing[i] != eannNO)
1302 "Current ref_t for group %s: %8.1f\n",
1303 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1312 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1314 if (ebin_->nsum_sim <= 0)
1318 fprintf(log, "Not enough data recorded to report energy averages\n");
1325 char buf1[22], buf2[22];
1327 fprintf(log, "\t<====== ############### ==>\n");
1328 fprintf(log, "\t<==== A V E R A G E S ====>\n");
1329 fprintf(log, "\t<== ############### ======>\n\n");
1332 "\tStatistics over %s steps using %s frames\n",
1333 gmx_step_str(ebin_->nsteps_sim, buf1),
1334 gmx_step_str(ebin_->nsum_sim, buf2));
1337 fprintf(log, " Energies (%s)\n", unit_energy);
1338 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1343 pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1348 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1349 pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1351 fprintf(log, " Force Virial (%s)\n", unit_energy);
1352 pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1357 fprintf(log, " Total Virial (%s)\n", unit_energy);
1358 pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1360 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1361 pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1366 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1367 pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1373 int padding = 8 - strlen(unit_energy);
1374 fprintf(log, "%*sEpot (%s) ", padding, "", unit_energy);
1375 for (int i = 0; (i < egNR); i++)
1379 fprintf(log, "%12s ", egrp_nm[i]);
1385 for (int i = 0; (i < nEg_); i++)
1387 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1388 for (int j = i; (j < nEg_); j++)
1390 int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1392 14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1393 fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1394 pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1402 pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1408 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1410 const t_ebin* const ebin = ebin_;
1412 enerhist->nsteps = ebin->nsteps;
1413 enerhist->nsum = ebin->nsum;
1414 enerhist->nsteps_sim = ebin->nsteps_sim;
1415 enerhist->nsum_sim = ebin->nsum_sim;
1419 /* This will only actually resize the first time */
1420 enerhist->ener_ave.resize(ebin->nener);
1421 enerhist->ener_sum.resize(ebin->nener);
1423 for (int i = 0; i < ebin->nener; i++)
1425 enerhist->ener_ave[i] = ebin->e[i].eav;
1426 enerhist->ener_sum[i] = ebin->e[i].esum;
1430 if (ebin->nsum_sim > 0)
1432 /* This will only actually resize the first time */
1433 enerhist->ener_sum_sim.resize(ebin->nener);
1435 for (int i = 0; i < ebin->nener; i++)
1437 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1442 mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1446 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1448 unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1450 if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1451 || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1454 "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1457 enerhist.ener_sum.size(),
1458 enerhist.ener_sum_sim.size());
1461 ebin_->nsteps = enerhist.nsteps;
1462 ebin_->nsum = enerhist.nsum;
1463 ebin_->nsteps_sim = enerhist.nsteps_sim;
1464 ebin_->nsum_sim = enerhist.nsum_sim;
1466 for (int i = 0; i < ebin_->nener; i++)
1468 ebin_->e[i].eav = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1469 ebin_->e[i].esum = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1470 ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1474 mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1478 int EnergyOutput::numEnergyTerms() const
1480 return ebin_->nener;
1483 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1485 if (fplog == nullptr)
1490 if (conservedEnergyTracker_)
1492 std::string partName = formatString("simulation part #%d", simulationPart);
1493 fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1495 else if (usingMdIntegrator)
1498 "\nCannot report drift of the conserved energy quantity because simulations share "