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