b66cb0a933dc37ca77ce0ed38491e8b9a7e56526
[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/fatalerror.h"
84 #include "gromacs/utility/mdmodulenotification.h"
85 #include "gromacs/utility/smalloc.h"
86 #include "gromacs/utility/stringutil.h"
87
88 #include "energydrifttracker.h"
89
90 //! Labels for energy file quantities
91 //! \{
92 static const char* conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
93
94 static std::array<const char*, 3> boxs_nm = { "Box-X", "Box-Y", "Box-Z" };
95
96 static std::array<const char*, 6> tricl_boxs_nm = { "Box-XX", "Box-YY", "Box-ZZ",
97                                                     "Box-YX", "Box-ZX", "Box-ZY" };
98
99 static const char* vol_nm[] = { "Volume" };
100
101 static const char* dens_nm[] = { "Density" };
102
103 static const char* pv_nm[] = { "pV" };
104
105 static const char* enthalpy_nm[] = { "Enthalpy" };
106
107 static std::array<const char*, 6> boxvel_nm = { "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
108                                                 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY" };
109
110 const char* egrp_nm[egNR + 1] = { "Coul-SR", "LJ-SR", "Buck-SR", "Coul-14", "LJ-14", nullptr };
111 //! \}
112
113 namespace gmx
114 {
115
116 /*! \brief Energy output class
117  *
118  * This is the collection of energy averages collected during mdrun, and to
119  * be written out to the .edr file.
120  *
121  * \todo Use more std containers.
122  * \todo Remove GMX_CONSTRAINTVIR
123  * \todo Write free-energy output also to energy file (after adding more tests)
124  */
125 EnergyOutput::EnergyOutput(ener_file*               fp_ene,
126                            const gmx_mtop_t*        mtop,
127                            const t_inputrec*        ir,
128                            const pull_t*            pull_work,
129                            FILE*                    fp_dhdl,
130                            bool                     isRerun,
131                            const StartingBehavior   startingBehavior,
132                            const bool               simulationsShareState,
133                            const MdModulesNotifier& mdModulesNotifier)
134 {
135     const char*        ener_nm[F_NRE];
136     static const char* vir_nm[]   = { "Vir-XX", "Vir-XY", "Vir-XZ", "Vir-YX", "Vir-YY",
137                                     "Vir-YZ", "Vir-ZX", "Vir-ZY", "Vir-ZZ" };
138     static const char* sv_nm[]    = { "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
139                                    "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
140                                    "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ" };
141     static const char* fv_nm[]    = { "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
142                                    "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
143                                    "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ" };
144     static const char* pres_nm[]  = { "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
145                                      "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ" };
146     static const char* surft_nm[] = { "#Surf*SurfTen" };
147     static const char* mu_nm[]    = { "Mu-X", "Mu-Y", "Mu-Z" };
148     static const char* vcos_nm[]  = { "2CosZ*Vel-X" };
149     static const char* visc_nm[]  = { "1/Viscosity" };
150     static const char* baro_nm[]  = { "Barostat" };
151
152     const SimulationGroups* groups;
153     char**                  gnm;
154     char                    buf[256];
155     const char*             bufi;
156     int                     i, j, ni, nj, n, k, kk, ncon, nset;
157     bool                    bBHAM, b14;
158
159     if (EI_DYNAMICS(ir->eI))
160     {
161         delta_t_ = ir->delta_t;
162     }
163     else
164     {
165         delta_t_ = 0;
166     }
167
168     groups = &mtop->groups;
169
170     bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
171     b14   = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
172
173     ncon         = gmx_mtop_ftype_count(mtop, F_CONSTR);
174     nset         = gmx_mtop_ftype_count(mtop, F_SETTLE);
175     bool bConstr = (ncon > 0 || nset > 0) && !isRerun;
176     bConstrVir_  = false;
177     nCrmsd_      = 0;
178     if (bConstr)
179     {
180         if (ncon > 0 && ir->eConstrAlg == econtLINCS)
181         {
182             nCrmsd_ = 1;
183         }
184         bConstrVir_ = (getenv("GMX_CONSTRAINTVIR") != nullptr);
185     }
186     else
187     {
188         nCrmsd_ = 0;
189     }
190
191     /* Energy monitoring */
192     for (i = 0; i < egNR; i++)
193     {
194         bEInd_[i] = false;
195     }
196
197     // Setting true only to those energy terms, that have active interactions and
198     // are not vsite terms (not VSITE2, VSITE3, VSITE3FD, VSITE3FAD, VSITE3OUT, VSITE4FD, VSITE4FDN, or VSITEN)
199     for (i = 0; i < F_NRE; i++)
200     {
201         bEner_[i] = (gmx_mtop_ftype_count(mtop, i) > 0)
202                     && ((interaction_function[i].flags & IF_VSITE) == 0);
203     }
204
205     if (!isRerun)
206     {
207         bEner_[F_EKIN] = EI_DYNAMICS(ir->eI);
208         bEner_[F_ETOT] = EI_DYNAMICS(ir->eI);
209         bEner_[F_TEMP] = EI_DYNAMICS(ir->eI);
210
211         bEner_[F_ECONSERVED] = integratorHasConservedEnergyQuantity(ir);
212         bEner_[F_PDISPCORR]  = (ir->eDispCorr != edispcNO);
213         bEner_[F_PRES]       = true;
214     }
215
216     bEner_[F_LJ]           = !bBHAM;
217     bEner_[F_BHAM]         = bBHAM;
218     bEner_[F_EQM]          = ir->bQMMM;
219     bEner_[F_RF_EXCL]      = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
220     bEner_[F_COUL_RECIP]   = EEL_FULL(ir->coulombtype);
221     bEner_[F_LJ_RECIP]     = EVDW_PME(ir->vdwtype);
222     bEner_[F_LJ14]         = b14;
223     bEner_[F_COUL14]       = b14;
224     bEner_[F_LJC14_Q]      = false;
225     bEner_[F_LJC_PAIRS_NB] = false;
226
227
228     bEner_[F_DVDL_COUL]      = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptCOUL];
229     bEner_[F_DVDL_VDW]       = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptVDW];
230     bEner_[F_DVDL_BONDED]    = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptBONDED];
231     bEner_[F_DVDL_RESTRAINT] = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptRESTRAINT];
232     bEner_[F_DKDL]           = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptMASS];
233     bEner_[F_DVDL]           = (ir->efep != efepNO) && ir->fepvals->separate_dvdl[efptFEP];
234
235     bEner_[F_CONSTR]   = false;
236     bEner_[F_CONSTRNC] = false;
237     bEner_[F_SETTLE]   = false;
238
239     bEner_[F_COUL_SR] = true;
240     bEner_[F_EPOT]    = true;
241
242     bEner_[F_DISPCORR]   = (ir->eDispCorr != edispcNO);
243     bEner_[F_DISRESVIOL] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
244     bEner_[F_ORIRESDEV]  = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
245     bEner_[F_COM_PULL]   = ((ir->bPull && pull_have_potential(*pull_work)) || ir->bRot);
246
247     MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
248     mdModulesNotifier.simulationSetupNotifications_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
249
250     bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
251
252
253     // Counting the energy terms that will be printed and saving their names
254     f_nre_ = 0;
255     for (i = 0; i < F_NRE; i++)
256     {
257         if (bEner_[i])
258         {
259             ener_nm[f_nre_] = interaction_function[i].longname;
260             f_nre_++;
261         }
262     }
263
264     epc_            = isRerun ? PressureCoupling::No : ir->epc;
265     bDiagPres_      = !TRICLINIC(ir->ref_p) && !isRerun;
266     ref_p_          = (ir->ref_p[XX][XX] + ir->ref_p[YY][YY] + ir->ref_p[ZZ][ZZ]) / DIM;
267     bTricl_         = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
268     bDynBox_        = inputrecDynamicBox(ir);
269     etc_            = isRerun ? TemperatureCoupling::No : ir->etc;
270     bNHC_trotter_   = inputrecNvtTrotter(ir) && !isRerun;
271     bPrintNHChains_ = ir->bPrintNHChains && !isRerun;
272     bMTTK_          = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
273     bMu_            = inputrecNeedMutot(ir);
274     bPres_          = !isRerun;
275
276     ebin_ = mk_ebin();
277     /* Pass NULL for unit to let get_ebin_space determine the units
278      * for interaction_function[i].longname
279      */
280     ie_ = get_ebin_space(ebin_, f_nre_, ener_nm, nullptr);
281     if (nCrmsd_)
282     {
283         /* This should be called directly after the call for ie_,
284          * such that iconrmsd_ follows directly in the list.
285          */
286         iconrmsd_ = get_ebin_space(ebin_, nCrmsd_, conrmsd_nm, "");
287     }
288     if (bDynBox_)
289     {
290         ib_    = get_ebin_space(ebin_,
291                              bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(),
292                              bTricl_ ? tricl_boxs_nm.data() : boxs_nm.data(),
293                              unit_length);
294         ivol_  = get_ebin_space(ebin_, 1, vol_nm, unit_volume);
295         idens_ = get_ebin_space(ebin_, 1, dens_nm, unit_density_SI);
296         if (bDiagPres_)
297         {
298             ipv_       = get_ebin_space(ebin_, 1, pv_nm, unit_energy);
299             ienthalpy_ = get_ebin_space(ebin_, 1, enthalpy_nm, unit_energy);
300         }
301     }
302     if (bConstrVir_)
303     {
304         isvir_ = get_ebin_space(ebin_, asize(sv_nm), sv_nm, unit_energy);
305         ifvir_ = get_ebin_space(ebin_, asize(fv_nm), fv_nm, unit_energy);
306     }
307     if (bPres_)
308     {
309         ivir_   = get_ebin_space(ebin_, asize(vir_nm), vir_nm, unit_energy);
310         ipres_  = get_ebin_space(ebin_, asize(pres_nm), pres_nm, unit_pres_bar);
311         isurft_ = get_ebin_space(ebin_, asize(surft_nm), surft_nm, unit_surft_bar);
312     }
313     if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
314     {
315         ipc_ = get_ebin_space(ebin_, bTricl_ ? boxvel_nm.size() : DIM, boxvel_nm.data(), unit_vel);
316     }
317     if (bMu_)
318     {
319         imu_ = get_ebin_space(ebin_, asize(mu_nm), mu_nm, unit_dipole_D);
320     }
321     if (ir->cos_accel != 0)
322     {
323         ivcos_ = get_ebin_space(ebin_, asize(vcos_nm), vcos_nm, unit_vel);
324         ivisc_ = get_ebin_space(ebin_, asize(visc_nm), visc_nm, unit_invvisc_SI);
325     }
326
327     /* Energy monitoring */
328     for (i = 0; i < egNR; i++)
329     {
330         bEInd_[i] = false;
331     }
332     bEInd_[egCOULSR] = true;
333     bEInd_[egLJSR]   = true;
334
335     if (bBHAM)
336     {
337         bEInd_[egLJSR]   = false;
338         bEInd_[egBHAMSR] = true;
339     }
340     if (b14)
341     {
342         bEInd_[egLJ14]   = true;
343         bEInd_[egCOUL14] = true;
344     }
345     nEc_ = 0;
346     for (i = 0; (i < egNR); i++)
347     {
348         if (bEInd_[i])
349         {
350             nEc_++;
351         }
352     }
353     n    = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
354     nEg_ = n;
355     nE_  = (n * (n + 1)) / 2;
356
357     snew(igrp_, nE_);
358     if (nE_ > 1)
359     {
360         n = 0;
361         snew(gnm, nEc_);
362         for (k = 0; (k < nEc_); k++)
363         {
364             snew(gnm[k], STRLEN);
365         }
366         for (i = 0; (i < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); i++)
367         {
368             ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
369             for (j = i; (j < gmx::ssize(groups->groups[SimulationAtomGroupType::EnergyOutput])); j++)
370             {
371                 nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
372                 for (k = kk = 0; (k < egNR); k++)
373                 {
374                     if (bEInd_[k])
375                     {
376                         sprintf(gnm[kk],
377                                 "%s:%s-%s",
378                                 egrp_nm[k],
379                                 *(groups->groupNames[ni]),
380                                 *(groups->groupNames[nj]));
381                         kk++;
382                     }
383                 }
384                 igrp_[n] = get_ebin_space(ebin_, nEc_, gnm, unit_energy);
385                 n++;
386             }
387         }
388         for (k = 0; (k < nEc_); k++)
389         {
390             sfree(gnm[k]);
391         }
392         sfree(gnm);
393
394         if (n != nE_)
395         {
396             gmx_incons("Number of energy terms wrong");
397         }
398     }
399
400     nTC_  = isRerun ? 0 : groups->groups[SimulationAtomGroupType::TemperatureCoupling].size();
401     nNHC_ = ir->opts.nhchainlength; /* shorthand for number of NH chains */
402     if (bMTTK_)
403     {
404         nTCP_ = 1; /* assume only one possible coupling system for barostat
405                            for now */
406     }
407     else
408     {
409         nTCP_ = 0;
410     }
411     if (etc_ == TemperatureCoupling::NoseHoover)
412     {
413         if (bNHC_trotter_)
414         {
415             mde_n_ = 2 * nNHC_ * nTC_;
416         }
417         else
418         {
419             mde_n_ = 2 * nTC_;
420         }
421         if (epc_ == PressureCoupling::Mttk)
422         {
423             mdeb_n_ = 2 * nNHC_ * nTCP_;
424         }
425     }
426     else
427     {
428         mde_n_  = nTC_;
429         mdeb_n_ = 0;
430     }
431
432     snew(tmp_r_, mde_n_);
433     // TODO redo the group name memory management to make it more clear
434     char** grpnms;
435     snew(grpnms, std::max(mde_n_, mdeb_n_)); // Just in case mdeb_n_ > mde_n_
436
437     for (i = 0; (i < nTC_); i++)
438     {
439         ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
440         sprintf(buf, "T-%s", *(groups->groupNames[ni]));
441         grpnms[i] = gmx_strdup(buf);
442     }
443     itemp_ = get_ebin_space(ebin_, nTC_, grpnms, unit_temp_K);
444     for (i = 0; i < nTC_; i++)
445     {
446         sfree(grpnms[i]);
447     }
448
449     int allocated = 0;
450     if (etc_ == TemperatureCoupling::NoseHoover)
451     {
452         if (bPrintNHChains_)
453         {
454             if (bNHC_trotter_)
455             {
456                 for (i = 0; (i < nTC_); i++)
457                 {
458                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
459                     bufi = *(groups->groupNames[ni]);
460                     for (j = 0; (j < nNHC_); j++)
461                     {
462                         sprintf(buf, "Xi-%d-%s", j, bufi);
463                         grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
464                         sprintf(buf, "vXi-%d-%s", j, bufi);
465                         grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
466                     }
467                 }
468                 itc_      = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
469                 allocated = mde_n_;
470                 if (bMTTK_)
471                 {
472                     for (i = 0; (i < nTCP_); i++)
473                     {
474                         bufi = baro_nm[0]; /* All barostat DOF's together for now. */
475                         for (j = 0; (j < nNHC_); j++)
476                         {
477                             sprintf(buf, "Xi-%d-%s", j, bufi);
478                             grpnms[2 * (i * nNHC_ + j)] = gmx_strdup(buf);
479                             sprintf(buf, "vXi-%d-%s", j, bufi);
480                             grpnms[2 * (i * nNHC_ + j) + 1] = gmx_strdup(buf);
481                         }
482                     }
483                     itcb_     = get_ebin_space(ebin_, mdeb_n_, grpnms, unit_invtime);
484                     allocated = mdeb_n_;
485                 }
486             }
487             else
488             {
489                 for (i = 0; (i < nTC_); i++)
490                 {
491                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
492                     bufi = *(groups->groupNames[ni]);
493                     sprintf(buf, "Xi-%s", bufi);
494                     grpnms[2 * i] = gmx_strdup(buf);
495                     sprintf(buf, "vXi-%s", bufi);
496                     grpnms[2 * i + 1] = gmx_strdup(buf);
497                 }
498                 itc_      = get_ebin_space(ebin_, mde_n_, grpnms, unit_invtime);
499                 allocated = mde_n_;
500             }
501         }
502     }
503     else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
504              || etc_ == TemperatureCoupling::VRescale)
505     {
506         for (i = 0; (i < nTC_); i++)
507         {
508             ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
509             sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
510             grpnms[i] = gmx_strdup(buf);
511         }
512         itc_      = get_ebin_space(ebin_, mde_n_, grpnms, "");
513         allocated = mde_n_;
514     }
515
516     for (i = 0; i < allocated; i++)
517     {
518         sfree(grpnms[i]);
519     }
520     sfree(grpnms);
521
522     /* Note that fp_ene should be valid on the master rank and null otherwise */
523     if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
524     {
525         do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
526     }
527
528     /* check whether we're going to write dh histograms */
529     dhc_ = nullptr;
530     if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
531     {
532         /* Currently dh histograms are only written with dynamics */
533         if (EI_DYNAMICS(ir->eI))
534         {
535             snew(dhc_, 1);
536
537             mde_delta_h_coll_init(dhc_, ir);
538         }
539         fp_dhdl_ = nullptr;
540         snew(dE_, ir->fepvals->n_lambda);
541     }
542     else
543     {
544         fp_dhdl_ = fp_dhdl;
545         snew(dE_, ir->fepvals->n_lambda);
546     }
547     if (ir->bSimTemp)
548     {
549         int i;
550         snew(temperatures_, ir->fepvals->n_lambda);
551         numTemperatures_ = ir->fepvals->n_lambda;
552         for (i = 0; i < ir->fepvals->n_lambda; i++)
553         {
554             temperatures_[i] = ir->simtempvals->temperatures[i];
555         }
556     }
557     else
558     {
559         numTemperatures_ = 0;
560     }
561
562     if (EI_MD(ir->eI) && !simulationsShareState)
563     {
564         conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop->natoms);
565     }
566 }
567
568 EnergyOutput::~EnergyOutput()
569 {
570     sfree(igrp_);
571     sfree(tmp_r_);
572     sfree(tmp_v_);
573     done_ebin(ebin_);
574     done_mde_delta_h_coll(dhc_);
575     sfree(dE_);
576     if (numTemperatures_ > 0)
577     {
578         sfree(temperatures_);
579     }
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 j, k = 0;
595     int Nsep = 0;
596
597     for (j = 0; j < efptNR; j++)
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 (j = 0; j < efptNR; j++)
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", efpt_singular_names[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; /* for simplicity */
651     t_expanded* expand         = ir->expandedvals;
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 (i = 0; i < efptNR; i++)
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 != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
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 == edhdlderivativesYES)
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 (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
719     {
720         nsets += 1; /*add fep state for expanded ensemble */
721     }
722
723     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
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 (expand->elmcmove > elmcmoveNO)
741     {
742         /* state for the fep_vals, if we have alchemical sampling */
743         setname[s++] = "Thermodynamic state";
744     }
745
746     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
747     {
748         std::string energy;
749         switch (fep->edHdLPrintEnergy)
750         {
751             case edHdLPrintEnergyPOTENTIAL:
752                 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
753                 break;
754             case edHdLPrintEnergyTOTAL:
755             case edHdLPrintEnergyYES:
756             default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
757         }
758         setname[s++] = energy;
759     }
760
761     if (fep->dhdl_derivatives == edhdlderivativesYES)
762     {
763         for (i = 0; i < efptNR; i++)
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, efpt_singular_names[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 (expand->elmcmove > elmcmoveNO)
794         {
795             nsetsbegin = 1; /* for including the expanded ensemble */
796         }
797         else
798         {
799             nsetsbegin = 0;
800         }
801
802         if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
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 t_expanded*       expand,
851                                        const matrix            box,
852                                        PTCouplingArrays        ptCouplingArrays,
853                                        int                     fep_state,
854                                        const tensor            svir,
855                                        const tensor            fvir,
856                                        const tensor            vir,
857                                        const tensor            pres,
858                                        const gmx_ekindata_t*   ekind,
859                                        const rvec              mu_tot,
860                                        const gmx::Constraints* constr)
861 {
862     int    j, k, kk, n, gid;
863     real   crmsd[2], tmp6[6];
864     real   bs[tricl_boxs_nm.size()], vol, dens, pv, enthalpy;
865     real   eee[egNR];
866     double store_dhdl[efptNR];
867     real   store_energy = 0;
868     real   tmp;
869
870     /* Do NOT use the box in the state variable, but the separate box provided
871      * as an argument. This is because we sometimes need to write the box from
872      * the last timestep to match the trajectory frames.
873      */
874     add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
875     if (nCrmsd_)
876     {
877         crmsd[0] = constr->rmsd();
878         add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
879     }
880     if (bDynBox_)
881     {
882         int nboxs;
883         if (bTricl_)
884         {
885             bs[0] = box[XX][XX];
886             bs[1] = box[YY][YY];
887             bs[2] = box[ZZ][ZZ];
888             bs[3] = box[YY][XX];
889             bs[4] = box[ZZ][XX];
890             bs[5] = box[ZZ][YY];
891             nboxs = tricl_boxs_nm.size();
892         }
893         else
894         {
895             bs[0] = box[XX][XX];
896             bs[1] = box[YY][YY];
897             bs[2] = box[ZZ][ZZ];
898             nboxs = boxs_nm.size();
899         }
900         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
901         dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
902         add_ebin(ebin_, ib_, nboxs, bs, bSum);
903         add_ebin(ebin_, ivol_, 1, &vol, bSum);
904         add_ebin(ebin_, idens_, 1, &dens, bSum);
905
906         if (bDiagPres_)
907         {
908             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
909                not the instantaneous pressure */
910             pv = vol * ref_p_ / PRESFAC;
911
912             add_ebin(ebin_, ipv_, 1, &pv, bSum);
913             enthalpy = pv + enerd->term[F_ETOT];
914             add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
915         }
916     }
917     if (bConstrVir_)
918     {
919         add_ebin(ebin_, isvir_, 9, svir[0], bSum);
920         add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
921     }
922     if (bPres_)
923     {
924         add_ebin(ebin_, ivir_, 9, vir[0], bSum);
925         add_ebin(ebin_, ipres_, 9, pres[0], bSum);
926         tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
927         add_ebin(ebin_, isurft_, 1, &tmp, bSum);
928     }
929     if (epc_ == PressureCoupling::ParrinelloRahman || epc_ == PressureCoupling::Mttk)
930     {
931         tmp6[0] = ptCouplingArrays.boxv[XX][XX];
932         tmp6[1] = ptCouplingArrays.boxv[YY][YY];
933         tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
934         tmp6[3] = ptCouplingArrays.boxv[YY][XX];
935         tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
936         tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
937         add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
938     }
939     if (bMu_)
940     {
941         add_ebin(ebin_, imu_, 3, mu_tot, bSum);
942     }
943     if (ekind && ekind->cosacc.cos_accel != 0)
944     {
945         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
946         dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
947         add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
948         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
949         tmp = 1
950               / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * PICO) * dens
951                  * gmx::square(box[ZZ][ZZ] * NANO / (2 * M_PI)));
952         add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
953     }
954     if (nE_ > 1)
955     {
956         n = 0;
957         for (int i = 0; (i < nEg_); i++)
958         {
959             for (j = i; (j < nEg_); j++)
960             {
961                 gid = GID(i, j, nEg_);
962                 for (k = kk = 0; (k < egNR); k++)
963                 {
964                     if (bEInd_[k])
965                     {
966                         eee[kk++] = enerd->grpp.ener[k][gid];
967                     }
968                 }
969                 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
970                 n++;
971             }
972         }
973     }
974
975     if (ekind)
976     {
977         for (int i = 0; (i < nTC_); i++)
978         {
979             tmp_r_[i] = ekind->tcstat[i].T;
980         }
981         add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
982
983         if (etc_ == TemperatureCoupling::NoseHoover)
984         {
985             /* whether to print Nose-Hoover chains: */
986             if (bPrintNHChains_)
987             {
988                 if (bNHC_trotter_)
989                 {
990                     for (int i = 0; (i < nTC_); i++)
991                     {
992                         for (j = 0; j < nNHC_; j++)
993                         {
994                             k                 = i * nNHC_ + j;
995                             tmp_r_[2 * k]     = ptCouplingArrays.nosehoover_xi[k];
996                             tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
997                         }
998                     }
999                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1000
1001                     if (bMTTK_)
1002                     {
1003                         for (int i = 0; (i < nTCP_); i++)
1004                         {
1005                             for (j = 0; j < nNHC_; j++)
1006                             {
1007                                 k                 = i * nNHC_ + j;
1008                                 tmp_r_[2 * k]     = ptCouplingArrays.nhpres_xi[k];
1009                                 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
1010                             }
1011                         }
1012                         add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1013                     }
1014                 }
1015                 else
1016                 {
1017                     for (int i = 0; (i < nTC_); i++)
1018                     {
1019                         tmp_r_[2 * i]     = ptCouplingArrays.nosehoover_xi[i];
1020                         tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1021                     }
1022                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1023                 }
1024             }
1025         }
1026         else if (etc_ == TemperatureCoupling::Berendsen || etc_ == TemperatureCoupling::Yes
1027                  || etc_ == TemperatureCoupling::VRescale)
1028         {
1029             for (int i = 0; (i < nTC_); i++)
1030             {
1031                 tmp_r_[i] = ekind->tcstat[i].lambda;
1032             }
1033             add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1034         }
1035     }
1036
1037     ebin_increase_count(1, ebin_, bSum);
1038
1039     // BAR + thermodynamic integration values
1040     if ((fp_dhdl_ || dhc_) && bDoDHDL)
1041     {
1042         const auto& foreignTerms = enerd->foreignLambdaTerms;
1043         for (int i = 0; i < foreignTerms.numLambdas(); i++)
1044         {
1045             /* zero for simulated tempering */
1046             dE_[i] = foreignTerms.deltaH(i);
1047             if (numTemperatures_ > 0)
1048             {
1049                 GMX_RELEASE_ASSERT(numTemperatures_ > fep_state,
1050                                    "Number of lambdas in state is bigger then in input record");
1051                 GMX_RELEASE_ASSERT(
1052                         numTemperatures_ >= foreignTerms.numLambdas(),
1053                         "Number of lambdas in energy data is bigger then in input record");
1054                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1055                 /* is this even useful to have at all? */
1056                 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1057             }
1058         }
1059
1060         if (fp_dhdl_)
1061         {
1062             fprintf(fp_dhdl_, "%.4f", time);
1063             /* the current free energy state */
1064
1065             /* print the current state if we are doing expanded ensemble */
1066             if (expand->elmcmove > elmcmoveNO)
1067             {
1068                 fprintf(fp_dhdl_, " %4d", fep_state);
1069             }
1070             /* total energy (for if the temperature changes */
1071
1072             if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1073             {
1074                 switch (fep->edHdLPrintEnergy)
1075                 {
1076                     case edHdLPrintEnergyPOTENTIAL: store_energy = enerd->term[F_EPOT]; break;
1077                     case edHdLPrintEnergyTOTAL:
1078                     case edHdLPrintEnergyYES:
1079                     default: store_energy = enerd->term[F_ETOT];
1080                 }
1081                 fprintf(fp_dhdl_, " %#.8g", store_energy);
1082             }
1083
1084             if (fep->dhdl_derivatives == edhdlderivativesYES)
1085             {
1086                 for (int i = 0; i < efptNR; i++)
1087                 {
1088                     if (fep->separate_dvdl[i])
1089                     {
1090                         /* assumes F_DVDL is first */
1091                         fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + i]);
1092                     }
1093                 }
1094             }
1095             for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1096             {
1097                 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1098             }
1099             if (bDynBox_ && bDiagPres_ && (epc_ != PressureCoupling::No)
1100                 && foreignTerms.numLambdas() > 0 && (fep->init_lambda < 0))
1101             {
1102                 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1103                                                          there are alternate state
1104                                                          lambda and we're not in
1105                                                          compatibility mode */
1106             }
1107             fprintf(fp_dhdl_, "\n");
1108             /* and the binary free energy output */
1109         }
1110         if (dhc_ && bDoDHDL)
1111         {
1112             int idhdl = 0;
1113             for (int i = 0; i < efptNR; i++)
1114             {
1115                 if (fep->separate_dvdl[i])
1116                 {
1117                     /* assumes F_DVDL is first */
1118                     store_dhdl[idhdl] = enerd->term[F_DVDL + i];
1119                     idhdl += 1;
1120                 }
1121             }
1122             store_energy = enerd->term[F_ETOT];
1123             /* store_dh is dE */
1124             mde_delta_h_coll_add_dh(
1125                     dhc_, static_cast<double>(fep_state), store_energy, pv, store_dhdl, dE_ + fep->lambda_start_n, 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->nr > 0)
1186     {
1187         t_oriresdata& orires = *fcd->orires;
1188         diagonalize_orires_tensors(&orires);
1189         nr[enxOR]     = orires.nr;
1190         block[enxOR]  = orires.otav;
1191         id[enxOR]     = enxOR;
1192         nr[enxORI]    = (orires.oinsl != orires.otav) ? orires.nr : 0;
1193         block[enxORI] = orires.oinsl;
1194         id[enxORI]    = enxORI;
1195         nr[enxORT]    = orires.nex * 12;
1196         block[enxORT] = orires.eig;
1197         id[enxORT]    = enxORT;
1198     }
1199
1200     /* whether we are going to write anything out: */
1201     if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1202     {
1203         /* the old-style blocks go first */
1204         fr.nblock = 0;
1205         for (int i = 0; i < enxNR; i++)
1206         {
1207             if (nr[i] > 0)
1208             {
1209                 fr.nblock = i + 1;
1210             }
1211         }
1212         add_blocks_enxframe(&fr, fr.nblock);
1213         for (int b = 0; b < fr.nblock; b++)
1214         {
1215             add_subblocks_enxblock(&(fr.block[b]), 1);
1216             fr.block[b].id        = id[b];
1217             fr.block[b].sub[0].nr = nr[b];
1218 #if !GMX_DOUBLE
1219             fr.block[b].sub[0].type = xdr_datatype_float;
1220             fr.block[b].sub[0].fval = block[b];
1221 #else
1222             fr.block[b].sub[0].type  = xdr_datatype_double;
1223             fr.block[b].sub[0].dval  = block[b];
1224 #endif
1225         }
1226
1227         /* check for disre block & fill it. */
1228         if (ndisre > 0)
1229         {
1230             int db = fr.nblock;
1231             fr.nblock += 1;
1232             add_blocks_enxframe(&fr, fr.nblock);
1233
1234             add_subblocks_enxblock(&(fr.block[db]), 2);
1235             const t_disresdata& disres = *fcd->disres;
1236             fr.block[db].id            = enxDISRE;
1237             fr.block[db].sub[0].nr     = ndisre;
1238             fr.block[db].sub[1].nr     = ndisre;
1239 #if !GMX_DOUBLE
1240             fr.block[db].sub[0].type = xdr_datatype_float;
1241             fr.block[db].sub[1].type = xdr_datatype_float;
1242             fr.block[db].sub[0].fval = disres.rt;
1243             fr.block[db].sub[1].fval = disres.rm3tav;
1244 #else
1245             fr.block[db].sub[0].type = xdr_datatype_double;
1246             fr.block[db].sub[1].type = xdr_datatype_double;
1247             fr.block[db].sub[0].dval = disres.rt;
1248             fr.block[db].sub[1].dval = disres.rm3tav;
1249 #endif
1250         }
1251         /* here we can put new-style blocks */
1252
1253         /* Free energy perturbation blocks */
1254         if (dhc_)
1255         {
1256             mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1257         }
1258
1259         /* we can now free & reset the data in the blocks */
1260         if (dhc_)
1261         {
1262             mde_delta_h_coll_reset(dhc_);
1263         }
1264
1265         /* AWH bias blocks. */
1266         if (awh != nullptr) // TODO: add boolean flag.
1267         {
1268             awh->writeToEnergyFrame(step, &fr);
1269         }
1270
1271         /* do the actual I/O */
1272         do_enx(fp_ene, &fr);
1273         if (fr.nre)
1274         {
1275             /* We have stored the sums, so reset the sum history */
1276             reset_ebin_sums(ebin_);
1277         }
1278     }
1279     free_enxframe(&fr);
1280     if (log)
1281     {
1282         if (bOR && fcd->orires->nr > 0)
1283         {
1284             print_orires_log(log, fcd->orires);
1285         }
1286
1287         fprintf(log, "   Energies (%s)\n", unit_energy);
1288         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1289         fprintf(log, "\n");
1290     }
1291 }
1292
1293 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts)
1294 {
1295     if (log)
1296     {
1297         if (opts)
1298         {
1299             for (int i = 0; i < opts->ngtc; i++)
1300             {
1301                 if (opts->annealing[i] != eannNO)
1302                 {
1303                     fprintf(log,
1304                             "Current ref_t for group %s: %8.1f\n",
1305                             *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1306                             opts->ref_t[i]);
1307                 }
1308             }
1309             fprintf(log, "\n");
1310         }
1311     }
1312 }
1313
1314 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1315 {
1316     if (ebin_->nsum_sim <= 0)
1317     {
1318         if (log)
1319         {
1320             fprintf(log, "Not enough data recorded to report energy averages\n");
1321         }
1322         return;
1323     }
1324     if (log)
1325     {
1326
1327         char buf1[22], buf2[22];
1328
1329         fprintf(log, "\t<======  ###############  ==>\n");
1330         fprintf(log, "\t<====  A V E R A G E S  ====>\n");
1331         fprintf(log, "\t<==  ###############  ======>\n\n");
1332
1333         fprintf(log,
1334                 "\tStatistics over %s steps using %s frames\n",
1335                 gmx_step_str(ebin_->nsteps_sim, buf1),
1336                 gmx_step_str(ebin_->nsum_sim, buf2));
1337         fprintf(log, "\n");
1338
1339         fprintf(log, "   Energies (%s)\n", unit_energy);
1340         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1341         fprintf(log, "\n");
1342
1343         if (bDynBox_)
1344         {
1345             pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1346             fprintf(log, "\n");
1347         }
1348         if (bConstrVir_)
1349         {
1350             fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1351             pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1352             fprintf(log, "\n");
1353             fprintf(log, "   Force Virial (%s)\n", unit_energy);
1354             pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1355             fprintf(log, "\n");
1356         }
1357         if (bPres_)
1358         {
1359             fprintf(log, "   Total Virial (%s)\n", unit_energy);
1360             pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1361             fprintf(log, "\n");
1362             fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1363             pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1364             fprintf(log, "\n");
1365         }
1366         if (bMu_)
1367         {
1368             fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1369             pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1370             fprintf(log, "\n");
1371         }
1372
1373         if (nE_ > 1)
1374         {
1375             int padding = 8 - strlen(unit_energy);
1376             fprintf(log, "%*sEpot (%s)   ", padding, "", unit_energy);
1377             for (int i = 0; (i < egNR); i++)
1378             {
1379                 if (bEInd_[i])
1380                 {
1381                     fprintf(log, "%12s   ", egrp_nm[i]);
1382                 }
1383             }
1384             fprintf(log, "\n");
1385
1386             int n = 0;
1387             for (int i = 0; (i < nEg_); i++)
1388             {
1389                 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1390                 for (int j = i; (j < nEg_); j++)
1391                 {
1392                     int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1393                     int padding =
1394                             14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1395                     fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1396                     pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1397                     n++;
1398                 }
1399             }
1400             fprintf(log, "\n");
1401         }
1402         if (nTC_ > 1)
1403         {
1404             pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1405             fprintf(log, "\n");
1406         }
1407     }
1408 }
1409
1410 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1411 {
1412     const t_ebin* const ebin = ebin_;
1413
1414     enerhist->nsteps     = ebin->nsteps;
1415     enerhist->nsum       = ebin->nsum;
1416     enerhist->nsteps_sim = ebin->nsteps_sim;
1417     enerhist->nsum_sim   = ebin->nsum_sim;
1418
1419     if (ebin->nsum > 0)
1420     {
1421         /* This will only actually resize the first time */
1422         enerhist->ener_ave.resize(ebin->nener);
1423         enerhist->ener_sum.resize(ebin->nener);
1424
1425         for (int i = 0; i < ebin->nener; i++)
1426         {
1427             enerhist->ener_ave[i] = ebin->e[i].eav;
1428             enerhist->ener_sum[i] = ebin->e[i].esum;
1429         }
1430     }
1431
1432     if (ebin->nsum_sim > 0)
1433     {
1434         /* This will only actually resize the first time */
1435         enerhist->ener_sum_sim.resize(ebin->nener);
1436
1437         for (int i = 0; i < ebin->nener; i++)
1438         {
1439             enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1440         }
1441     }
1442     if (dhc_)
1443     {
1444         mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1445     }
1446 }
1447
1448 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1449 {
1450     unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1451
1452     if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1453         || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1454     {
1455         gmx_fatal(FARGS,
1456                   "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1457                   "or %zu).",
1458                   nener,
1459                   enerhist.ener_sum.size(),
1460                   enerhist.ener_sum_sim.size());
1461     }
1462
1463     ebin_->nsteps     = enerhist.nsteps;
1464     ebin_->nsum       = enerhist.nsum;
1465     ebin_->nsteps_sim = enerhist.nsteps_sim;
1466     ebin_->nsum_sim   = enerhist.nsum_sim;
1467
1468     for (int i = 0; i < ebin_->nener; i++)
1469     {
1470         ebin_->e[i].eav      = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1471         ebin_->e[i].esum     = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1472         ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1473     }
1474     if (dhc_)
1475     {
1476         mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1477     }
1478 }
1479
1480 int EnergyOutput::numEnergyTerms() const
1481 {
1482     return ebin_->nener;
1483 }
1484
1485 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1486 {
1487     if (fplog == nullptr)
1488     {
1489         return;
1490     }
1491
1492     if (conservedEnergyTracker_)
1493     {
1494         std::string partName = formatString("simulation part #%d", simulationPart);
1495         fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1496     }
1497     else if (usingMdIntegrator)
1498     {
1499         fprintf(fplog,
1500                 "\nCannot report drift of the conserved energy quantity because simulations share "
1501                 "state\n\n");
1502     }
1503 }
1504
1505 } // namespace gmx