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