Make non bonded energy terms enum class
[alexxy/gromacs.git] / src / gromacs / mdlib / energyoutput.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \internal \file
39  * \brief Defines code that writes energy-like quantities.
40  *
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>
44  *
45  * \ingroup module_mdlib
46  */
47 #include "gmxpre.h"
48
49 #include "energyoutput.h"
50
51 #include <cfloat>
52 #include <cstdlib>
53 #include <cstring>
54
55 #include <array>
56 #include <string>
57
58 #include "gromacs/applied_forces/awh/awh.h"
59 #include "gromacs/fileio/enxio.h"
60 #include "gromacs/fileio/gmxfio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/gmxlib/network.h"
63 #include "gromacs/listed_forces/disre.h"
64 #include "gromacs/listed_forces/orires.h"
65 #include "gromacs/math/functions.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/constr.h"
69 #include "gromacs/mdlib/ebin.h"
70 #include "gromacs/mdlib/mdebin_bar.h"
71 #include "gromacs/mdrunutility/handlerestart.h"
72 #include "gromacs/mdtypes/energyhistory.h"
73 #include "gromacs/mdtypes/fcdata.h"
74 #include "gromacs/mdtypes/group.h"
75 #include "gromacs/mdtypes/inputrec.h"
76 #include "gromacs/mdtypes/md_enums.h"
77 #include "gromacs/mdtypes/state.h"
78 #include "gromacs/pbcutil/pbc.h"
79 #include "gromacs/pulling/pull.h"
80 #include "gromacs/topology/mtop_util.h"
81 #include "gromacs/trajectory/energyframe.h"
82 #include "gromacs/utility/arraysize.h"
83 #include "gromacs/utility/enumerationhelpers.h"
84 #include "gromacs/utility/fatalerror.h"
85 #include "gromacs/utility/mdmodulenotification.h"
86 #include "gromacs/utility/smalloc.h"
87 #include "gromacs/utility/stringutil.h"
88
89 #include "energydrifttracker.h"
90
91 //! Labels for energy file quantities
92 //! \{
93 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
94
95 static std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
96
97 static std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
98                                                     "Box-YX", "Box-ZX", "Box-ZY" };
99
100 static const char* vol_nm[] = { "Volume" };
101
102 static const char* dens_nm[] = { "Density" };
103
104 static const char* pv_nm[] = { "pV" };
105
106 static const char* enthalpy_nm[] = { "Enthalpy" };
107
108 static std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
109                                                 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
110
111 const char* enumValueToString(NonBondedEnergyTerms enumValue)
112 {
113     static constexpr gmx::EnumerationArray<NonBondedEnergyTerms, const char*> nonBondedEnergyTermTypeNames = {
114         "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14"
115     };
116     return nonBondedEnergyTermTypeNames[enumValue];
117 }
118
119 //! \}
120
121 namespace gmx
122 {
123
124 /*! \brief Energy output class
125  *
126  * This is the collection of energy averages collected during mdrun, and to
127  * be written out to the .edr file.
128  *
129  * \todo Use more std containers.
130  * \todo Remove GMX_CONSTRAINTVIR
131  * \todo Write free-energy output also to energy file (after adding more tests)
132  */
133 EnergyOutput::EnergyOutput(ener_file*               fp_ene,
134                            const gmx_mtop_t&        mtop,
135                            const t_inputrec&        inputrec,
136                            const pull_t*            pull_work,
137                            FILE*                    fp_dhdl,
138                            bool                     isRerun,
139                            const StartingBehavior   startingBehavior,
140                            const bool               simulationsShareState,
141                            const MdModulesNotifier& mdModulesNotifier)
142 {
143     const char*        ener_nm[F_NRE];
144     static const char* vir_nm[]   = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
145                                     "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
146     static const char* sv_nm[]    = { "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
147                                    "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
148                                    "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ" };
149     static const char* fv_nm[]    = { "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
150                                    "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
151                                    "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ" };
152     static const char* pres_nm[]  = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
153                                      "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
154     static const char* surft_nm[] = { "#Surf*SurfTen" };
155     static const char* mu_nm[]    = { "Mu-X", "Mu-Y", "Mu-Z" };
156     static const char* vcos_nm[]  = { "2CosZ*Vel-X" };
157     static const char* visc_nm[]  = { "1/Viscosity" };
158     static const char* baro_nm[]  = { "Barostat" };
159
160     const SimulationGroups* groups;
161     char**                  gnm;
162     char                    buf[256];
163     const char*             bufi;
164     int                     i, j, ni, nj, n, ncon, nset;
165     bool                    bBHAM, b14;
166
167     if (EI_DYNAMICS(inputrec.eI))
168     {
169         delta_t_ = inputrec.delta_t;
170     }
171     else
172     {
173         delta_t_ = 0;
174     }
175
176     groups = &mtop.groups;
177
178     bBHAM = (mtop.ffparams.numTypes() > 0) && (mtop.ffparams.functype[0] == F_BHAM);
179     b14   = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
180
181     ncon         = gmx_mtop_ftype_count(mtop, F_CONSTR);
182     nset         = gmx_mtop_ftype_count(mtop, F_SETTLE);
183     bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
184     bConstrVir_  = false;
185     nCrmsd_      = 0;
186     if (bConstr)
187     {
188         if (ncon > 0 && inputrec.eConstrAlg == ConstraintAlgorithm::Lincs)
189         {
190             nCrmsd_ = 1;
191         }
192         bConstrVir_ = (getenv("GMX_CONSTRAINTVIR") != nullptr);
193     }
194     else
195     {
196         nCrmsd_ = 0;
197     }
198
199     /* Energy monitoring */
200     for (auto& term : bEInd_)
201     {
202         term = false;
203     }
204
205     // Setting true only to those energy terms, that have active interactions and
206     // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
207     for (i = 0; i < F_NRE; i++)
208     {
209         bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
210                     && ((interaction_function[i].flags & IF_VSITE) == 0);
211     }
212
213     if (!isRerun)
214     {
215         bEner_[F_EKIN] = EI_DYNAMICS(inputrec.eI);
216         bEner_[F_ETOT] = EI_DYNAMICS(inputrec.eI);
217         bEner_[F_TEMP] = EI_DYNAMICS(inputrec.eI);
218
219         bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(&inputrec);
220         bEner_[F_PDISPCORR]  = (inputrec.eDispCorr != DispersionCorrectionType::No);
221         bEner_[F_PRES]       = true;
222     }
223
224     bEner_[F_LJ]   = !bBHAM;
225     bEner_[F_BHAM] = bBHAM;
226     bEner_[F_EQM]  = inputrec.bQMMM;
227     bEner_[F_RF_EXCL] = (EEL_RF(inputrec.coulombtype) && inputrec.cutoff_scheme == CutoffScheme::Group);
228     bEner_[F_COUL_RECIP]   = EEL_FULL(inputrec.coulombtype);
229     bEner_[F_LJ_RECIP]     = EVDW_PME(inputrec.vdwtype);
230     bEner_[F_LJ14]         = b14;
231     bEner_[F_COUL14]       = b14;
232     bEner_[F_LJC14_Q]      = false;
233     bEner_[F_LJC_PAIRS_NB] = false;
234
235
236     bEner_[F_DVDL_COUL] = (inputrec.efep != FreeEnergyPerturbationType::No)
237                           && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Coul];
238     bEner_[F_DVDL_VDW] = (inputrec.efep != FreeEnergyPerturbationType::No)
239                          && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Vdw];
240     bEner_[F_DVDL_BONDED] = (inputrec.efep != FreeEnergyPerturbationType::No)
241                             && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Bonded];
242     bEner_[F_DVDL_RESTRAINT] =
243             (inputrec.efep != FreeEnergyPerturbationType::No)
244             && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Restraint];
245     bEner_[F_DKDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
246                      && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Mass];
247     bEner_[F_DVDL] = (inputrec.efep != FreeEnergyPerturbationType::No)
248                      && inputrec.fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep];
249
250     bEner_[F_CONSTR]   = false;
251     bEner_[F_CONSTRNC] = false;
252     bEner_[F_SETTLE]   = false;
253
254     bEner_[F_COUL_SR] = true;
255     bEner_[F_EPOT]    = true;
256
257     bEner_[F_DISPCORR]   = (inputrec.eDispCorr != DispersionCorrectionType::No);
258     bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
259     bEner_[F_ORIRESDEV]  = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
260     bEner_[F_COM_PULL]   = ((inputrec.bPull && pull_have_potential(*pull_work)) || inputrec.bRot);
261
262     MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
263     mdModulesNotifier.simulationSetupNotifications_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
264
265     bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
266
267
268     // Counting the energy terms that will be printed and saving their names
269     f_nre_ = 0;
270     for (i = 0; i < F_NRE; i++)
271     {
272         if (bEner_[i])
273         {
274             ener_nm[f_nre_] = interaction_function[i].longname;
275             f_nre_++;
276         }
277     }
278
279     epc_       = isRerun ? PressureCoupling::No : inputrec.epc;
280     bDiagPres_ = !TRICLINIC(inputrec.ref_p) && !isRerun;
281     ref_p_     = (inputrec.ref_p[XX][XX] + inputrec.ref_p[YY][YY] + inputrec.ref_p[ZZ][ZZ]) / DIM;
282     bTricl_    = TRICLINIC(inputrec.compress) || TRICLINIC(inputrec.deform);
283     bDynBox_   = inputrecDynamicBox(&inputrec);
284     etc_       = isRerun ? TemperatureCoupling::No : inputrec.etc;
285     bNHC_trotter_   = inputrecNvtTrotter(&inputrec) && !isRerun;
286     bPrintNHChains_ = inputrec.bPrintNHChains && !isRerun;
287     bMTTK_          = (inputrecNptTrotter(&inputrec) || inputrecNphTrotter(&inputrec)) && !isRerun;
288     bMu_            = inputrecNeedMutot(&inputrec);
289     bPres_          = !isRerun;
290
291     ebin_ = mk_ebin();
292     /* Pass NULL for unit to let get_ebin_space determine the units
293      * for interaction_function[i].longname
294      */
295     ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
296     if (nCrmsd_)
297     {
298         /* This should be called directly after the call for ie_,
299          * such that iconrmsd_ follows directly in the list.
300          */
301         iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
302     }
303     if (bDynBox_)
304     {
305         ib_    = get_ebin_space(ebin_,
306                              bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
307                              bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
308                              unit_length);
309         ivol_  = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
310         idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
311         if (bDiagPres_)
312         {
313             ipv_       = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
314             ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
315         }
316     }
317     if (bConstrVir_)
318     {
319         isvir_ = get_ebin_space(ebin_, asize(sv_nm), sv_nm, unit_energy);
320         ifvir_ = get_ebin_space(ebin_, asize(fv_nm), fv_nm, unit_energy);
321     }
322     if (bPres_)
323     {
324         ivir_   = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
325         ipres_  = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
326         isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
327     }
328     if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
329     {
330         ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
331     }
332     if (bMu_)
333     {
334         imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
335     }
336     if (inputrec.cos_accel != 0)
337     {
338         ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
339         ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
340     }
341
342     /* Energy monitoring */
343     for (auto& term : bEInd_)
344     {
345         term = false;
346     }
347     bEInd_[NonBondedEnergyTerms::CoulombSR] = true;
348     bEInd_[NonBondedEnergyTerms::LJSR]      = true;
349
350     if (bBHAM)
351     {
352         bEInd_[NonBondedEnergyTerms::LJSR]         = false;
353         bEInd_[NonBondedEnergyTerms::BuckinghamSR] = true;
354     }
355     if (b14)
356     {
357         bEInd_[NonBondedEnergyTerms::LJ14]      = true;
358         bEInd_[NonBondedEnergyTerms::Coulomb14] = true;
359     }
360     nEc_ = 0;
361     for (auto term : bEInd_)
362     {
363         if (term)
364         {
365             nEc_++;
366         }
367     }
368     n    = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
369     nEg_ = n;
370     nE_  = (n * (n + 1)) / 2;
371
372     snew(igrp_, nE_);
373     if (nE_ > 1)
374     {
375         n = 0;
376         snew(gnm, nEc_);
377         for (int k = 0; (k < nEc_); k++)
378         {
379             snew(gnm[k], STRLEN);
380         }
381         for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
382         {
383             ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
384             for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
385             {
386                 nj    = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
387                 int k = 0;
388                 for (auto key : keysOf(bEInd_))
389                 {
390                     if (bEInd_[key])
391                     {
392                         sprintf(gnm[k],
393                                 "%s:%s-%s",
394                                 enumValueToString(key),
395                                 *(groups->groupNames[ni]),
396                                 *(groups->groupNames[nj]));
397                         k++;
398                     }
399                 }
400                 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
401                 n++;
402             }
403         }
404         for (int k = 0; (k < nEc_); k++)
405         {
406             sfree(gnm[k]);
407         }
408         sfree(gnm);
409
410         if (n != nE_)
411         {
412             gmx_incons("Number of energy terms wrong");
413         }
414     }
415
416     nTC_  = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
417     nNHC_ = inputrec.opts.nhchainlength; /* shorthand for number of NH chains */
418     if (bMTTK_)
419     {
420         nTCP_ = 1; /* assume only one possible coupling system for barostat
421                            for now */
422     }
423     else
424     {
425         nTCP_ = 0;
426     }
427     if (etc_ == TemperatureCoupling::NoseHoover)
428     {
429         if (bNHC_trotter_)
430         {
431             mde_n_ = 2 * nNHC_ * nTC_;
432         }
433         else
434         {
435             mde_n_ = 2 * nTC_;
436         }
437         if (epc_ == PressureCoupling::Mttk)
438         {
439             mdeb_n_ = 2 * nNHC_ * nTCP_;
440         }
441     }
442     else
443     {
444         mde_n_  = nTC_;
445         mdeb_n_ = 0;
446     }
447
448     snew(tmp_r_, mde_n_);
449     // TODO redo the group name memory management to make it more clear
450     char** grpnms;
451     snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
452
453     for (i = 0; (i < nTC_); i++)
454     {
455         ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
456         sprintf(buf, "T-%s", *(groups->groupNames[ni]));
457         grpnms[i] = gmx_strdup(buf);
458     }
459     itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
460     for (i = 0; i < nTC_; i++)
461     {
462         sfree(grpnms[i]);
463     }
464
465     int allocated = 0;
466     if (etc_ == TemperatureCoupling::NoseHoover)
467     {
468         if (bPrintNHChains_)
469         {
470             if (bNHC_trotter_)
471             {
472                 for (i = 0; (i < nTC_); i++)
473                 {
474                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
475                     bufi = *(groups->groupNames[ni]);
476                     for (j = 0; (j < nNHC_); j++)
477                     {
478                         sprintf(buf, "Xi-%d-%s", j, bufi);
479                         grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
480                         sprintf(buf, "vXi-%d-%s", j, bufi);
481                         grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
482                     }
483                 }
484                 itc_      = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
485                 allocated = mde_n_;
486                 if (bMTTK_)
487                 {
488                     for (i = 0; (i < nTCP_); i++)
489                     {
490                         bufi = baro_nm[0]; /* All barostat DOF's together for now. */
491                         for (j = 0; (j < nNHC_); j++)
492                         {
493                             sprintf(buf, "Xi-%d-%s", j, bufi);
494                             grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
495                             sprintf(buf, "vXi-%d-%s", j, bufi);
496                             grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
497                         }
498                     }
499                     itcb_     = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
500                     allocated = mdeb_n_;
501                 }
502             }
503             else
504             {
505                 for (i = 0; (i < nTC_); i++)
506                 {
507                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
508                     bufi = *(groups->groupNames[ni]);
509                     sprintf(buf, "Xi-%s", bufi);
510                     grpnms[2 * i] = gmx_strdup(buf);
511                     sprintf(buf, "vXi-%s", bufi);
512                     grpnms[2 * i + 1] = gmx_strdup(buf);
513                 }
514                 itc_      = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
515                 allocated = mde_n_;
516             }
517         }
518     }
519     else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
520              || etc_ == TemperatureCoupling::VRescale)
521     {
522         for (i = 0; (i < nTC_); i++)
523         {
524             ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
525             sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
526             grpnms[i] = gmx_strdup(buf);
527         }
528         itc_      = get_ebin_space(ebin_, mde_n_, grpnms, "");
529         allocated = mde_n_;
530     }
531
532     for (i = 0; i < allocated; i++)
533     {
534         sfree(grpnms[i]);
535     }
536     sfree(grpnms);
537
538     /* Note that fp_ene should be valid on the master rank and null otherwise */
539     if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
540     {
541         do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
542     }
543
544     /* check whether we're going to write dh histograms */
545     dhc_ = nullptr;
546     if (inputrec.fepvals->separate_dhdl_file == SeparateDhdlFile::No)
547     {
548         /* Currently dh histograms are only written with dynamics */
549         if (EI_DYNAMICS(inputrec.eI))
550         {
551             snew(dhc_, 1);
552
553             mde_delta_h_coll_init(dhc_, inputrec);
554         }
555         fp_dhdl_ = nullptr;
556         snew(dE_, inputrec.fepvals->n_lambda);
557     }
558     else
559     {
560         fp_dhdl_ = fp_dhdl;
561         snew(dE_, inputrec.fepvals->n_lambda);
562     }
563     if (inputrec.bSimTemp)
564     {
565         int i;
566         snew(temperatures_, inputrec.fepvals->n_lambda);
567         numTemperatures_ = inputrec.fepvals->n_lambda;
568         for (i = 0; i < inputrec.fepvals->n_lambda; i++)
569         {
570             temperatures_[i] = inputrec.simtempvals->temperatures[i];
571         }
572     }
573     else
574     {
575         numTemperatures_ = 0;
576     }
577
578     if (EI_MD(inputrec.eI) && !simulationsShareState)
579     {
580         conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop.natoms);
581     }
582 }
583
584 EnergyOutput::~EnergyOutput()
585 {
586     sfree(igrp_);
587     sfree(tmp_r_);
588     sfree(tmp_v_);
589     done_ebin(ebin_);
590     done_mde_delta_h_coll(dhc_);
591     sfree(dE_);
592     if (numTemperatures_ > 0)
593     {
594         sfree(temperatures_);
595     }
596 }
597
598 } // namespace gmx
599
600 /*! \brief Print a lambda vector to a string
601  *
602  * \param[in] fep                The inputrec's FEP input data
603  * \param[in] i                  The index of the lambda vector
604  * \param[in] get_native_lambda  Whether to print the native lambda
605  * \param[in] get_names          Whether to print the names rather than the values
606  * \param[in,out] str            The pre-allocated string buffer to print to.
607  */
608 static void print_lambda_vector(t_lambda* fep, int i, bool get_native_lambda, bool get_names, char* str)
609 {
610     int k    = 0;
611     int Nsep = 0;
612
613     for (auto j : keysOf(fep->separate_dvdl))
614     {
615         if (fep->separate_dvdl[j])
616         {
617             Nsep++;
618         }
619     }
620     str[0] = 0; /* reset the string */
621     if (Nsep > 1)
622     {
623         str += sprintf(str, "("); /* set the opening parenthesis*/
624     }
625     for (auto j : keysOf(fep->separate_dvdl))
626     {
627         if (fep->separate_dvdl[j])
628         {
629             if (!get_names)
630             {
631                 if (get_native_lambda && fep->init_lambda >= 0)
632                 {
633                     str += sprintf(str, "%.4f", fep->init_lambda);
634                 }
635                 else
636                 {
637                     str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
638                 }
639             }
640             else
641             {
642                 str += sprintf(str, "%s", enumValueToStringSingular(j));
643             }
644             /* print comma for the next item */
645             if (k < Nsep - 1)
646             {
647                 str += sprintf(str, ", ");
648             }
649             k++;
650         }
651     }
652     if (Nsep > 1)
653     {
654         /* and add the closing parenthesis */
655         sprintf(str, ")");
656     }
657 }
658
659 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv)
660 {
661     FILE*       fp;
662     const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
663                *lambdastate = "\\lambda state";
664     int         i, nsets, nsets_de, nsetsbegin;
665     int         n_lambda_terms = 0;
666     t_lambda*   fep            = ir->fepvals.get(); /* for simplicity */
667     t_expanded* expand         = ir->expandedvals.get();
668     char        lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
669
670     int  nsets_dhdl = 0;
671     int  s          = 0;
672     int  nsetsextend;
673     bool write_pV = false;
674
675     /* count the number of different lambda terms */
676     for (auto i : keysOf(fep->separate_dvdl))
677     {
678         if (fep->separate_dvdl[i])
679         {
680             n_lambda_terms++;
681         }
682     }
683
684     std::string title, label_x, label_y;
685     if (fep->n_lambda == 0)
686     {
687         title   = gmx::formatString("%s", dhdl);
688         label_x = gmx::formatString("Time (ps)");
689         label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
690     }
691     else
692     {
693         title   = gmx::formatString("%s and %s", dhdl, deltag);
694         label_x = gmx::formatString("Time (ps)");
695         label_y = gmx::formatString(
696                 "%s and %s (%s %s)", dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
697     }
698     fp = gmx_fio_fopen(filename, "w+");
699     xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
700
701     std::string buf;
702     if (!(ir->bSimTemp))
703     {
704         buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
705     }
706     if ((ir->efep != FreeEnergyPerturbationType::SlowGrowth)
707         && (ir->efep != FreeEnergyPerturbationType::Expanded))
708     {
709         if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
710         {
711             /* compatibility output */
712             buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
713         }
714         else
715         {
716             print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
717             print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
718             buf += gmx::formatString(
719                     "%s %d: %s = %s", lambdastate, fep->init_fep_state, lambda_name_str, lambda_vec_str);
720         }
721     }
722     xvgr_subtitle(fp, buf.c_str(), oenv);
723
724
725     nsets_dhdl = 0;
726     if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
727     {
728         nsets_dhdl = n_lambda_terms;
729     }
730     /* count the number of delta_g states */
731     nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
732
733     nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
734
735     if (fep->n_lambda > 0 && (expand->elmcmove > LambdaMoveCalculation::No))
736     {
737         nsets += 1; /*add fep state for expanded ensemble */
738     }
739
740     if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
741     {
742         nsets += 1; /* add energy to the dhdl as well */
743     }
744
745     nsetsextend = nsets;
746     if ((ir->epc != PressureCoupling::No) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
747     {
748         nsetsextend += 1; /* for PV term, other terms possible if required for
749                              the reduced potential (only needed with foreign
750                              lambda, and only output when init_lambda is not
751                              set in order to maintain compatibility of the
752                              dhdl.xvg file) */
753         write_pV = true;
754     }
755     std::vector<std::string> setname(nsetsextend);
756
757     if (expand->elmcmove > LambdaMoveCalculation::No)
758     {
759         /* state for the fep_vals, if we have alchemical sampling */
760         setname[s++] = "Thermodynamic state";
761     }
762
763     if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
764     {
765         std::string energy;
766         switch (fep->edHdLPrintEnergy)
767         {
768             case FreeEnergyPrintEnergy::Potential:
769                 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
770                 break;
771             case FreeEnergyPrintEnergy::Total:
772             case FreeEnergyPrintEnergy::Yes:
773             default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
774         }
775         setname[s++] = energy;
776     }
777
778     if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
779     {
780         for (auto i : keysOf(fep->separate_dvdl))
781         {
782             if (fep->separate_dvdl[i])
783             {
784                 std::string derivative;
785                 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
786                 {
787                     /* compatibility output */
788                     derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
789                 }
790                 else
791                 {
792                     double lam = fep->init_lambda;
793                     if (fep->init_lambda < 0)
794                     {
795                         lam = fep->all_lambda[i][fep->init_fep_state];
796                     }
797                     derivative = gmx::formatString("%s %s = %.4f", dhdl, enumValueToStringSingular(i), lam);
798                 }
799                 setname[s++] = derivative;
800             }
801         }
802     }
803
804     if (fep->n_lambda > 0)
805     {
806         /* g_bar has to determine the lambda values used in this simulation
807          * from this xvg legend.
808          */
809
810         if (expand->elmcmove > LambdaMoveCalculation::No)
811         {
812             nsetsbegin = 1; /* for including the expanded ensemble */
813         }
814         else
815         {
816             nsetsbegin = 0;
817         }
818
819         if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
820         {
821             nsetsbegin += 1;
822         }
823         nsetsbegin += nsets_dhdl;
824
825         for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
826         {
827             print_lambda_vector(fep, i, false, false, lambda_vec_str);
828             std::string buf;
829             if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
830             {
831                 /* for compatible dhdl.xvg files */
832                 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
833             }
834             else
835             {
836                 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
837             }
838
839             if (ir->bSimTemp)
840             {
841                 /* print the temperature for this state if doing simulated annealing */
842                 buf += gmx::formatString(
843                         "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
844             }
845             setname[s++] = buf;
846         }
847         if (write_pV)
848         {
849             setname[s++] = gmx::formatString("pV (%s)", unit_energy);
850         }
851
852         xvgrLegend(fp, setname, oenv);
853     }
854
855     return fp;
856 }
857
858 namespace gmx
859 {
860
861 void EnergyOutput::addDataAtEnergyStep(bool                    bDoDHDL,
862                                        bool                    bSum,
863                                        double                  time,
864                                        real                    tmass,
865                                        const gmx_enerdata_t*   enerd,
866                                        const t_lambda*         fep,
867                                        const t_expanded*       expand,
868                                        const matrix            box,
869                                        PTCouplingArrays        ptCouplingArrays,
870                                        int                     fep_state,
871                                        const tensor            svir,
872                                        const tensor            fvir,
873                                        const tensor            vir,
874                                        const tensor            pres,
875                                        const gmx_ekindata_t*   ekind,
876                                        const rvec              mu_tot,
877                                        const gmx::Constraints* constr)
878 {
879     int  j, k, kk, n, gid;
880     real crmsd[2], tmp6[6];
881     real bs[tricl_boxs_nm.size()], vol, dens, enthalpy;
882     real eee[static_cast<int>(NonBondedEnergyTerms::Count)];
883     gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, double> store_dhdl;
884     real                                                              store_energy = 0;
885     real                                                              tmp;
886     real pv = 0.0; // static analyzer warns about uninitialized variable warnings here.
887
888     /* Do NOT use the box in the state variable, but the separate box provided
889      * as an argument. This is because we sometimes need to write the box from
890      * the last timestep to match the trajectory frames.
891      */
892     add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
893     if (nCrmsd_)
894     {
895         crmsd[0] = constr->rmsd();
896         add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
897     }
898     if (bDynBox_)
899     {
900         int nboxs;
901         if (bTricl_)
902         {
903             bs[0] = box[XX][XX];
904             bs[1] = box[YY][YY];
905             bs[2] = box[ZZ][ZZ];
906             bs[3] = box[YY][XX];
907             bs[4] = box[ZZ][XX];
908             bs[5] = box[ZZ][YY];
909             nboxs = tricl_boxs_nm.size();
910         }
911         else
912         {
913             bs[0] = box[XX][XX];
914             bs[1] = box[YY][YY];
915             bs[2] = box[ZZ][ZZ];
916             nboxs = boxs_nm.size();
917         }
918         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
919         dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
920         add_ebin(ebin_, ib_, nboxs, bs, bSum);
921         add_ebin(ebin_, ivol_, 1, &vol, bSum);
922         add_ebin(ebin_, idens_, 1, &dens, bSum);
923
924         if (bDiagPres_)
925         {
926             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
927                not the instantaneous pressure */
928             pv = vol * ref_p_ / gmx::c_presfac;
929
930             add_ebin(ebin_, ipv_, 1, &pv, bSum);
931             enthalpy = pv + enerd->term[F_ETOT];
932             add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
933         }
934     }
935     if (bConstrVir_)
936     {
937         add_ebin(ebin_, isvir_, 9, svir[0], bSum);
938         add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
939     }
940     if (bPres_)
941     {
942         add_ebin(ebin_, ivir_, 9, vir[0], bSum);
943         add_ebin(ebin_, ipres_, 9, pres[0], bSum);
944         tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
945         add_ebin(ebin_, isurft_, 1, &tmp, bSum);
946     }
947     if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
948     {
949         tmp6[0] = ptCouplingArrays.boxv[XX][XX];
950         tmp6[1] = ptCouplingArrays.boxv[YY][YY];
951         tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
952         tmp6[3] = ptCouplingArrays.boxv[YY][XX];
953         tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
954         tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
955         add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
956     }
957     if (bMu_)
958     {
959         add_ebin(ebin_, imu_, 3, mu_tot, bSum);
960     }
961     if (ekind && ekind->cosacc.cos_accel != 0)
962     {
963         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
964         dens = (tmass * gmx::c_amu) / (vol * gmx::c_nano * gmx::c_nano * gmx::c_nano);
965         add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
966         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
967         tmp = 1
968               / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * gmx::c_pico) * dens
969                  * gmx::square(box[ZZ][ZZ] * gmx::c_nano / (2 * M_PI)));
970         add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
971     }
972     if (nE_ > 1)
973     {
974         n = 0;
975         for (int i = 0; (i < nEg_); i++)
976         {
977             for (j = i; (j < nEg_); j++)
978             {
979                 gid = GID(i, j, nEg_);
980                 for (k = kk = 0; (k < static_cast<int>(NonBondedEnergyTerms::Count)); k++)
981                 {
982                     if (bEInd_[k])
983                     {
984                         eee[kk++] = enerd->grpp.energyGroupPairTerms[k][gid];
985                     }
986                 }
987                 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
988                 n++;
989             }
990         }
991     }
992
993     if (ekind)
994     {
995         for (int i = 0; (i < nTC_); i++)
996         {
997             tmp_r_[i] = ekind->tcstat[i].T;
998         }
999         add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
1000
1001         if (etc_ == TemperatureCoupling::NoseHoover)
1002         {
1003             /* whether to print Nose-Hoover chains: */
1004             if (bPrintNHChains_)
1005             {
1006                 if (bNHC_trotter_)
1007                 {
1008                     for (int i = 0; (i < nTC_); i++)
1009                     {
1010                         for (j = 0; j < nNHC_; j++)
1011                         {
1012                             k                 = i * nNHC_ + j;
1013                             tmp_r_[2 * k]     = ptCouplingArrays.nosehoover_xi[k];
1014                             tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
1015                         }
1016                     }
1017                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1018
1019                     if (bMTTK_)
1020                     {
1021                         for (int i = 0; (i < nTCP_); i++)
1022                         {
1023                             for (j = 0; j < nNHC_; j++)
1024                             {
1025                                 k                 = i * nNHC_ + j;
1026                                 tmp_r_[2 * k]     = ptCouplingArrays.nhpres_xi[k];
1027                                 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
1028                             }
1029                         }
1030                         add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1031                     }
1032                 }
1033                 else
1034                 {
1035                     for (int i = 0; (i < nTC_); i++)
1036                     {
1037                         tmp_r_[2 * i]     = ptCouplingArrays.nosehoover_xi[i];
1038                         tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1039                     }
1040                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1041                 }
1042             }
1043         }
1044         else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
1045                  || etc_ == TemperatureCoupling::VRescale)
1046         {
1047             for (int i = 0; (i < nTC_); i++)
1048             {
1049                 tmp_r_[i] = ekind->tcstat[i].lambda;
1050             }
1051             add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1052         }
1053     }
1054
1055     ebin_increase_count(1, ebin_, bSum);
1056
1057     // BAR + thermodynamic integration values
1058     if ((fp_dhdl_ || dhc_) && bDoDHDL)
1059     {
1060         const auto& foreignTerms = enerd->foreignLambdaTerms;
1061         for (int i = 0; i < foreignTerms.numLambdas(); i++)
1062         {
1063             /* zero for simulated tempering */
1064             dE_[i] = foreignTerms.deltaH(i);
1065             if (numTemperatures_ > 0)
1066             {
1067                 GMX_RELEASE_ASSERT(numTemperatures_ > fep_state,
1068                                    "Number of lambdas in state is bigger then in input record");
1069                 GMX_RELEASE_ASSERT(
1070                         numTemperatures_ >= foreignTerms.numLambdas(),
1071                         "Number of lambdas in energy data is bigger then in input record");
1072                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1073                 /* is this even useful to have at all? */
1074                 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1075             }
1076         }
1077
1078         if (fp_dhdl_)
1079         {
1080             fprintf(fp_dhdl_, "%.4f", time);
1081             /* the current free energy state */
1082
1083             /* print the current state if we are doing expanded ensemble */
1084             if (expand->elmcmove > LambdaMoveCalculation::No)
1085             {
1086                 fprintf(fp_dhdl_, " %4d", fep_state);
1087             }
1088             /* total energy (for if the temperature changes */
1089
1090             if (fep->edHdLPrintEnergy != FreeEnergyPrintEnergy::No)
1091             {
1092                 switch (fep->edHdLPrintEnergy)
1093                 {
1094                     case FreeEnergyPrintEnergy::Potential:
1095                         store_energy = enerd->term[F_EPOT];
1096                         break;
1097                     case FreeEnergyPrintEnergy::Total:
1098                     case FreeEnergyPrintEnergy::Yes:
1099                     default: store_energy = enerd->term[F_ETOT];
1100                 }
1101                 fprintf(fp_dhdl_, " %#.8g", store_energy);
1102             }
1103
1104             if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
1105             {
1106                 for (auto i : keysOf(fep->separate_dvdl))
1107                 {
1108                     if (fep->separate_dvdl[i])
1109                     {
1110                         /* assumes F_DVDL is first */
1111                         fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + static_cast<int>(i)]);
1112                     }
1113                 }
1114             }
1115             for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1116             {
1117                 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1118             }
1119             if (bDynBox_ && bDiagPres_ && (epc_ != PressureCoupling::No)
1120                 && foreignTerms.numLambdas() > 0 && (fep->init_lambda < 0))
1121             {
1122                 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1123                                                          there are alternate state
1124                                                          lambda and we're not in
1125                                                          compatibility mode */
1126             }
1127             fprintf(fp_dhdl_, "\n");
1128             /* and the binary free energy output */
1129         }
1130         if (dhc_ && bDoDHDL)
1131         {
1132             int idhdl = 0;
1133             for (auto i : keysOf(fep->separate_dvdl))
1134             {
1135                 if (fep->separate_dvdl[i])
1136                 {
1137                     /* assumes F_DVDL is first */
1138                     store_dhdl[idhdl] = enerd->term[F_DVDL + static_cast<int>(i)];
1139                     idhdl += 1;
1140                 }
1141             }
1142             store_energy = enerd->term[F_ETOT];
1143             /* store_dh is dE */
1144             mde_delta_h_coll_add_dh(
1145                     dhc_, static_cast<double>(fep_state), store_energy, pv, store_dhdl, dE_ + fep->lambda_start_n, time);
1146         }
1147     }
1148
1149     if (conservedEnergyTracker_)
1150     {
1151         conservedEnergyTracker_->addPoint(
1152                 time, bEner_[F_ECONSERVED] ? enerd->term[F_ECONSERVED] : enerd->term[F_ETOT]);
1153     }
1154 }
1155
1156 void EnergyOutput::recordNonEnergyStep()
1157 {
1158     ebin_increase_count(1, ebin_, false);
1159 }
1160
1161 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1162 {
1163     char buf[22];
1164
1165     fprintf(log,
1166             "   %12s   %12s\n"
1167             "   %12s   %12.5f\n\n",
1168             "Step",
1169             "Time",
1170             gmx_step_str(steps, buf),
1171             time);
1172 }
1173
1174 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1175                                          bool       bEne,
1176                                          bool       bDR,
1177                                          bool       bOR,
1178                                          FILE*      log,
1179                                          int64_t    step,
1180                                          double     time,
1181                                          t_fcdata*  fcd,
1182                                          gmx::Awh*  awh)
1183 {
1184     t_enxframe fr;
1185     init_enxframe(&fr);
1186     fr.t       = time;
1187     fr.step    = step;
1188     fr.nsteps  = ebin_->nsteps;
1189     fr.dt      = delta_t_;
1190     fr.nsum    = ebin_->nsum;
1191     fr.nre     = (bEne) ? ebin_->nener : 0;
1192     fr.ener    = ebin_->e;
1193     int ndisre = bDR ? fcd->disres->npair : 0;
1194     /* these are for the old-style blocks (1 subblock, only reals), because
1195        there can be only one per ID for these */
1196     int   nr[enxNR];
1197     int   id[enxNR];
1198     real* block[enxNR];
1199     /* Optional additional old-style (real-only) blocks. */
1200     for (int i = 0; i < enxNR; i++)
1201     {
1202         nr[i] = 0;
1203     }
1204
1205     if (bOR && fcd->orires->nr > 0)
1206     {
1207         t_oriresdata& orires = *fcd->orires;
1208         diagonalize_orires_tensors(&orires);
1209         nr[enxOR]     = orires.nr;
1210         block[enxOR]  = orires.otav;
1211         id[enxOR]     = enxOR;
1212         nr[enxORI]    = (orires.oinsl != orires.otav) ? orires.nr : 0;
1213         block[enxORI] = orires.oinsl;
1214         id[enxORI]    = enxORI;
1215         nr[enxORT]    = orires.nex * 12;
1216         block[enxORT] = orires.eig;
1217         id[enxORT]    = enxORT;
1218     }
1219
1220     /* whether we are going to write anything out: */
1221     if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1222     {
1223         /* the old-style blocks go first */
1224         fr.nblock = 0;
1225         for (int i = 0; i < enxNR; i++)
1226         {
1227             if (nr[i] > 0)
1228             {
1229                 fr.nblock = i + 1;
1230             }
1231         }
1232         add_blocks_enxframe(&fr, fr.nblock);
1233         for (int b = 0; b < fr.nblock; b++)
1234         {
1235             add_subblocks_enxblock(&(fr.block[b]), 1);
1236             fr.block[b].id        = id[b];
1237             fr.block[b].sub[0].nr = nr[b];
1238 #if !GMX_DOUBLE
1239             fr.block[b].sub[0].type = xdr_datatype_float;
1240             fr.block[b].sub[0].fval = block[b];
1241 #else
1242             fr.block[b].sub[0].type  = xdr_datatype_double;
1243             fr.block[b].sub[0].dval  = block[b];
1244 #endif
1245         }
1246
1247         /* check for disre block & fill it. */
1248         if (ndisre > 0)
1249         {
1250             int db = fr.nblock;
1251             fr.nblock += 1;
1252             add_blocks_enxframe(&fr, fr.nblock);
1253
1254             add_subblocks_enxblock(&(fr.block[db]), 2);
1255             const t_disresdata& disres = *fcd->disres;
1256             fr.block[db].id            = enxDISRE;
1257             fr.block[db].sub[0].nr     = ndisre;
1258             fr.block[db].sub[1].nr     = ndisre;
1259 #if !GMX_DOUBLE
1260             fr.block[db].sub[0].type = xdr_datatype_float;
1261             fr.block[db].sub[1].type = xdr_datatype_float;
1262             fr.block[db].sub[0].fval = disres.rt;
1263             fr.block[db].sub[1].fval = disres.rm3tav;
1264 #else
1265             fr.block[db].sub[0].type = xdr_datatype_double;
1266             fr.block[db].sub[1].type = xdr_datatype_double;
1267             fr.block[db].sub[0].dval = disres.rt;
1268             fr.block[db].sub[1].dval = disres.rm3tav;
1269 #endif
1270         }
1271         /* here we can put new-style blocks */
1272
1273         /* Free energy perturbation blocks */
1274         if (dhc_)
1275         {
1276             mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1277         }
1278
1279         /* we can now free & reset the data in the blocks */
1280         if (dhc_)
1281         {
1282             mde_delta_h_coll_reset(dhc_);
1283         }
1284
1285         /* AWH bias blocks. */
1286         if (awh != nullptr) // TODO: add boolean flag.
1287         {
1288             awh->writeToEnergyFrame(step, &fr);
1289         }
1290
1291         /* do the actual I/O */
1292         do_enx(fp_ene, &fr);
1293         if (fr.nre)
1294         {
1295             /* We have stored the sums, so reset the sum history */
1296             reset_ebin_sums(ebin_);
1297         }
1298     }
1299     free_enxframe(&fr);
1300     if (log)
1301     {
1302         if (bOR && fcd->orires->nr > 0)
1303         {
1304             print_orires_log(log, fcd->orires);
1305         }
1306
1307         fprintf(log, "   Energies (%s)\n", unit_energy);
1308         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1309         fprintf(log, "\n");
1310     }
1311 }
1312
1313 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts)
1314 {
1315     if (log)
1316     {
1317         if (opts)
1318         {
1319             for (int i = 0; i < opts->ngtc; i++)
1320             {
1321                 if (opts->annealing[i] != SimulatedAnnealing::No)
1322                 {
1323                     fprintf(log,
1324                             "Current ref_t for group %s: %8.1f\n",
1325                             *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1326                             opts->ref_t[i]);
1327                 }
1328             }
1329             fprintf(log, "\n");
1330         }
1331     }
1332 }
1333
1334 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1335 {
1336     if (ebin_->nsum_sim <= 0)
1337     {
1338         if (log)
1339         {
1340             fprintf(log, "Not enough data recorded to report energy averages\n");
1341         }
1342         return;
1343     }
1344     if (log)
1345     {
1346
1347         char buf1[22], buf2[22];
1348
1349         fprintf(log, "\t<======  ###############  ==>\n");
1350         fprintf(log, "\t<====  A V E R A G E S  ====>\n");
1351         fprintf(log, "\t<==  ###############  ======>\n\n");
1352
1353         fprintf(log,
1354                 "\tStatistics over %s steps using %s frames\n",
1355                 gmx_step_str(ebin_->nsteps_sim, buf1),
1356                 gmx_step_str(ebin_->nsum_sim, buf2));
1357         fprintf(log, "\n");
1358
1359         fprintf(log, "   Energies (%s)\n", unit_energy);
1360         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1361         fprintf(log, "\n");
1362
1363         if (bDynBox_)
1364         {
1365             pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1366             fprintf(log, "\n");
1367         }
1368         if (bConstrVir_)
1369         {
1370             fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1371             pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1372             fprintf(log, "\n");
1373             fprintf(log, "   Force Virial (%s)\n", unit_energy);
1374             pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1375             fprintf(log, "\n");
1376         }
1377         if (bPres_)
1378         {
1379             fprintf(log, "   Total Virial (%s)\n", unit_energy);
1380             pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1381             fprintf(log, "\n");
1382             fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1383             pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1384             fprintf(log, "\n");
1385         }
1386         if (bMu_)
1387         {
1388             fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1389             pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1390             fprintf(log, "\n");
1391         }
1392
1393         if (nE_ > 1)
1394         {
1395             int padding = 8 - strlen(unit_energy);
1396             fprintf(log, "%*sEpot (%s)   ", padding, "", unit_energy);
1397             for (auto key : keysOf(bEInd_))
1398             {
1399                 if (bEInd_[key])
1400                 {
1401                     fprintf(log, "%12s   ", enumValueToString(key));
1402                 }
1403             }
1404             fprintf(log, "\n");
1405
1406             int n = 0;
1407             for (int i = 0; (i < nEg_); i++)
1408             {
1409                 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1410                 for (int j = i; (j < nEg_); j++)
1411                 {
1412                     int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1413                     int padding =
1414                             14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1415                     fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1416                     pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1417                     n++;
1418                 }
1419             }
1420             fprintf(log, "\n");
1421         }
1422         if (nTC_ > 1)
1423         {
1424             pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1425             fprintf(log, "\n");
1426         }
1427     }
1428 }
1429
1430 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1431 {
1432     const t_ebin* const ebin = ebin_;
1433
1434     enerhist->nsteps     = ebin->nsteps;
1435     enerhist->nsum       = ebin->nsum;
1436     enerhist->nsteps_sim = ebin->nsteps_sim;
1437     enerhist->nsum_sim   = ebin->nsum_sim;
1438
1439     if (ebin->nsum > 0)
1440     {
1441         /* This will only actually resize the first time */
1442         enerhist->ener_ave.resize(ebin->nener);
1443         enerhist->ener_sum.resize(ebin->nener);
1444
1445         for (int i = 0; i < ebin->nener; i++)
1446         {
1447             enerhist->ener_ave[i] = ebin->e[i].eav;
1448             enerhist->ener_sum[i] = ebin->e[i].esum;
1449         }
1450     }
1451
1452     if (ebin->nsum_sim > 0)
1453     {
1454         /* This will only actually resize the first time */
1455         enerhist->ener_sum_sim.resize(ebin->nener);
1456
1457         for (int i = 0; i < ebin->nener; i++)
1458         {
1459             enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1460         }
1461     }
1462     if (dhc_)
1463     {
1464         mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1465     }
1466 }
1467
1468 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1469 {
1470     unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1471
1472     if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1473         || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1474     {
1475         gmx_fatal(FARGS,
1476                   "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1477                   "or %zu).",
1478                   nener,
1479                   enerhist.ener_sum.size(),
1480                   enerhist.ener_sum_sim.size());
1481     }
1482
1483     ebin_->nsteps     = enerhist.nsteps;
1484     ebin_->nsum       = enerhist.nsum;
1485     ebin_->nsteps_sim = enerhist.nsteps_sim;
1486     ebin_->nsum_sim   = enerhist.nsum_sim;
1487
1488     for (int i = 0; i < ebin_->nener; i++)
1489     {
1490         ebin_->e[i].eav      = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1491         ebin_->e[i].esum     = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1492         ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1493     }
1494     if (dhc_)
1495     {
1496         mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1497     }
1498 }
1499
1500 int EnergyOutput::numEnergyTerms() const
1501 {
1502     return ebin_->nener;
1503 }
1504
1505 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1506 {
1507     if (fplog == nullptr)
1508     {
1509         return;
1510     }
1511
1512     if (conservedEnergyTracker_)
1513     {
1514         std::string partName = formatString("simulation part #%d", simulationPart);
1515         fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1516     }
1517     else if (usingMdIntegrator)
1518     {
1519         fprintf(fplog,
1520                 "\nCannot report drift of the conserved energy quantity because simulations share "
1521                 "state\n\n");
1522     }
1523 }
1524
1525 } // namespace gmx