2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
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/mdtypes/energyhistory.h"
71 #include "gromacs/mdtypes/fcdata.h"
72 #include "gromacs/mdtypes/group.h"
73 #include "gromacs/mdtypes/inputrec.h"
74 #include "gromacs/mdtypes/md_enums.h"
75 #include "gromacs/mdtypes/state.h"
76 #include "gromacs/pbcutil/pbc.h"
77 #include "gromacs/pulling/pull.h"
78 #include "gromacs/topology/mtop_util.h"
79 #include "gromacs/trajectory/energyframe.h"
80 #include "gromacs/utility/arraysize.h"
81 #include "gromacs/utility/fatalerror.h"
82 #include "gromacs/utility/mdmodulenotification.h"
83 #include "gromacs/utility/smalloc.h"
84 #include "gromacs/utility/stringutil.h"
86 //! Labels for energy file quantities
88 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
90 static std::array<const char *, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
92 static std::array<const char *, 6> tricl_boxs_nm = {
93 "Box-XX", "Box-YY", "Box-ZZ",
94 "Box-YX", "Box-ZX", "Box-ZY"
97 static const char *vol_nm[] = { "Volume" };
99 static const char *dens_nm[] = {"Density" };
101 static const char *pv_nm[] = {"pV" };
103 static const char *enthalpy_nm[] = {"Enthalpy" };
105 static std::array<const char *, 6> boxvel_nm = {
106 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
107 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
110 const char *egrp_nm[egNR+1] = {
111 "Coul-SR", "LJ-SR", "Buck-SR",
112 "Coul-14", "LJ-14", nullptr
119 /*! \brief Energy output class
121 * This is the collection of energy averages collected during mdrun, and to
122 * be written out to the .edr file.
124 * \todo Use more std containers.
125 * \todo Remove GMX_CONSTRAINTVIR
126 * \todo Write free-energy output also to energy file (after adding more tests)
128 EnergyOutput::EnergyOutput(ener_file * fp_ene,
129 const gmx_mtop_t * mtop,
130 const t_inputrec * ir,
131 const pull_t * pull_work,
134 const MdModulesNotifier &mdModulesNotifier)
136 const char *ener_nm[F_NRE];
137 static const char *vir_nm[] = {
138 "Vir-XX", "Vir-XY", "Vir-XZ",
139 "Vir-YX", "Vir-YY", "Vir-YZ",
140 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
142 static const char *sv_nm[] = {
143 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
144 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
145 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
147 static const char *fv_nm[] = {
148 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
149 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
150 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
152 static const char *pres_nm[] = {
153 "Pres-XX", "Pres-XY", "Pres-XZ",
154 "Pres-YX", "Pres-YY", "Pres-YZ",
155 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
157 static const char *surft_nm[] = {
160 static const char *mu_nm[] = {
161 "Mu-X", "Mu-Y", "Mu-Z"
163 static const char *vcos_nm[] = {
166 static const char *visc_nm[] = {
169 static const char *baro_nm[] = {
173 const SimulationGroups *groups;
177 int i, j, ni, nj, n, k, kk, ncon, nset;
180 if (EI_DYNAMICS(ir->eI))
182 delta_t_ = ir->delta_t;
189 groups = &mtop->groups;
191 bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
192 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
193 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
195 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
196 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
197 bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
202 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
206 bConstrVir_ = (getenv("GMX_CONSTRAINTVIR") != nullptr);
213 /* Energy monitoring */
214 for (i = 0; i < egNR; i++)
219 // Setting true only to those energy terms, that have active interactions and
220 // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
221 for (i = 0; i < F_NRE; i++)
223 bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0) &&
224 ((interaction_function[i].flags & IF_VSITE) == 0);
229 bEner_[F_EKIN] = EI_DYNAMICS(ir->eI);
230 bEner_[F_ETOT] = EI_DYNAMICS(ir->eI);
231 bEner_[F_TEMP] = EI_DYNAMICS(ir->eI);
233 bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(ir);
234 bEner_[F_PDISPCORR] = (ir->eDispCorr != edispcNO);
235 bEner_[F_PRES] = true;
238 bEner_[F_LJ] = !bBHAM;
239 bEner_[F_BHAM] = bBHAM;
240 bEner_[F_EQM] = ir->bQMMM;
241 bEner_[F_RF_EXCL] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
242 bEner_[F_COUL_RECIP] = EEL_FULL(ir->coulombtype);
243 bEner_[F_LJ_RECIP] = EVDW_PME(ir->vdwtype);
244 bEner_[F_LJ14] = b14;
245 bEner_[F_COUL14] = b14;
246 bEner_[F_LJC14_Q] = false;
247 bEner_[F_LJC_PAIRS_NB] = false;
250 bEner_[F_DVDL_COUL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptCOUL];
251 bEner_[F_DVDL_VDW] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptVDW];
252 bEner_[F_DVDL_BONDED] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptBONDED];
253 bEner_[F_DVDL_RESTRAINT] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptRESTRAINT];
254 bEner_[F_DKDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptMASS];
255 bEner_[F_DVDL] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptFEP];
257 bEner_[F_CONSTR] = false;
258 bEner_[F_CONSTRNC] = false;
259 bEner_[F_SETTLE] = false;
261 bEner_[F_COUL_SR] = true;
262 bEner_[F_EPOT] = true;
264 bEner_[F_DISPCORR] = (ir->eDispCorr != edispcNO);
265 bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
266 bEner_[F_ORIRESDEV] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
267 bEner_[F_COM_PULL] = ((ir->bPull && pull_have_potential(pull_work)) || ir->bRot);
269 MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
270 mdModulesNotifier.notifier_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
272 bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
275 // Counting the energy terms that will be printed and saving their names
277 for (i = 0; i < F_NRE; i++)
281 ener_nm[f_nre_] = interaction_function[i].longname;
286 epc_ = isRerun ? epcNO : ir->epc;
287 bDiagPres_ = !TRICLINIC(ir->ref_p) && !isRerun;
288 ref_p_ = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
289 bTricl_ = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
290 bDynBox_ = inputrecDynamicBox(ir);
291 etc_ = isRerun ? etcNO : ir->etc;
292 bNHC_trotter_ = inputrecNvtTrotter(ir) && !isRerun;
293 bPrintNHChains_ = ir->bPrintNHChains && !isRerun;
294 bMTTK_ = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
295 bMu_ = inputrecNeedMutot(ir);
299 /* Pass NULL for unit to let get_ebin_space determine the units
300 * for interaction_function[i].longname
302 ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
305 /* This should be called directly after the call for ie_,
306 * such that iconrmsd_ follows directly in the list.
308 iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
312 ib_ = get_ebin_space(ebin_,
313 bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
314 bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
316 ivol_ = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
317 idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
320 ipv_ = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
321 ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
326 isvir_ = get_ebin_space(ebin_, asize(sv_nm), sv_nm, unit_energy);
327 ifvir_ = get_ebin_space(ebin_, asize(fv_nm), fv_nm, unit_energy);
331 ivir_ = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
332 ipres_ = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
333 isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm,
336 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
338 ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM,
339 boxvel_nm.data(), unit_vel);
343 imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
345 if (ir->cos_accel != 0)
347 ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
348 ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm,
352 /* Energy monitoring */
353 for (i = 0; i < egNR; i++)
357 bEInd_[egCOULSR] = true;
358 bEInd_[egLJSR ] = true;
362 bEInd_[egLJSR] = false;
363 bEInd_[egBHAMSR] = true;
367 bEInd_[egLJ14] = true;
368 bEInd_[egCOUL14] = true;
371 for (i = 0; (i < egNR); i++)
378 n = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
387 for (k = 0; (k < nEc_); k++)
389 snew(gnm[k], STRLEN);
391 for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
393 ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
394 for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
396 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
397 for (k = kk = 0; (k < egNR); k++)
401 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
402 *(groups->groupNames[ni]), *(groups->groupNames[nj]));
406 igrp_[n] = get_ebin_space(ebin_, nEc_,
411 for (k = 0; (k < nEc_); k++)
419 gmx_incons("Number of energy terms wrong");
423 nTC_ = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
424 nNHC_ = ir->opts.nhchainlength; /* shorthand for number of NH chains */
427 nTCP_ = 1; /* assume only one possible coupling system for barostat
434 if (etc_ == etcNOSEHOOVER)
438 mde_n_ = 2*nNHC_*nTC_;
446 mdeb_n_ = 2*nNHC_*nTCP_;
455 snew(tmp_r_, mde_n_);
456 // TODO redo the group name memory management to make it more clear
458 snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
460 for (i = 0; (i < nTC_); i++)
462 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
463 sprintf(buf, "T-%s", *(groups->groupNames[ni]));
464 grpnms[i] = gmx_strdup(buf);
466 itemp_ = get_ebin_space(ebin_, nTC_, grpnms,
468 for (i = 0; i < nTC_; i++)
474 if (etc_ == etcNOSEHOOVER)
480 for (i = 0; (i < nTC_); i++)
482 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
483 bufi = *(groups->groupNames[ni]);
484 for (j = 0; (j < nNHC_); j++)
486 sprintf(buf, "Xi-%d-%s", j, bufi);
487 grpnms[2*(i*nNHC_+j)] = gmx_strdup(buf);
488 sprintf(buf, "vXi-%d-%s", j, bufi);
489 grpnms[2*(i*nNHC_+j)+1] = gmx_strdup(buf);
492 itc_ = get_ebin_space(ebin_, mde_n_,
493 grpnms, unit_invtime);
497 for (i = 0; (i < nTCP_); i++)
499 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
500 for (j = 0; (j < nNHC_); j++)
502 sprintf(buf, "Xi-%d-%s", j, bufi);
503 grpnms[2*(i*nNHC_+j)] = gmx_strdup(buf);
504 sprintf(buf, "vXi-%d-%s", j, bufi);
505 grpnms[2*(i*nNHC_+j)+1] = gmx_strdup(buf);
508 itcb_ = get_ebin_space(ebin_, mdeb_n_,
509 grpnms, unit_invtime);
515 for (i = 0; (i < nTC_); i++)
517 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
518 bufi = *(groups->groupNames[ni]);
519 sprintf(buf, "Xi-%s", bufi);
520 grpnms[2*i] = gmx_strdup(buf);
521 sprintf(buf, "vXi-%s", bufi);
522 grpnms[2*i+1] = gmx_strdup(buf);
524 itc_ = get_ebin_space(ebin_, mde_n_,
525 grpnms, unit_invtime);
530 else if (etc_ == etcBERENDSEN || etc_ == etcYES ||
533 for (i = 0; (i < nTC_); i++)
535 ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
536 sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
537 grpnms[i] = gmx_strdup(buf);
539 itc_ = get_ebin_space(ebin_, mde_n_, grpnms, "");
543 for (i = 0; i < allocated; i++)
549 nU_ = groups->groups[SimulationAtomGroupType::Acceleration].size();
554 for (i = 0; (i < nU_); i++)
556 ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
557 sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
558 grpnms[3*i+XX] = gmx_strdup(buf);
559 sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
560 grpnms[3*i+YY] = gmx_strdup(buf);
561 sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
562 grpnms[3*i+ZZ] = gmx_strdup(buf);
564 iu_ = get_ebin_space(ebin_, 3*nU_, grpnms, unit_vel);
565 for (i = 0; i < 3*nU_; i++)
574 do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
577 /* check whether we're going to write dh histograms */
579 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
581 /* Currently dh histograms are only written with dynamics */
582 if (EI_DYNAMICS(ir->eI))
586 mde_delta_h_coll_init(dhc_, ir);
589 snew(dE_, ir->fepvals->n_lambda);
594 snew(dE_, ir->fepvals->n_lambda);
599 snew(temperatures_, ir->fepvals->n_lambda);
600 numTemperatures_ = ir->fepvals->n_lambda;
601 for (i = 0; i < ir->fepvals->n_lambda; i++)
603 temperatures_[i] = ir->simtempvals->temperatures[i];
608 numTemperatures_ = 0;
613 EnergyOutput::~EnergyOutput()
619 done_mde_delta_h_coll(dhc_);
621 if (numTemperatures_ > 0)
623 sfree(temperatures_);
629 /*! \brief Print a lambda vector to a string
631 * \param[in] fep The inputrec's FEP input data
632 * \param[in] i The index of the lambda vector
633 * \param[in] get_native_lambda Whether to print the native lambda
634 * \param[in] get_names Whether to print the names rather than the values
635 * \param[in,out] str The pre-allocated string buffer to print to.
637 static void print_lambda_vector(t_lambda *fep, int i,
638 bool get_native_lambda, bool get_names,
644 for (j = 0; j < efptNR; j++)
646 if (fep->separate_dvdl[j])
651 str[0] = 0; /* reset the string */
654 str += sprintf(str, "("); /* set the opening parenthesis*/
656 for (j = 0; j < efptNR; j++)
658 if (fep->separate_dvdl[j])
662 if (get_native_lambda && fep->init_lambda >= 0)
664 str += sprintf(str, "%.4f", fep->init_lambda);
668 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
673 str += sprintf(str, "%s", efpt_singular_names[j]);
675 /* print comma for the next item */
678 str += sprintf(str, ", ");
685 /* and add the closing parenthesis */
690 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
691 const gmx_output_env_t *oenv)
694 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
695 *lambdastate = "\\lambda state";
696 int i, nsets, nsets_de, nsetsbegin;
697 int n_lambda_terms = 0;
698 t_lambda *fep = ir->fepvals; /* for simplicity */
699 t_expanded *expand = ir->expandedvals;
700 char lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
705 bool write_pV = false;
707 /* count the number of different lambda terms */
708 for (i = 0; i < efptNR; i++)
710 if (fep->separate_dvdl[i])
716 std::string title, label_x, label_y;
717 if (fep->n_lambda == 0)
719 title = gmx::formatString("%s", dhdl);
720 label_x = gmx::formatString("Time (ps)");
721 label_y = gmx::formatString("%s (%s %s)",
722 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
726 title = gmx::formatString("%s and %s", dhdl, deltag);
727 label_x = gmx::formatString("Time (ps)");
728 label_y = gmx::formatString("%s and %s (%s %s)",
729 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
731 fp = gmx_fio_fopen(filename, "w+");
732 xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
737 buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
739 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
741 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
743 /* compatibility output */
744 buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
748 print_lambda_vector(fep, fep->init_fep_state, true, false,
750 print_lambda_vector(fep, fep->init_fep_state, true, true,
752 buf += gmx::formatString("%s %d: %s = %s",
753 lambdastate, fep->init_fep_state,
754 lambda_name_str, lambda_vec_str);
757 xvgr_subtitle(fp, buf.c_str(), oenv);
761 if (fep->dhdl_derivatives == edhdlderivativesYES)
763 nsets_dhdl = n_lambda_terms;
765 /* count the number of delta_g states */
766 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
768 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
770 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
772 nsets += 1; /*add fep state for expanded ensemble */
775 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
777 nsets += 1; /* add energy to the dhdl as well */
781 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
783 nsetsextend += 1; /* for PV term, other terms possible if required for
784 the reduced potential (only needed with foreign
785 lambda, and only output when init_lambda is not
786 set in order to maintain compatibility of the
790 std::vector<std::string> setname(nsetsextend);
792 if (expand->elmcmove > elmcmoveNO)
794 /* state for the fep_vals, if we have alchemical sampling */
795 setname[s++] = "Thermodynamic state";
798 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
801 switch (fep->edHdLPrintEnergy)
803 case edHdLPrintEnergyPOTENTIAL:
804 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
806 case edHdLPrintEnergyTOTAL:
807 case edHdLPrintEnergyYES:
809 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
811 setname[s++] = energy;
814 if (fep->dhdl_derivatives == edhdlderivativesYES)
816 for (i = 0; i < efptNR; i++)
818 if (fep->separate_dvdl[i])
820 std::string derivative;
821 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
823 /* compatibility output */
824 derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
828 double lam = fep->init_lambda;
829 if (fep->init_lambda < 0)
831 lam = fep->all_lambda[i][fep->init_fep_state];
833 derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
836 setname[s++] = derivative;
841 if (fep->n_lambda > 0)
843 /* g_bar has to determine the lambda values used in this simulation
844 * from this xvg legend.
847 if (expand->elmcmove > elmcmoveNO)
849 nsetsbegin = 1; /* for including the expanded ensemble */
856 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
860 nsetsbegin += nsets_dhdl;
862 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
864 print_lambda_vector(fep, i, false, false, lambda_vec_str);
866 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
868 /* for compatible dhdl.xvg files */
869 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
873 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
878 /* print the temperature for this state if doing simulated annealing */
879 buf += gmx::formatString("T = %g (%s)",
880 ir->simtempvals->temperatures[s-(nsetsbegin)],
887 setname[s++] = gmx::formatString("pV (%s)", unit_energy);
890 xvgrLegend(fp, setname, oenv);
899 void EnergyOutput::addDataAtEnergyStep(bool bDoDHDL,
903 const gmx_enerdata_t *enerd,
904 const t_state *state,
906 const t_expanded *expand,
912 const gmx_ekindata_t *ekind,
914 const gmx::Constraints *constr)
916 int j, k, kk, n, gid;
917 real crmsd[2], tmp6[6];
918 real bs[tricl_boxs_nm.size()], vol, dens, pv, enthalpy;
920 double store_dhdl[efptNR];
921 real store_energy = 0;
924 /* Do NOT use the box in the state variable, but the separate box provided
925 * as an argument. This is because we sometimes need to write the box from
926 * the last timestep to match the trajectory frames.
928 add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
931 crmsd[0] = constr->rmsd();
932 add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
945 nboxs = tricl_boxs_nm.size();
952 nboxs = boxs_nm.size();
954 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
955 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
956 add_ebin(ebin_, ib_, nboxs, bs, bSum);
957 add_ebin(ebin_, ivol_, 1, &vol, bSum);
958 add_ebin(ebin_, idens_, 1, &dens, bSum);
962 /* This is pV (in kJ/mol). The pressure is the reference pressure,
963 not the instantaneous pressure */
964 pv = vol*ref_p_/PRESFAC;
966 add_ebin(ebin_, ipv_, 1, &pv, bSum);
967 enthalpy = pv + enerd->term[F_ETOT];
968 add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
973 add_ebin(ebin_, isvir_, 9, svir[0], bSum);
974 add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
978 add_ebin(ebin_, ivir_, 9, vir[0], bSum);
979 add_ebin(ebin_, ipres_, 9, pres[0], bSum);
980 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
981 add_ebin(ebin_, isurft_, 1, &tmp, bSum);
983 if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
985 tmp6[0] = state->boxv[XX][XX];
986 tmp6[1] = state->boxv[YY][YY];
987 tmp6[2] = state->boxv[ZZ][ZZ];
988 tmp6[3] = state->boxv[YY][XX];
989 tmp6[4] = state->boxv[ZZ][XX];
990 tmp6[5] = state->boxv[ZZ][YY];
991 add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
995 add_ebin(ebin_, imu_, 3, mu_tot, bSum);
997 if (ekind && ekind->cosacc.cos_accel != 0)
999 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1000 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1001 add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
1002 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1003 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1004 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1005 add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
1010 for (int i = 0; (i < nEg_); i++)
1012 for (j = i; (j < nEg_); j++)
1014 gid = GID(i, j, nEg_);
1015 for (k = kk = 0; (k < egNR); k++)
1019 eee[kk++] = enerd->grpp.ener[k][gid];
1022 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
1030 for (int i = 0; (i < nTC_); i++)
1032 tmp_r_[i] = ekind->tcstat[i].T;
1034 add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
1036 if (etc_ == etcNOSEHOOVER)
1038 /* whether to print Nose-Hoover chains: */
1039 if (bPrintNHChains_)
1043 for (int i = 0; (i < nTC_); i++)
1045 for (j = 0; j < nNHC_; j++)
1048 tmp_r_[2*k] = state->nosehoover_xi[k];
1049 tmp_r_[2*k+1] = state->nosehoover_vxi[k];
1052 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1056 for (int i = 0; (i < nTCP_); i++)
1058 for (j = 0; j < nNHC_; j++)
1061 tmp_r_[2*k] = state->nhpres_xi[k];
1062 tmp_r_[2*k+1] = state->nhpres_vxi[k];
1065 add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1070 for (int i = 0; (i < nTC_); i++)
1072 tmp_r_[2*i] = state->nosehoover_xi[i];
1073 tmp_r_[2*i+1] = state->nosehoover_vxi[i];
1075 add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1079 else if (etc_ == etcBERENDSEN || etc_ == etcYES ||
1080 etc_ == etcVRESCALE)
1082 for (int i = 0; (i < nTC_); i++)
1084 tmp_r_[i] = ekind->tcstat[i].lambda;
1086 add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1090 if (ekind && nU_ > 1)
1092 for (int i = 0; (i < nU_); i++)
1094 copy_rvec(ekind->grpstat[i].u, tmp_v_[i]);
1096 add_ebin(ebin_, iu_, 3*nU_, tmp_v_[0], bSum);
1099 ebin_increase_count(1, ebin_, bSum);
1101 // BAR + thermodynamic integration values
1102 if ((fp_dhdl_ || dhc_) && bDoDHDL)
1104 for (gmx::index i = 0; i < static_cast<gmx::index>(enerd->enerpart_lambda.size()) - 1; i++)
1106 /* zero for simulated tempering */
1107 dE_[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1108 if (numTemperatures_ > 0)
1110 GMX_RELEASE_ASSERT(numTemperatures_ > state->fep_state, "Number of lambdas in state is bigger then in input record");
1111 GMX_RELEASE_ASSERT(numTemperatures_ >= static_cast<gmx::index>(enerd->enerpart_lambda.size()) - 1, "Number of lambdas in energy data is bigger then in input record");
1112 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1113 /* is this even useful to have at all? */
1114 dE_[i] += (temperatures_[i]/
1115 temperatures_[state->fep_state]-1.0)*
1116 enerd->term[F_EKIN];
1122 fprintf(fp_dhdl_, "%.4f", time);
1123 /* the current free energy state */
1125 /* print the current state if we are doing expanded ensemble */
1126 if (expand->elmcmove > elmcmoveNO)
1128 fprintf(fp_dhdl_, " %4d", state->fep_state);
1130 /* total energy (for if the temperature changes */
1132 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1134 switch (fep->edHdLPrintEnergy)
1136 case edHdLPrintEnergyPOTENTIAL:
1137 store_energy = enerd->term[F_EPOT];
1139 case edHdLPrintEnergyTOTAL:
1140 case edHdLPrintEnergyYES:
1142 store_energy = enerd->term[F_ETOT];
1144 fprintf(fp_dhdl_, " %#.8g", store_energy);
1147 if (fep->dhdl_derivatives == edhdlderivativesYES)
1149 for (int i = 0; i < efptNR; i++)
1151 if (fep->separate_dvdl[i])
1153 /* assumes F_DVDL is first */
1154 fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL+i]);
1158 for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1160 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1165 !enerd->enerpart_lambda.empty() &&
1166 (fep->init_lambda < 0))
1168 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1169 there are alternate state
1170 lambda and we're not in
1171 compatibility mode */
1173 fprintf(fp_dhdl_, "\n");
1174 /* and the binary free energy output */
1176 if (dhc_ && bDoDHDL)
1179 for (int i = 0; i < efptNR; i++)
1181 if (fep->separate_dvdl[i])
1183 /* assumes F_DVDL is first */
1184 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1188 store_energy = enerd->term[F_ETOT];
1189 /* store_dh is dE */
1190 mde_delta_h_coll_add_dh(dhc_,
1191 static_cast<double>(state->fep_state),
1195 dE_ + fep->lambda_start_n,
1202 void EnergyOutput::recordNonEnergyStep()
1204 ebin_increase_count(1, ebin_, false);
1207 void EnergyOutput::printHeader(FILE *log, int64_t steps, double time)
1211 fprintf(log, " %12s %12s\n"
1213 "Step", "Time", gmx_step_str(steps, buf), time);
1216 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1218 int64_t step, double time,
1226 fr.nsteps = ebin_->nsteps;
1228 fr.nsum = ebin_->nsum;
1229 fr.nre = (bEne) ? ebin_->nener : 0;
1231 int ndisre = bDR ? fcd->disres.npair : 0;
1232 /* these are for the old-style blocks (1 subblock, only reals), because
1233 there can be only one per ID for these */
1237 /* Optional additional old-style (real-only) blocks. */
1238 for (int i = 0; i < enxNR; i++)
1243 if (bOR && fcd->orires.nr > 0)
1245 diagonalize_orires_tensors(&(fcd->orires));
1246 nr[enxOR] = fcd->orires.nr;
1247 block[enxOR] = fcd->orires.otav;
1249 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1251 block[enxORI] = fcd->orires.oinsl;
1252 id[enxORI] = enxORI;
1253 nr[enxORT] = fcd->orires.nex*12;
1254 block[enxORT] = fcd->orires.eig;
1255 id[enxORT] = enxORT;
1258 /* whether we are going to write anything out: */
1259 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1261 /* the old-style blocks go first */
1263 for (int i = 0; i < enxNR; i++)
1270 add_blocks_enxframe(&fr, fr.nblock);
1271 for (int b = 0; b < fr.nblock; b++)
1273 add_subblocks_enxblock(&(fr.block[b]), 1);
1274 fr.block[b].id = id[b];
1275 fr.block[b].sub[0].nr = nr[b];
1277 fr.block[b].sub[0].type = xdr_datatype_float;
1278 fr.block[b].sub[0].fval = block[b];
1280 fr.block[b].sub[0].type = xdr_datatype_double;
1281 fr.block[b].sub[0].dval = block[b];
1285 /* check for disre block & fill it. */
1290 add_blocks_enxframe(&fr, fr.nblock);
1292 add_subblocks_enxblock(&(fr.block[db]), 2);
1293 fr.block[db].id = enxDISRE;
1294 fr.block[db].sub[0].nr = ndisre;
1295 fr.block[db].sub[1].nr = ndisre;
1297 fr.block[db].sub[0].type = xdr_datatype_float;
1298 fr.block[db].sub[1].type = xdr_datatype_float;
1299 fr.block[db].sub[0].fval = fcd->disres.rt;
1300 fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1302 fr.block[db].sub[0].type = xdr_datatype_double;
1303 fr.block[db].sub[1].type = xdr_datatype_double;
1304 fr.block[db].sub[0].dval = fcd->disres.rt;
1305 fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1308 /* here we can put new-style blocks */
1310 /* Free energy perturbation blocks */
1313 mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1316 /* we can now free & reset the data in the blocks */
1319 mde_delta_h_coll_reset(dhc_);
1322 /* AWH bias blocks. */
1323 if (awh != nullptr) // TODO: add boolean flag.
1325 awh->writeToEnergyFrame(step, &fr);
1328 /* do the actual I/O */
1329 do_enx(fp_ene, &fr);
1332 /* We have stored the sums, so reset the sum history */
1333 reset_ebin_sums(ebin_);
1339 if (bOR && fcd->orires.nr > 0)
1341 print_orires_log(log, &(fcd->orires));
1344 fprintf(log, " Energies (%s)\n", unit_energy);
1345 pr_ebin(log, ebin_, ie_, f_nre_+nCrmsd_, 5, eprNORMAL, true);
1350 void EnergyOutput::printAnnealingTemperatures(FILE *log, SimulationGroups *groups, t_grpopts *opts)
1356 for (int i = 0; i < opts->ngtc; i++)
1358 if (opts->annealing[i] != eannNO)
1360 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1361 *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1370 void EnergyOutput::printAverages(FILE *log, const SimulationGroups *groups)
1372 if (ebin_->nsum_sim <= 0)
1376 fprintf(log, "Not enough data recorded to report energy averages\n");
1383 char buf1[22], buf2[22];
1385 fprintf(log, "\t<====== ############### ==>\n");
1386 fprintf(log, "\t<==== A V E R A G E S ====>\n");
1387 fprintf(log, "\t<== ############### ======>\n\n");
1389 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1390 gmx_step_str(ebin_->nsteps_sim, buf1),
1391 gmx_step_str(ebin_->nsum_sim, buf2));
1394 fprintf(log, " Energies (%s)\n", unit_energy);
1395 pr_ebin(log, ebin_, ie_, f_nre_+nCrmsd_, 5, eprAVER, true);
1400 pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5,
1406 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1407 pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1409 fprintf(log, " Force Virial (%s)\n", unit_energy);
1410 pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1415 fprintf(log, " Total Virial (%s)\n", unit_energy);
1416 pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1418 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1419 pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1424 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1425 pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1431 int padding = 8 - strlen(unit_energy);
1432 fprintf(log, "%*sEpot (%s) ", padding, "", unit_energy);
1433 for (int i = 0; (i < egNR); i++)
1437 fprintf(log, "%12s ", egrp_nm[i]);
1443 for (int i = 0; (i < nEg_); i++)
1445 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1446 for (int j = i; (j < nEg_); j++)
1448 int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1449 int padding = 14 - (strlen(*(groups->groupNames[ni])) +
1450 strlen(*(groups->groupNames[nj])));
1451 fprintf(log, "%*s%s-%s", padding, "",
1452 *(groups->groupNames[ni]),
1453 *(groups->groupNames[nj]));
1454 pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER,
1463 pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1468 fprintf(log, "%15s %12s %12s %12s\n",
1469 "Group", "Ux", "Uy", "Uz");
1470 for (int i = 0; (i < nU_); i++)
1472 int ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
1473 fprintf(log, "%15s", *groups->groupNames[ni]);
1474 pr_ebin(log, ebin_, iu_+3*i, 3, 3, eprAVER, false);
1481 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1483 const t_ebin * const ebin = ebin_;
1485 enerhist->nsteps = ebin->nsteps;
1486 enerhist->nsum = ebin->nsum;
1487 enerhist->nsteps_sim = ebin->nsteps_sim;
1488 enerhist->nsum_sim = ebin->nsum_sim;
1492 /* This will only actually resize the first time */
1493 enerhist->ener_ave.resize(ebin->nener);
1494 enerhist->ener_sum.resize(ebin->nener);
1496 for (int i = 0; i < ebin->nener; i++)
1498 enerhist->ener_ave[i] = ebin->e[i].eav;
1499 enerhist->ener_sum[i] = ebin->e[i].esum;
1503 if (ebin->nsum_sim > 0)
1505 /* This will only actually resize the first time */
1506 enerhist->ener_sum_sim.resize(ebin->nener);
1508 for (int i = 0; i < ebin->nener; i++)
1510 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1515 mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1519 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1521 unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1523 if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size()) ||
1524 (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1526 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1527 nener, enerhist.ener_sum.size(), enerhist.ener_sum_sim.size());
1530 ebin_->nsteps = enerhist.nsteps;
1531 ebin_->nsum = enerhist.nsum;
1532 ebin_->nsteps_sim = enerhist.nsteps_sim;
1533 ebin_->nsum_sim = enerhist.nsum_sim;
1535 for (int i = 0; i < ebin_->nener; i++)
1538 (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1540 (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1541 ebin_->e_sim[i].esum =
1542 (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1546 mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1550 int EnergyOutput::numEnergyTerms() const
1552 return ebin_->nener;