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