2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Defines code that writes energy-like quantities.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \author Paul Bauer <paul.bauer.q@gmail.com>
42 * \author Artem Zhmurov <zhmurov@gmail.com>
44 * \ingroup module_mdlib
48 #include "energyoutput.h"
57 #include "gromacs/awh/awh.h"
58 #include "gromacs/fileio/enxio.h"
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/listed_forces/disre.h"
63 #include "gromacs/listed_forces/orires.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/constr.h"
68 #include "gromacs/mdlib/ebin.h"
69 #include "gromacs/mdlib/mdebin_bar.h"
70 #include "gromacs/mdrunutility/handlerestart.h"
71 #include "gromacs/mdtypes/energyhistory.h"
72 #include "gromacs/mdtypes/fcdata.h"
73 #include "gromacs/mdtypes/group.h"
74 #include "gromacs/mdtypes/inputrec.h"
75 #include "gromacs/mdtypes/md_enums.h"
76 #include "gromacs/mdtypes/state.h"
77 #include "gromacs/pbcutil/pbc.h"
78 #include "gromacs/pulling/pull.h"
79 #include "gromacs/topology/mtop_util.h"
80 #include "gromacs/trajectory/energyframe.h"
81 #include "gromacs/utility/arraysize.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/mdmodulenotification.h"
84 #include "gromacs/utility/smalloc.h"
85 #include "gromacs/utility/stringutil.h"
87 //! Labels for energy file quantities
89 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
91 static std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
93 static std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
94 "Box-YX", "Box-ZX", "Box-ZY" };
96 static const char* vol_nm[] = { "Volume" };
98 static const char* dens_nm[] = { "Density" };
100 static const char* pv_nm[] = { "pV" };
102 static const char* enthalpy_nm[] = { "Enthalpy" };
104 static std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
105 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
107 const char* egrp_nm[egNR + 1] = { "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14", nullptr };
113 /*! \brief Energy output class
115 * This is the collection of energy averages collected during mdrun, and to
116 * be written out to the .edr file.
118 * \todo Use more std containers.
119 * \todo Remove GMX_CONSTRAINTVIR
120 * \todo Write free-energy output also to energy file (after adding more tests)
122 EnergyOutput::EnergyOutput(ener_file* fp_ene,
123 const gmx_mtop_t* mtop,
124 const t_inputrec* ir,
125 const pull_t* pull_work,
128 const StartingBehavior startingBehavior,
129 const MdModulesNotifier& mdModulesNotifier)
131 const char* ener_nm[F_NRE];
132 static const char* vir_nm[] = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
133 "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
134 static const char* sv_nm[] = { "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
135 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
136 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ" };
137 static const char* fv_nm[] = { "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
138 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
139 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ" };
140 static const char* pres_nm[] = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
141 "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
142 static const char* surft_nm[] = { "#Surf*SurfTen" };
143 static const char* mu_nm[] = { "Mu-X", "Mu-Y", "Mu-Z" };
144 static const char* vcos_nm[] = { "2CosZ*Vel-X" };
145 static const char* visc_nm[] = { "1/Viscosity" };
146 static const char* baro_nm[] = { "Barostat" };
148 const SimulationGroups* groups;
152 int i, j, ni, nj, n, k, kk, ncon, nset;
155 if (EI_DYNAMICS(ir->eI))
157 delta_t_ = ir->delta_t;
164 groups = &mtop->groups;
166 bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
167 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
169 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
170 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
171 bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
176 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
180 bConstrVir_ = (getenv("GMX_CONSTRAINTVIR") != nullptr);
187 /* Energy monitoring */
188 for (i = 0; i < egNR; i++)
193 // Setting true only to those energy terms, that have active interactions and
194 // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
195 for (i = 0; i < F_NRE; i++)
197 bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
198 && ((interaction_function[i].flags & IF_VSITE) == 0);
203 bEner_[F_EKIN] = EI_DYNAMICS(ir->eI);
204 bEner_[F_ETOT] = EI_DYNAMICS(ir->eI);
205 bEner_[F_TEMP] = EI_DYNAMICS(ir->eI);
207 bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(ir);
208 bEner_[F_PDISPCORR] = (ir->eDispCorr != edispcNO);
209 bEner_[F_PRES] = true;
212 bEner_[F_LJ] = !bBHAM;
213 bEner_[F_BHAM] = bBHAM;
214 bEner_[F_EQM] = ir->bQMMM;
215 bEner_[F_RF_EXCL] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
216 bEner_[F_COUL_RECIP] = EEL_FULL(ir->coulombtype);
217 bEner_[F_LJ_RECIP] = EVDW_PME(ir->vdwtype);
218 bEner_[F_LJ14] = b14;
219 bEner_[F_COUL14] = b14;
220 bEner_[F_LJC14_Q] = false;
221 bEner_[F_LJC_PAIRS_NB] = false;
224 bEner_[F_DVDL_COUL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptCOUL];
225 bEner_[F_DVDL_VDW] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptVDW];
226 bEner_[F_DVDL_BONDED] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptBONDED];
227 bEner_[F_DVDL_RESTRAINT] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptRESTRAINT];
228 bEner_[F_DKDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptMASS];
229 bEner_[F_DVDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptFEP];
231 bEner_[F_CONSTR] = false;
232 bEner_[F_CONSTRNC] = false;
233 bEner_[F_SETTLE] = false;
235 bEner_[F_COUL_SR] = true;
236 bEner_[F_EPOT] = true;
238 bEner_[F_DISPCORR] = (ir->eDispCorr != edispcNO);
239 bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
240 bEner_[F_ORIRESDEV] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
241 bEner_[F_COM_PULL] = ((ir->bPull && pull_have_potential(pull_work)) || ir->bRot);
243 MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
244 mdModulesNotifier.notifier_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
246 bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
249 // Counting the energy terms that will be printed and saving their names
251 for (i = 0; i < F_NRE; i++)
255 ener_nm[f_nre_] = interaction_function[i].longname;
260 epc_ = isRerun ? epcNO : ir->epc;
261 bDiagPres_ = !TRICLINIC(ir->ref_p) && !isRerun;
262 ref_p_ = (ir->ref_p[XX][XX] + ir->ref_p[YY][YY] + ir->ref_p[ZZ][ZZ]) / DIM;
263 bTricl_ = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
264 bDynBox_ = inputrecDynamicBox(ir);
265 etc_ = isRerun ? etcNO : ir->etc;
266 bNHC_trotter_ = inputrecNvtTrotter(ir) && !isRerun;
267 bPrintNHChains_ = ir->bPrintNHChains && !isRerun;
268 bMTTK_ = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
269 bMu_ = inputrecNeedMutot(ir);
273 /* Pass NULL for unit to let get_ebin_space determine the units
274 * for interaction_function[i].longname
276 ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
279 /* This should be called directly after the call for ie_,
280 * such that iconrmsd_ follows directly in the list.
282 iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
286 ib_ = get_ebin_space(ebin_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
287 bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(), unit_length);
288 ivol_ = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
289 idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
292 ipv_ = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
293 ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
298 isvir_ = get_ebin_space(ebin_, asize(sv_nm), sv_nm, unit_energy);
299 ifvir_ = get_ebin_space(ebin_, asize(fv_nm), fv_nm, unit_energy);
303 ivir_ = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
304 ipres_ = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
305 isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
307 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
309 ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
313 imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
315 if (ir->cos_accel != 0)
317 ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
318 ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
321 /* Energy monitoring */
322 for (i = 0; i < egNR; i++)
326 bEInd_[egCOULSR] = true;
327 bEInd_[egLJSR] = true;
331 bEInd_[egLJSR] = false;
332 bEInd_[egBHAMSR] = true;
336 bEInd_[egLJ14] = true;
337 bEInd_[egCOUL14] = true;
340 for (i = 0; (i < egNR); i++)
347 n = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
349 nE_ = (n * (n + 1)) / 2;
356 for (k = 0; (k < nEc_); k++)
358 snew(gnm[k], STRLEN);
360 for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
362 ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
363 for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
365 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
366 for (k = kk = 0; (k < egNR); k++)
370 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k], *(groups->groupNames[ni]),
371 *(groups->groupNames[nj]));
375 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
379 for (k = 0; (k < nEc_); k++)
387 gmx_incons("Number of energy terms wrong");
391 nTC_ = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
392 nNHC_ = ir->opts.nhchainlength; /* shorthand for number of NH chains */
395 nTCP_ = 1; /* assume only one possible coupling system for barostat
402 if (etc_ == etcNOSEHOOVER)
406 mde_n_ = 2 * nNHC_ * nTC_;
414 mdeb_n_ = 2 * nNHC_ * nTCP_;
423 snew(tmp_r_, mde_n_);
424 // TODO redo the group name memory management to make it more clear
426 snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
428 for (i = 0; (i < nTC_); i++)
430 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
431 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
432 grpnms[i] = gmx_strdup(buf);
434 itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
435 for (i = 0; i < nTC_; i++)
441 if (etc_ == etcNOSEHOOVER)
447 for (i = 0; (i < nTC_); i++)
449 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
450 bufi = *(groups->groupNames[ni]);
451 for (j = 0; (j < nNHC_); j++)
453 sprintf(buf, "Xi-%d-%s", j, bufi);
454 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
455 sprintf(buf, "vXi-%d-%s", j, bufi);
456 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
459 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
463 for (i = 0; (i < nTCP_); i++)
465 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
466 for (j = 0; (j < nNHC_); j++)
468 sprintf(buf, "Xi-%d-%s", j, bufi);
469 grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
470 sprintf(buf, "vXi-%d-%s", j, bufi);
471 grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
474 itcb_ = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
480 for (i = 0; (i < nTC_); i++)
482 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
483 bufi = *(groups->groupNames[ni]);
484 sprintf(buf, "Xi-%s", bufi);
485 grpnms[2 * i] = gmx_strdup(buf);
486 sprintf(buf, "vXi-%s", bufi);
487 grpnms[2 * i + 1] = gmx_strdup(buf);
489 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
494 else if (etc_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
496 for (i = 0; (i < nTC_); i++)
498 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
499 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
500 grpnms[i] = gmx_strdup(buf);
502 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, "");
506 for (i = 0; i < allocated; i++)
512 nU_ = groups->groups[SimulationAtomGroupType::Acceleration].size();
516 snew(grpnms, 3 * nU_);
517 for (i = 0; (i < nU_); i++)
519 ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
520 sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
521 grpnms[3 * i + XX] = gmx_strdup(buf);
522 sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
523 grpnms[3 * i + YY] = gmx_strdup(buf);
524 sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
525 grpnms[3 * i + ZZ] = gmx_strdup(buf);
527 iu_ = get_ebin_space(ebin_, 3 * nU_, grpnms, unit_vel);
528 for (i = 0; i < 3 * nU_; i++)
535 /* Note that fp_ene should be valid on the master rank and null otherwise */
536 if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
538 do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
541 /* check whether we're going to write dh histograms */
543 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
545 /* Currently dh histograms are only written with dynamics */
546 if (EI_DYNAMICS(ir->eI))
550 mde_delta_h_coll_init(dhc_, ir);
553 snew(dE_, ir->fepvals->n_lambda);
558 snew(dE_, ir->fepvals->n_lambda);
563 snew(temperatures_, ir->fepvals->n_lambda);
564 numTemperatures_ = ir->fepvals->n_lambda;
565 for (i = 0; i < ir->fepvals->n_lambda; i++)
567 temperatures_[i] = ir->simtempvals->temperatures[i];
572 numTemperatures_ = 0;
576 EnergyOutput::~EnergyOutput()
582 done_mde_delta_h_coll(dhc_);
584 if (numTemperatures_ > 0)
586 sfree(temperatures_);
592 /*! \brief Print a lambda vector to a string
594 * \param[in] fep The inputrec's FEP input data
595 * \param[in] i The index of the lambda vector
596 * \param[in] get_native_lambda Whether to print the native lambda
597 * \param[in] get_names Whether to print the names rather than the values
598 * \param[in,out] str The pre-allocated string buffer to print to.
600 static void print_lambda_vector(t_lambda* fep, int i, bool get_native_lambda, bool get_names, char* str)
605 for (j = 0; j < efptNR; j++)
607 if (fep->separate_dvdl[j])
612 str[0] = 0; /* reset the string */
615 str += sprintf(str, "("); /* set the opening parenthesis*/
617 for (j = 0; j < efptNR; j++)
619 if (fep->separate_dvdl[j])
623 if (get_native_lambda && fep->init_lambda >= 0)
625 str += sprintf(str, "%.4f", fep->init_lambda);
629 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
634 str += sprintf(str, "%s", efpt_singular_names[j]);
636 /* print comma for the next item */
639 str += sprintf(str, ", ");
646 /* and add the closing parenthesis */
651 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv)
654 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
655 *lambdastate = "\\lambda state";
656 int i, nsets, nsets_de, nsetsbegin;
657 int n_lambda_terms = 0;
658 t_lambda* fep = ir->fepvals; /* for simplicity */
659 t_expanded* expand = ir->expandedvals;
660 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
665 bool write_pV = false;
667 /* count the number of different lambda terms */
668 for (i = 0; i < efptNR; i++)
670 if (fep->separate_dvdl[i])
676 std::string title, label_x, label_y;
677 if (fep->n_lambda == 0)
679 title = gmx::formatString("%s", dhdl);
680 label_x = gmx::formatString("Time (ps)");
681 label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
685 title = gmx::formatString("%s and %s", dhdl, deltag);
686 label_x = gmx::formatString("Time (ps)");
687 label_y = gmx::formatString("%s and %s (%s %s)", dhdl, deltag, unit_energy,
688 "[\\8l\\4]\\S-1\\N");
690 fp = gmx_fio_fopen(filename, "w+");
691 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
696 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
698 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
700 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
702 /* compatibility output */
703 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
707 print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
708 print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
709 buf += gmx::formatString("%s %d: %s = %s", lambdastate, fep->init_fep_state,
710 lambda_name_str, lambda_vec_str);
713 xvgr_subtitle(fp, buf.c_str(), oenv);
717 if (fep->dhdl_derivatives == edhdlderivativesYES)
719 nsets_dhdl = n_lambda_terms;
721 /* count the number of delta_g states */
722 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
724 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
726 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
728 nsets += 1; /*add fep state for expanded ensemble */
731 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
733 nsets += 1; /* add energy to the dhdl as well */
737 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
739 nsetsextend += 1; /* for PV term, other terms possible if required for
740 the reduced potential (only needed with foreign
741 lambda, and only output when init_lambda is not
742 set in order to maintain compatibility of the
746 std::vector<std::string> setname(nsetsextend);
748 if (expand->elmcmove > elmcmoveNO)
750 /* state for the fep_vals, if we have alchemical sampling */
751 setname[s++] = "Thermodynamic state";
754 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
757 switch (fep->edHdLPrintEnergy)
759 case edHdLPrintEnergyPOTENTIAL:
760 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
762 case edHdLPrintEnergyTOTAL:
763 case edHdLPrintEnergyYES:
764 default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
766 setname[s++] = energy;
769 if (fep->dhdl_derivatives == edhdlderivativesYES)
771 for (i = 0; i < efptNR; i++)
773 if (fep->separate_dvdl[i])
775 std::string derivative;
776 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
778 /* compatibility output */
779 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
783 double lam = fep->init_lambda;
784 if (fep->init_lambda < 0)
786 lam = fep->all_lambda[i][fep->init_fep_state];
788 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i], lam);
790 setname[s++] = derivative;
795 if (fep->n_lambda > 0)
797 /* g_bar has to determine the lambda values used in this simulation
798 * from this xvg legend.
801 if (expand->elmcmove > elmcmoveNO)
803 nsetsbegin = 1; /* for including the expanded ensemble */
810 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
814 nsetsbegin += nsets_dhdl;
816 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
818 print_lambda_vector(fep, i, false, false, lambda_vec_str);
820 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
822 /* for compatible dhdl.xvg files */
823 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
827 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
832 /* print the temperature for this state if doing simulated annealing */
833 buf += gmx::formatString(
834 "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
840 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
843 xvgrLegend(fp, setname, oenv);
852 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
856 const gmx_enerdata_t* enerd,
857 const t_state* state,
859 const t_expanded* expand,
865 const gmx_ekindata_t* ekind,
867 const gmx::Constraints* constr)
869 int j, k, kk, n, gid;
870 real crmsd[2], tmp6[6];
871 real bs[tricl_boxs_nm.size()], vol, dens, pv, enthalpy;
873 double store_dhdl[efptNR];
874 real store_energy = 0;
877 /* Do NOT use the box in the state variable, but the separate box provided
878 * as an argument. This is because we sometimes need to write the box from
879 * the last timestep to match the trajectory frames.
881 add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
884 crmsd[0] = constr->rmsd();
885 add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
898 nboxs = tricl_boxs_nm.size();
905 nboxs = boxs_nm.size();
907 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
908 dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
909 add_ebin(ebin_, ib_, nboxs, bs, bSum);
910 add_ebin(ebin_, ivol_, 1, &vol, bSum);
911 add_ebin(ebin_, idens_, 1, &dens, bSum);
915 /* This is pV (in kJ/mol). The pressure is the reference pressure,
916 not the instantaneous pressure */
917 pv = vol * ref_p_ / PRESFAC;
919 add_ebin(ebin_, ipv_, 1, &pv, bSum);
920 enthalpy = pv + enerd->term[F_ETOT];
921 add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
926 add_ebin(ebin_, isvir_, 9, svir[0], bSum);
927 add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
931 add_ebin(ebin_, ivir_, 9, vir[0], bSum);
932 add_ebin(ebin_, ipres_, 9, pres[0], bSum);
933 tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
934 add_ebin(ebin_, isurft_, 1, &tmp, bSum);
936 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
938 tmp6[0] = state->boxv[XX][XX];
939 tmp6[1] = state->boxv[YY][YY];
940 tmp6[2] = state->boxv[ZZ][ZZ];
941 tmp6[3] = state->boxv[YY][XX];
942 tmp6[4] = state->boxv[ZZ][XX];
943 tmp6[5] = state->boxv[ZZ][YY];
944 add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
948 add_ebin(ebin_, imu_, 3, mu_tot, bSum);
950 if (ekind && ekind->cosacc.cos_accel != 0)
952 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
953 dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
954 add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
955 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
957 / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * PICO) * dens
958 * gmx::square(box[ZZ][ZZ] * NANO / (2 * M_PI)));
959 add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
964 for (int i = 0; (i < nEg_); i++)
966 for (j = i; (j < nEg_); j++)
968 gid = GID(i, j, nEg_);
969 for (k = kk = 0; (k < egNR); k++)
973 eee[kk++] = enerd->grpp.ener[k][gid];
976 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
984 for (int i = 0; (i < nTC_); i++)
986 tmp_r_[i] = ekind->tcstat[i].T;
988 add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
990 if (etc_ == etcNOSEHOOVER)
992 /* whether to print Nose-Hoover chains: */
997 for (int i = 0; (i < nTC_); i++)
999 for (j = 0; j < nNHC_; j++)
1002 tmp_r_[2 * k] = state->nosehoover_xi[k];
1003 tmp_r_[2 * k + 1] = state->nosehoover_vxi[k];
1006 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1010 for (int i = 0; (i < nTCP_); i++)
1012 for (j = 0; j < nNHC_; j++)
1015 tmp_r_[2 * k] = state->nhpres_xi[k];
1016 tmp_r_[2 * k + 1] = state->nhpres_vxi[k];
1019 add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1024 for (int i = 0; (i < nTC_); i++)
1026 tmp_r_[2 * i] = state->nosehoover_xi[i];
1027 tmp_r_[2 * i + 1] = state->nosehoover_vxi[i];
1029 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1033 else if (etc_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
1035 for (int i = 0; (i < nTC_); i++)
1037 tmp_r_[i] = ekind->tcstat[i].lambda;
1039 add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1043 if (ekind && nU_ > 1)
1045 for (int i = 0; (i < nU_); i++)
1047 copy_rvec(ekind->grpstat[i].u, tmp_v_[i]);
1049 add_ebin(ebin_, iu_, 3 * nU_, tmp_v_[0], bSum);
1052 ebin_increase_count(1, ebin_, bSum);
1054 // BAR + thermodynamic integration values
1055 if ((fp_dhdl_ || dhc_) && bDoDHDL)
1057 for (gmx::index i = 0; i < static_cast<gmx::index>(enerd->enerpart_lambda.size()) - 1; i++)
1059 /* zero for simulated tempering */
1060 dE_[i] = enerd->enerpart_lambda[i + 1] - enerd->enerpart_lambda[0];
1061 if (numTemperatures_ > 0)
1063 GMX_RELEASE_ASSERT(numTemperatures_ > state->fep_state,
1064 "Number of lambdas in state is bigger then in input record");
1066 numTemperatures_ >= static_cast<gmx::index>(enerd->enerpart_lambda.size()) - 1,
1067 "Number of lambdas in energy data is bigger then in input record");
1068 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1069 /* is this even useful to have at all? */
1070 dE_[i] += (temperatures_[i] / temperatures_[state->fep_state] - 1.0) * enerd->term[F_EKIN];
1076 fprintf(fp_dhdl_, "%.4f", time);
1077 /* the current free energy state */
1079 /* print the current state if we are doing expanded ensemble */
1080 if (expand->elmcmove > elmcmoveNO)
1082 fprintf(fp_dhdl_, " %4d", state->fep_state);
1084 /* total energy (for if the temperature changes */
1086 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1088 switch (fep->edHdLPrintEnergy)
1090 case edHdLPrintEnergyPOTENTIAL: store_energy = enerd->term[F_EPOT]; break;
1091 case edHdLPrintEnergyTOTAL:
1092 case edHdLPrintEnergyYES:
1093 default: store_energy = enerd->term[F_ETOT];
1095 fprintf(fp_dhdl_, " %#.8g", store_energy);
1098 if (fep->dhdl_derivatives == edhdlderivativesYES)
1100 for (int i = 0; i < efptNR; i++)
1102 if (fep->separate_dvdl[i])
1104 /* assumes F_DVDL is first */
1105 fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + i]);
1109 for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1111 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1113 if (bDynBox_ && bDiagPres_ && (epc_ != epcNO) && !enerd->enerpart_lambda.empty()
1114 && (fep->init_lambda < 0))
1116 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1117 there are alternate state
1118 lambda and we're not in
1119 compatibility mode */
1121 fprintf(fp_dhdl_, "\n");
1122 /* and the binary free energy output */
1124 if (dhc_ && bDoDHDL)
1127 for (int i = 0; i < efptNR; i++)
1129 if (fep->separate_dvdl[i])
1131 /* assumes F_DVDL is first */
1132 store_dhdl[idhdl] = enerd->term[F_DVDL + i];
1136 store_energy = enerd->term[F_ETOT];
1137 /* store_dh is dE */
1138 mde_delta_h_coll_add_dh(dhc_, static_cast<double>(state->fep_state), store_energy, pv,
1139 store_dhdl, dE_ + fep->lambda_start_n, time);
1144 void EnergyOutput::recordNonEnergyStep()
1146 ebin_increase_count(1, ebin_, false);
1149 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1156 "Step", "Time", gmx_step_str(steps, buf), time);
1159 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1173 fr.nsteps = ebin_->nsteps;
1175 fr.nsum = ebin_->nsum;
1176 fr.nre = (bEne) ? ebin_->nener : 0;
1178 int ndisre = bDR ? fcd->disres.npair : 0;
1179 /* these are for the old-style blocks (1 subblock, only reals), because
1180 there can be only one per ID for these */
1184 /* Optional additional old-style (real-only) blocks. */
1185 for (int i = 0; i < enxNR; i++)
1190 if (bOR && fcd->orires.nr > 0)
1192 diagonalize_orires_tensors(&(fcd->orires));
1193 nr[enxOR] = fcd->orires.nr;
1194 block[enxOR] = fcd->orires.otav;
1196 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ? fcd->orires.nr : 0;
1197 block[enxORI] = fcd->orires.oinsl;
1198 id[enxORI] = enxORI;
1199 nr[enxORT] = fcd->orires.nex * 12;
1200 block[enxORT] = fcd->orires.eig;
1201 id[enxORT] = enxORT;
1204 /* whether we are going to write anything out: */
1205 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1207 /* the old-style blocks go first */
1209 for (int i = 0; i < enxNR; i++)
1216 add_blocks_enxframe(&fr, fr.nblock);
1217 for (int b = 0; b < fr.nblock; b++)
1219 add_subblocks_enxblock(&(fr.block[b]), 1);
1220 fr.block[b].id = id[b];
1221 fr.block[b].sub[0].nr = nr[b];
1223 fr.block[b].sub[0].type = xdr_datatype_float;
1224 fr.block[b].sub[0].fval = block[b];
1226 fr.block[b].sub[0].type = xdr_datatype_double;
1227 fr.block[b].sub[0].dval = block[b];
1231 /* check for disre block & fill it. */
1236 add_blocks_enxframe(&fr, fr.nblock);
1238 add_subblocks_enxblock(&(fr.block[db]), 2);
1239 fr.block[db].id = enxDISRE;
1240 fr.block[db].sub[0].nr = ndisre;
1241 fr.block[db].sub[1].nr = ndisre;
1243 fr.block[db].sub[0].type = xdr_datatype_float;
1244 fr.block[db].sub[1].type = xdr_datatype_float;
1245 fr.block[db].sub[0].fval = fcd->disres.rt;
1246 fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1248 fr.block[db].sub[0].type = xdr_datatype_double;
1249 fr.block[db].sub[1].type = xdr_datatype_double;
1250 fr.block[db].sub[0].dval = fcd->disres.rt;
1251 fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1254 /* here we can put new-style blocks */
1256 /* Free energy perturbation blocks */
1259 mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1262 /* we can now free & reset the data in the blocks */
1265 mde_delta_h_coll_reset(dhc_);
1268 /* AWH bias blocks. */
1269 if (awh != nullptr) // TODO: add boolean flag.
1271 awh->writeToEnergyFrame(step, &fr);
1274 /* do the actual I/O */
1275 do_enx(fp_ene, &fr);
1278 /* We have stored the sums, so reset the sum history */
1279 reset_ebin_sums(ebin_);
1285 if (bOR && fcd->orires.nr > 0)
1287 print_orires_log(log, &(fcd->orires));
1290 fprintf(log, " Energies (%s)\n", unit_energy);
1291 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1296 void EnergyOutput::printAnnealingTemperatures(FILE* log, SimulationGroups* groups, t_grpopts* opts)
1302 for (int i = 0; i < opts->ngtc; i++)
1304 if (opts->annealing[i] != eannNO)
1306 fprintf(log, "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");
1335 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1336 gmx_step_str(ebin_->nsteps_sim, buf1), gmx_step_str(ebin_->nsum_sim, buf2));
1339 fprintf(log, " Energies (%s)\n", unit_energy);
1340 pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1345 pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1350 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1351 pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1353 fprintf(log, " Force Virial (%s)\n", unit_energy);
1354 pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1359 fprintf(log, " Total Virial (%s)\n", unit_energy);
1360 pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1362 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1363 pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1368 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1369 pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1375 int padding = 8 - strlen(unit_energy);
1376 fprintf(log, "%*sEpot (%s) ", padding, "", unit_energy);
1377 for (int i = 0; (i < egNR); i++)
1381 fprintf(log, "%12s ", egrp_nm[i]);
1387 for (int i = 0; (i < nEg_); i++)
1389 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1390 for (int j = i; (j < nEg_); j++)
1392 int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1394 14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1395 fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]),
1396 *(groups->groupNames[nj]));
1397 pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1405 pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1410 fprintf(log, "%15s %12s %12s %12s\n", "Group", "Ux", "Uy", "Uz");
1411 for (int i = 0; (i < nU_); i++)
1413 int ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
1414 fprintf(log, "%15s", *groups->groupNames[ni]);
1415 pr_ebin(log, ebin_, iu_ + 3 * i, 3, 3, eprAVER, false);
1422 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1424 const t_ebin* const ebin = ebin_;
1426 enerhist->nsteps = ebin->nsteps;
1427 enerhist->nsum = ebin->nsum;
1428 enerhist->nsteps_sim = ebin->nsteps_sim;
1429 enerhist->nsum_sim = ebin->nsum_sim;
1433 /* This will only actually resize the first time */
1434 enerhist->ener_ave.resize(ebin->nener);
1435 enerhist->ener_sum.resize(ebin->nener);
1437 for (int i = 0; i < ebin->nener; i++)
1439 enerhist->ener_ave[i] = ebin->e[i].eav;
1440 enerhist->ener_sum[i] = ebin->e[i].esum;
1444 if (ebin->nsum_sim > 0)
1446 /* This will only actually resize the first time */
1447 enerhist->ener_sum_sim.resize(ebin->nener);
1449 for (int i = 0; i < ebin->nener; i++)
1451 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1456 mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1460 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1462 unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1464 if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1465 || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1468 "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1470 nener, enerhist.ener_sum.size(), enerhist.ener_sum_sim.size());
1473 ebin_->nsteps = enerhist.nsteps;
1474 ebin_->nsum = enerhist.nsum;
1475 ebin_->nsteps_sim = enerhist.nsteps_sim;
1476 ebin_->nsum_sim = enerhist.nsum_sim;
1478 for (int i = 0; i < ebin_->nener; i++)
1480 ebin_->e[i].eav = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1481 ebin_->e[i].esum = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1482 ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1486 mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1490 int EnergyOutput::numEnergyTerms() const
1492 return ebin_->nener;