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