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