Apply re-formatting to C++ in src/ tree.
[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/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 ? epcNO : 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 ? etcNO : 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_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
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_ == etcNOSEHOOVER)
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_ == epcMTTK)
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_ == etcNOSEHOOVER)
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_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
504     {
505         for (i = 0; (i < nTC_); i++)
506         {
507             ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling][i];
508             sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
509             grpnms[i] = gmx_strdup(buf);
510         }
511         itc_      = get_ebin_space(ebin_, mde_n_, grpnms, "");
512         allocated = mde_n_;
513     }
514
515     for (i = 0; i < allocated; i++)
516     {
517         sfree(grpnms[i]);
518     }
519     sfree(grpnms);
520
521     nU_ = groups->groups[SimulationAtomGroupType::Acceleration].size();
522     snew(tmp_v_, nU_);
523     if (nU_ > 1)
524     {
525         snew(grpnms, 3 * nU_);
526         for (i = 0; (i < nU_); i++)
527         {
528             ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
529             sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
530             grpnms[3 * i + XX] = gmx_strdup(buf);
531             sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
532             grpnms[3 * i + YY] = gmx_strdup(buf);
533             sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
534             grpnms[3 * i + ZZ] = gmx_strdup(buf);
535         }
536         iu_ = get_ebin_space(ebin_, 3 * nU_, grpnms, unit_vel);
537         for (i = 0; i < 3 * nU_; i++)
538         {
539             sfree(grpnms[i]);
540         }
541         sfree(grpnms);
542     }
543
544     /* Note that fp_ene should be valid on the master rank and null otherwise */
545     if (fp_ene != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
546     {
547         do_enxnms(fp_ene, &ebin_->nener, &ebin_->enm);
548     }
549
550     /* check whether we're going to write dh histograms */
551     dhc_ = nullptr;
552     if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
553     {
554         /* Currently dh histograms are only written with dynamics */
555         if (EI_DYNAMICS(ir->eI))
556         {
557             snew(dhc_, 1);
558
559             mde_delta_h_coll_init(dhc_, ir);
560         }
561         fp_dhdl_ = nullptr;
562         snew(dE_, ir->fepvals->n_lambda);
563     }
564     else
565     {
566         fp_dhdl_ = fp_dhdl;
567         snew(dE_, ir->fepvals->n_lambda);
568     }
569     if (ir->bSimTemp)
570     {
571         int i;
572         snew(temperatures_, ir->fepvals->n_lambda);
573         numTemperatures_ = ir->fepvals->n_lambda;
574         for (i = 0; i < ir->fepvals->n_lambda; i++)
575         {
576             temperatures_[i] = ir->simtempvals->temperatures[i];
577         }
578     }
579     else
580     {
581         numTemperatures_ = 0;
582     }
583
584     if (EI_MD(ir->eI) && !simulationsShareState)
585     {
586         conservedEnergyTracker_ = std::make_unique<EnergyDriftTracker>(mtop->natoms);
587     }
588 }
589
590 EnergyOutput::~EnergyOutput()
591 {
592     sfree(igrp_);
593     sfree(tmp_r_);
594     sfree(tmp_v_);
595     done_ebin(ebin_);
596     done_mde_delta_h_coll(dhc_);
597     sfree(dE_);
598     if (numTemperatures_ > 0)
599     {
600         sfree(temperatures_);
601     }
602 }
603
604 } // namespace gmx
605
606 /*! \brief Print a lambda vector to a string
607  *
608  * \param[in] fep                The inputrec's FEP input data
609  * \param[in] i                  The index of the lambda vector
610  * \param[in] get_native_lambda  Whether to print the native lambda
611  * \param[in] get_names          Whether to print the names rather than the values
612  * \param[in,out] str            The pre-allocated string buffer to print to.
613  */
614 static void print_lambda_vector(t_lambda* fep, int i, bool get_native_lambda, bool get_names, char* str)
615 {
616     int j, k = 0;
617     int Nsep = 0;
618
619     for (j = 0; j < efptNR; j++)
620     {
621         if (fep->separate_dvdl[j])
622         {
623             Nsep++;
624         }
625     }
626     str[0] = 0; /* reset the string */
627     if (Nsep > 1)
628     {
629         str += sprintf(str, "("); /* set the opening parenthesis*/
630     }
631     for (j = 0; j < efptNR; j++)
632     {
633         if (fep->separate_dvdl[j])
634         {
635             if (!get_names)
636             {
637                 if (get_native_lambda && fep->init_lambda >= 0)
638                 {
639                     str += sprintf(str, "%.4f", fep->init_lambda);
640                 }
641                 else
642                 {
643                     str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
644                 }
645             }
646             else
647             {
648                 str += sprintf(str, "%s", efpt_singular_names[j]);
649             }
650             /* print comma for the next item */
651             if (k < Nsep - 1)
652             {
653                 str += sprintf(str, ", ");
654             }
655             k++;
656         }
657     }
658     if (Nsep > 1)
659     {
660         /* and add the closing parenthesis */
661         sprintf(str, ")");
662     }
663 }
664
665 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv)
666 {
667     FILE*       fp;
668     const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
669                *lambdastate = "\\lambda state";
670     int         i, nsets, nsets_de, nsetsbegin;
671     int         n_lambda_terms = 0;
672     t_lambda*   fep            = ir->fepvals; /* for simplicity */
673     t_expanded* expand         = ir->expandedvals;
674     char        lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
675
676     int  nsets_dhdl = 0;
677     int  s          = 0;
678     int  nsetsextend;
679     bool write_pV = false;
680
681     /* count the number of different lambda terms */
682     for (i = 0; i < efptNR; i++)
683     {
684         if (fep->separate_dvdl[i])
685         {
686             n_lambda_terms++;
687         }
688     }
689
690     std::string title, label_x, label_y;
691     if (fep->n_lambda == 0)
692     {
693         title   = gmx::formatString("%s", dhdl);
694         label_x = gmx::formatString("Time (ps)");
695         label_y = gmx::formatString("%s (%s %s)", dhdl, unit_energy, "[\\lambda]\\S-1\\N");
696     }
697     else
698     {
699         title   = gmx::formatString("%s and %s", dhdl, deltag);
700         label_x = gmx::formatString("Time (ps)");
701         label_y = gmx::formatString(
702                 "%s and %s (%s %s)", dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
703     }
704     fp = gmx_fio_fopen(filename, "w+");
705     xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
706
707     std::string buf;
708     if (!(ir->bSimTemp))
709     {
710         buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
711     }
712     if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
713     {
714         if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
715         {
716             /* compatibility output */
717             buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
718         }
719         else
720         {
721             print_lambda_vector(fep, fep->init_fep_state, true, false, lambda_vec_str);
722             print_lambda_vector(fep, fep->init_fep_state, true, true, lambda_name_str);
723             buf += gmx::formatString(
724                     "%s %d: %s = %s", lambdastate, fep->init_fep_state, lambda_name_str, lambda_vec_str);
725         }
726     }
727     xvgr_subtitle(fp, buf.c_str(), oenv);
728
729
730     nsets_dhdl = 0;
731     if (fep->dhdl_derivatives == edhdlderivativesYES)
732     {
733         nsets_dhdl = n_lambda_terms;
734     }
735     /* count the number of delta_g states */
736     nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
737
738     nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
739
740     if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
741     {
742         nsets += 1; /*add fep state for expanded ensemble */
743     }
744
745     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
746     {
747         nsets += 1; /* add energy to the dhdl as well */
748     }
749
750     nsetsextend = nsets;
751     if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
752     {
753         nsetsextend += 1; /* for PV term, other terms possible if required for
754                              the reduced potential (only needed with foreign
755                              lambda, and only output when init_lambda is not
756                              set in order to maintain compatibility of the
757                              dhdl.xvg file) */
758         write_pV = true;
759     }
760     std::vector<std::string> setname(nsetsextend);
761
762     if (expand->elmcmove > elmcmoveNO)
763     {
764         /* state for the fep_vals, if we have alchemical sampling */
765         setname[s++] = "Thermodynamic state";
766     }
767
768     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
769     {
770         std::string energy;
771         switch (fep->edHdLPrintEnergy)
772         {
773             case edHdLPrintEnergyPOTENTIAL:
774                 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
775                 break;
776             case edHdLPrintEnergyTOTAL:
777             case edHdLPrintEnergyYES:
778             default: energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
779         }
780         setname[s++] = energy;
781     }
782
783     if (fep->dhdl_derivatives == edhdlderivativesYES)
784     {
785         for (i = 0; i < efptNR; i++)
786         {
787             if (fep->separate_dvdl[i])
788             {
789                 std::string derivative;
790                 if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
791                 {
792                     /* compatibility output */
793                     derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
794                 }
795                 else
796                 {
797                     double lam = fep->init_lambda;
798                     if (fep->init_lambda < 0)
799                     {
800                         lam = fep->all_lambda[i][fep->init_fep_state];
801                     }
802                     derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i], lam);
803                 }
804                 setname[s++] = derivative;
805             }
806         }
807     }
808
809     if (fep->n_lambda > 0)
810     {
811         /* g_bar has to determine the lambda values used in this simulation
812          * from this xvg legend.
813          */
814
815         if (expand->elmcmove > elmcmoveNO)
816         {
817             nsetsbegin = 1; /* for including the expanded ensemble */
818         }
819         else
820         {
821             nsetsbegin = 0;
822         }
823
824         if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
825         {
826             nsetsbegin += 1;
827         }
828         nsetsbegin += nsets_dhdl;
829
830         for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
831         {
832             print_lambda_vector(fep, i, false, false, lambda_vec_str);
833             std::string buf;
834             if ((fep->init_lambda >= 0) && (n_lambda_terms == 1))
835             {
836                 /* for compatible dhdl.xvg files */
837                 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
838             }
839             else
840             {
841                 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
842             }
843
844             if (ir->bSimTemp)
845             {
846                 /* print the temperature for this state if doing simulated annealing */
847                 buf += gmx::formatString(
848                         "T = %g (%s)", ir->simtempvals->temperatures[s - (nsetsbegin)], unit_temp_K);
849             }
850             setname[s++] = buf;
851         }
852         if (write_pV)
853         {
854             setname[s++] = gmx::formatString("pV (%s)", unit_energy);
855         }
856
857         xvgrLegend(fp, setname, oenv);
858     }
859
860     return fp;
861 }
862
863 namespace gmx
864 {
865
866 void EnergyOutput::addDataAtEnergyStep(bool                    bDoDHDL,
867                                        bool                    bSum,
868                                        double                  time,
869                                        real                    tmass,
870                                        const gmx_enerdata_t*   enerd,
871                                        const t_lambda*         fep,
872                                        const t_expanded*       expand,
873                                        const matrix            box,
874                                        PTCouplingArrays        ptCouplingArrays,
875                                        int                     fep_state,
876                                        const tensor            svir,
877                                        const tensor            fvir,
878                                        const tensor            vir,
879                                        const tensor            pres,
880                                        const gmx_ekindata_t*   ekind,
881                                        const rvec              mu_tot,
882                                        const gmx::Constraints* constr)
883 {
884     int    j, k, kk, n, gid;
885     real   crmsd[2], tmp6[6];
886     real   bs[tricl_boxs_nm.size()], vol, dens, pv, enthalpy;
887     real   eee[egNR];
888     double store_dhdl[efptNR];
889     real   store_energy = 0;
890     real   tmp;
891
892     /* Do NOT use the box in the state variable, but the separate box provided
893      * as an argument. This is because we sometimes need to write the box from
894      * the last timestep to match the trajectory frames.
895      */
896     add_ebin_indexed(ebin_, ie_, gmx::ArrayRef<bool>(bEner_), enerd->term, bSum);
897     if (nCrmsd_)
898     {
899         crmsd[0] = constr->rmsd();
900         add_ebin(ebin_, iconrmsd_, nCrmsd_, crmsd, false);
901     }
902     if (bDynBox_)
903     {
904         int nboxs;
905         if (bTricl_)
906         {
907             bs[0] = box[XX][XX];
908             bs[1] = box[YY][YY];
909             bs[2] = box[ZZ][ZZ];
910             bs[3] = box[YY][XX];
911             bs[4] = box[ZZ][XX];
912             bs[5] = box[ZZ][YY];
913             nboxs = tricl_boxs_nm.size();
914         }
915         else
916         {
917             bs[0] = box[XX][XX];
918             bs[1] = box[YY][YY];
919             bs[2] = box[ZZ][ZZ];
920             nboxs = boxs_nm.size();
921         }
922         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
923         dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
924         add_ebin(ebin_, ib_, nboxs, bs, bSum);
925         add_ebin(ebin_, ivol_, 1, &vol, bSum);
926         add_ebin(ebin_, idens_, 1, &dens, bSum);
927
928         if (bDiagPres_)
929         {
930             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
931                not the instantaneous pressure */
932             pv = vol * ref_p_ / PRESFAC;
933
934             add_ebin(ebin_, ipv_, 1, &pv, bSum);
935             enthalpy = pv + enerd->term[F_ETOT];
936             add_ebin(ebin_, ienthalpy_, 1, &enthalpy, bSum);
937         }
938     }
939     if (bConstrVir_)
940     {
941         add_ebin(ebin_, isvir_, 9, svir[0], bSum);
942         add_ebin(ebin_, ifvir_, 9, fvir[0], bSum);
943     }
944     if (bPres_)
945     {
946         add_ebin(ebin_, ivir_, 9, vir[0], bSum);
947         add_ebin(ebin_, ipres_, 9, pres[0], bSum);
948         tmp = (pres[ZZ][ZZ] - (pres[XX][XX] + pres[YY][YY]) * 0.5) * box[ZZ][ZZ];
949         add_ebin(ebin_, isurft_, 1, &tmp, bSum);
950     }
951     if (epc_ == epcPARRINELLORAHMAN || epc_ == epcMTTK)
952     {
953         tmp6[0] = ptCouplingArrays.boxv[XX][XX];
954         tmp6[1] = ptCouplingArrays.boxv[YY][YY];
955         tmp6[2] = ptCouplingArrays.boxv[ZZ][ZZ];
956         tmp6[3] = ptCouplingArrays.boxv[YY][XX];
957         tmp6[4] = ptCouplingArrays.boxv[ZZ][XX];
958         tmp6[5] = ptCouplingArrays.boxv[ZZ][YY];
959         add_ebin(ebin_, ipc_, bTricl_ ? 6 : 3, tmp6, bSum);
960     }
961     if (bMu_)
962     {
963         add_ebin(ebin_, imu_, 3, mu_tot, bSum);
964     }
965     if (ekind && ekind->cosacc.cos_accel != 0)
966     {
967         vol  = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
968         dens = (tmass * AMU) / (vol * NANO * NANO * NANO);
969         add_ebin(ebin_, ivcos_, 1, &(ekind->cosacc.vcos), bSum);
970         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
971         tmp = 1
972               / (ekind->cosacc.cos_accel / (ekind->cosacc.vcos * PICO) * dens
973                  * gmx::square(box[ZZ][ZZ] * NANO / (2 * M_PI)));
974         add_ebin(ebin_, ivisc_, 1, &tmp, bSum);
975     }
976     if (nE_ > 1)
977     {
978         n = 0;
979         for (int i = 0; (i < nEg_); i++)
980         {
981             for (j = i; (j < nEg_); j++)
982             {
983                 gid = GID(i, j, nEg_);
984                 for (k = kk = 0; (k < egNR); k++)
985                 {
986                     if (bEInd_[k])
987                     {
988                         eee[kk++] = enerd->grpp.ener[k][gid];
989                     }
990                 }
991                 add_ebin(ebin_, igrp_[n], nEc_, eee, bSum);
992                 n++;
993             }
994         }
995     }
996
997     if (ekind)
998     {
999         for (int i = 0; (i < nTC_); i++)
1000         {
1001             tmp_r_[i] = ekind->tcstat[i].T;
1002         }
1003         add_ebin(ebin_, itemp_, nTC_, tmp_r_, bSum);
1004
1005         if (etc_ == etcNOSEHOOVER)
1006         {
1007             /* whether to print Nose-Hoover chains: */
1008             if (bPrintNHChains_)
1009             {
1010                 if (bNHC_trotter_)
1011                 {
1012                     for (int i = 0; (i < nTC_); i++)
1013                     {
1014                         for (j = 0; j < nNHC_; j++)
1015                         {
1016                             k                 = i * nNHC_ + j;
1017                             tmp_r_[2 * k]     = ptCouplingArrays.nosehoover_xi[k];
1018                             tmp_r_[2 * k + 1] = ptCouplingArrays.nosehoover_vxi[k];
1019                         }
1020                     }
1021                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1022
1023                     if (bMTTK_)
1024                     {
1025                         for (int i = 0; (i < nTCP_); i++)
1026                         {
1027                             for (j = 0; j < nNHC_; j++)
1028                             {
1029                                 k                 = i * nNHC_ + j;
1030                                 tmp_r_[2 * k]     = ptCouplingArrays.nhpres_xi[k];
1031                                 tmp_r_[2 * k + 1] = ptCouplingArrays.nhpres_vxi[k];
1032                             }
1033                         }
1034                         add_ebin(ebin_, itcb_, mdeb_n_, tmp_r_, bSum);
1035                     }
1036                 }
1037                 else
1038                 {
1039                     for (int i = 0; (i < nTC_); i++)
1040                     {
1041                         tmp_r_[2 * i]     = ptCouplingArrays.nosehoover_xi[i];
1042                         tmp_r_[2 * i + 1] = ptCouplingArrays.nosehoover_vxi[i];
1043                     }
1044                     add_ebin(ebin_, itc_, mde_n_, tmp_r_, bSum);
1045                 }
1046             }
1047         }
1048         else if (etc_ == etcBERENDSEN || etc_ == etcYES || etc_ == etcVRESCALE)
1049         {
1050             for (int i = 0; (i < nTC_); i++)
1051             {
1052                 tmp_r_[i] = ekind->tcstat[i].lambda;
1053             }
1054             add_ebin(ebin_, itc_, nTC_, tmp_r_, bSum);
1055         }
1056     }
1057
1058     if (ekind && nU_ > 1)
1059     {
1060         for (int i = 0; (i < nU_); i++)
1061         {
1062             copy_rvec(ekind->grpstat[i].u, tmp_v_[i]);
1063         }
1064         add_ebin(ebin_, iu_, 3 * nU_, tmp_v_[0], bSum);
1065     }
1066
1067     ebin_increase_count(1, ebin_, bSum);
1068
1069     // BAR + thermodynamic integration values
1070     if ((fp_dhdl_ || dhc_) && bDoDHDL)
1071     {
1072         const auto& foreignTerms = enerd->foreignLambdaTerms;
1073         for (int i = 0; i < foreignTerms.numLambdas(); i++)
1074         {
1075             /* zero for simulated tempering */
1076             dE_[i] = foreignTerms.deltaH(i);
1077             if (numTemperatures_ > 0)
1078             {
1079                 GMX_RELEASE_ASSERT(numTemperatures_ > fep_state,
1080                                    "Number of lambdas in state is bigger then in input record");
1081                 GMX_RELEASE_ASSERT(
1082                         numTemperatures_ >= foreignTerms.numLambdas(),
1083                         "Number of lambdas in energy data is bigger then in input record");
1084                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1085                 /* is this even useful to have at all? */
1086                 dE_[i] += (temperatures_[i] / temperatures_[fep_state] - 1.0) * enerd->term[F_EKIN];
1087             }
1088         }
1089
1090         if (fp_dhdl_)
1091         {
1092             fprintf(fp_dhdl_, "%.4f", time);
1093             /* the current free energy state */
1094
1095             /* print the current state if we are doing expanded ensemble */
1096             if (expand->elmcmove > elmcmoveNO)
1097             {
1098                 fprintf(fp_dhdl_, " %4d", fep_state);
1099             }
1100             /* total energy (for if the temperature changes */
1101
1102             if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1103             {
1104                 switch (fep->edHdLPrintEnergy)
1105                 {
1106                     case edHdLPrintEnergyPOTENTIAL: store_energy = enerd->term[F_EPOT]; break;
1107                     case edHdLPrintEnergyTOTAL:
1108                     case edHdLPrintEnergyYES:
1109                     default: store_energy = enerd->term[F_ETOT];
1110                 }
1111                 fprintf(fp_dhdl_, " %#.8g", store_energy);
1112             }
1113
1114             if (fep->dhdl_derivatives == edhdlderivativesYES)
1115             {
1116                 for (int i = 0; i < efptNR; i++)
1117                 {
1118                     if (fep->separate_dvdl[i])
1119                     {
1120                         /* assumes F_DVDL is first */
1121                         fprintf(fp_dhdl_, " %#.8g", enerd->term[F_DVDL + i]);
1122                     }
1123                 }
1124             }
1125             for (int i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1126             {
1127                 fprintf(fp_dhdl_, " %#.8g", dE_[i]);
1128             }
1129             if (bDynBox_ && bDiagPres_ && (epc_ != epcNO) && foreignTerms.numLambdas() > 0
1130                 && (fep->init_lambda < 0))
1131             {
1132                 fprintf(fp_dhdl_, " %#.8g", pv); /* PV term only needed when
1133                                                          there are alternate state
1134                                                          lambda and we're not in
1135                                                          compatibility mode */
1136             }
1137             fprintf(fp_dhdl_, "\n");
1138             /* and the binary free energy output */
1139         }
1140         if (dhc_ && bDoDHDL)
1141         {
1142             int idhdl = 0;
1143             for (int i = 0; i < efptNR; i++)
1144             {
1145                 if (fep->separate_dvdl[i])
1146                 {
1147                     /* assumes F_DVDL is first */
1148                     store_dhdl[idhdl] = enerd->term[F_DVDL + i];
1149                     idhdl += 1;
1150                 }
1151             }
1152             store_energy = enerd->term[F_ETOT];
1153             /* store_dh is dE */
1154             mde_delta_h_coll_add_dh(
1155                     dhc_, static_cast<double>(fep_state), store_energy, pv, store_dhdl, dE_ + fep->lambda_start_n, time);
1156         }
1157     }
1158
1159     if (conservedEnergyTracker_)
1160     {
1161         conservedEnergyTracker_->addPoint(
1162                 time, bEner_[F_ECONSERVED] ? enerd->term[F_ECONSERVED] : enerd->term[F_ETOT]);
1163     }
1164 }
1165
1166 void EnergyOutput::recordNonEnergyStep()
1167 {
1168     ebin_increase_count(1, ebin_, false);
1169 }
1170
1171 void EnergyOutput::printHeader(FILE* log, int64_t steps, double time)
1172 {
1173     char buf[22];
1174
1175     fprintf(log,
1176             "   %12s   %12s\n"
1177             "   %12s   %12.5f\n\n",
1178             "Step",
1179             "Time",
1180             gmx_step_str(steps, buf),
1181             time);
1182 }
1183
1184 void EnergyOutput::printStepToEnergyFile(ener_file* fp_ene,
1185                                          bool       bEne,
1186                                          bool       bDR,
1187                                          bool       bOR,
1188                                          FILE*      log,
1189                                          int64_t    step,
1190                                          double     time,
1191                                          t_fcdata*  fcd,
1192                                          gmx::Awh*  awh)
1193 {
1194     t_enxframe fr;
1195     init_enxframe(&fr);
1196     fr.t       = time;
1197     fr.step    = step;
1198     fr.nsteps  = ebin_->nsteps;
1199     fr.dt      = delta_t_;
1200     fr.nsum    = ebin_->nsum;
1201     fr.nre     = (bEne) ? ebin_->nener : 0;
1202     fr.ener    = ebin_->e;
1203     int ndisre = bDR ? fcd->disres->npair : 0;
1204     /* these are for the old-style blocks (1 subblock, only reals), because
1205        there can be only one per ID for these */
1206     int   nr[enxNR];
1207     int   id[enxNR];
1208     real* block[enxNR];
1209     /* Optional additional old-style (real-only) blocks. */
1210     for (int i = 0; i < enxNR; i++)
1211     {
1212         nr[i] = 0;
1213     }
1214
1215     if (bOR && fcd->orires->nr > 0)
1216     {
1217         t_oriresdata& orires = *fcd->orires;
1218         diagonalize_orires_tensors(&orires);
1219         nr[enxOR]     = orires.nr;
1220         block[enxOR]  = orires.otav;
1221         id[enxOR]     = enxOR;
1222         nr[enxORI]    = (orires.oinsl != orires.otav) ? orires.nr : 0;
1223         block[enxORI] = orires.oinsl;
1224         id[enxORI]    = enxORI;
1225         nr[enxORT]    = orires.nex * 12;
1226         block[enxORT] = orires.eig;
1227         id[enxORT]    = enxORT;
1228     }
1229
1230     /* whether we are going to write anything out: */
1231     if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1232     {
1233         /* the old-style blocks go first */
1234         fr.nblock = 0;
1235         for (int i = 0; i < enxNR; i++)
1236         {
1237             if (nr[i] > 0)
1238             {
1239                 fr.nblock = i + 1;
1240             }
1241         }
1242         add_blocks_enxframe(&fr, fr.nblock);
1243         for (int b = 0; b < fr.nblock; b++)
1244         {
1245             add_subblocks_enxblock(&(fr.block[b]), 1);
1246             fr.block[b].id        = id[b];
1247             fr.block[b].sub[0].nr = nr[b];
1248 #if !GMX_DOUBLE
1249             fr.block[b].sub[0].type = xdr_datatype_float;
1250             fr.block[b].sub[0].fval = block[b];
1251 #else
1252             fr.block[b].sub[0].type  = xdr_datatype_double;
1253             fr.block[b].sub[0].dval  = block[b];
1254 #endif
1255         }
1256
1257         /* check for disre block & fill it. */
1258         if (ndisre > 0)
1259         {
1260             int db = fr.nblock;
1261             fr.nblock += 1;
1262             add_blocks_enxframe(&fr, fr.nblock);
1263
1264             add_subblocks_enxblock(&(fr.block[db]), 2);
1265             const t_disresdata& disres = *fcd->disres;
1266             fr.block[db].id            = enxDISRE;
1267             fr.block[db].sub[0].nr     = ndisre;
1268             fr.block[db].sub[1].nr     = ndisre;
1269 #if !GMX_DOUBLE
1270             fr.block[db].sub[0].type = xdr_datatype_float;
1271             fr.block[db].sub[1].type = xdr_datatype_float;
1272             fr.block[db].sub[0].fval = disres.rt;
1273             fr.block[db].sub[1].fval = disres.rm3tav;
1274 #else
1275             fr.block[db].sub[0].type = xdr_datatype_double;
1276             fr.block[db].sub[1].type = xdr_datatype_double;
1277             fr.block[db].sub[0].dval = disres.rt;
1278             fr.block[db].sub[1].dval = disres.rm3tav;
1279 #endif
1280         }
1281         /* here we can put new-style blocks */
1282
1283         /* Free energy perturbation blocks */
1284         if (dhc_)
1285         {
1286             mde_delta_h_coll_handle_block(dhc_, &fr, fr.nblock);
1287         }
1288
1289         /* we can now free & reset the data in the blocks */
1290         if (dhc_)
1291         {
1292             mde_delta_h_coll_reset(dhc_);
1293         }
1294
1295         /* AWH bias blocks. */
1296         if (awh != nullptr) // TODO: add boolean flag.
1297         {
1298             awh->writeToEnergyFrame(step, &fr);
1299         }
1300
1301         /* do the actual I/O */
1302         do_enx(fp_ene, &fr);
1303         if (fr.nre)
1304         {
1305             /* We have stored the sums, so reset the sum history */
1306             reset_ebin_sums(ebin_);
1307         }
1308     }
1309     free_enxframe(&fr);
1310     if (log)
1311     {
1312         if (bOR && fcd->orires->nr > 0)
1313         {
1314             print_orires_log(log, fcd->orires);
1315         }
1316
1317         fprintf(log, "   Energies (%s)\n", unit_energy);
1318         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprNORMAL, true);
1319         fprintf(log, "\n");
1320     }
1321 }
1322
1323 void EnergyOutput::printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, t_grpopts* opts)
1324 {
1325     if (log)
1326     {
1327         if (opts)
1328         {
1329             for (int i = 0; i < opts->ngtc; i++)
1330             {
1331                 if (opts->annealing[i] != eannNO)
1332                 {
1333                     fprintf(log,
1334                             "Current ref_t for group %s: %8.1f\n",
1335                             *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling][i]]),
1336                             opts->ref_t[i]);
1337                 }
1338             }
1339             fprintf(log, "\n");
1340         }
1341     }
1342 }
1343
1344 void EnergyOutput::printAverages(FILE* log, const SimulationGroups* groups)
1345 {
1346     if (ebin_->nsum_sim <= 0)
1347     {
1348         if (log)
1349         {
1350             fprintf(log, "Not enough data recorded to report energy averages\n");
1351         }
1352         return;
1353     }
1354     if (log)
1355     {
1356
1357         char buf1[22], buf2[22];
1358
1359         fprintf(log, "\t<======  ###############  ==>\n");
1360         fprintf(log, "\t<====  A V E R A G E S  ====>\n");
1361         fprintf(log, "\t<==  ###############  ======>\n\n");
1362
1363         fprintf(log,
1364                 "\tStatistics over %s steps using %s frames\n",
1365                 gmx_step_str(ebin_->nsteps_sim, buf1),
1366                 gmx_step_str(ebin_->nsum_sim, buf2));
1367         fprintf(log, "\n");
1368
1369         fprintf(log, "   Energies (%s)\n", unit_energy);
1370         pr_ebin(log, ebin_, ie_, f_nre_ + nCrmsd_, 5, eprAVER, true);
1371         fprintf(log, "\n");
1372
1373         if (bDynBox_)
1374         {
1375             pr_ebin(log, ebin_, ib_, bTricl_ ? tricl_boxs_nm.size() : boxs_nm.size(), 5, eprAVER, true);
1376             fprintf(log, "\n");
1377         }
1378         if (bConstrVir_)
1379         {
1380             fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1381             pr_ebin(log, ebin_, isvir_, 9, 3, eprAVER, false);
1382             fprintf(log, "\n");
1383             fprintf(log, "   Force Virial (%s)\n", unit_energy);
1384             pr_ebin(log, ebin_, ifvir_, 9, 3, eprAVER, false);
1385             fprintf(log, "\n");
1386         }
1387         if (bPres_)
1388         {
1389             fprintf(log, "   Total Virial (%s)\n", unit_energy);
1390             pr_ebin(log, ebin_, ivir_, 9, 3, eprAVER, false);
1391             fprintf(log, "\n");
1392             fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1393             pr_ebin(log, ebin_, ipres_, 9, 3, eprAVER, false);
1394             fprintf(log, "\n");
1395         }
1396         if (bMu_)
1397         {
1398             fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1399             pr_ebin(log, ebin_, imu_, 3, 3, eprAVER, false);
1400             fprintf(log, "\n");
1401         }
1402
1403         if (nE_ > 1)
1404         {
1405             int padding = 8 - strlen(unit_energy);
1406             fprintf(log, "%*sEpot (%s)   ", padding, "", unit_energy);
1407             for (int i = 0; (i < egNR); i++)
1408             {
1409                 if (bEInd_[i])
1410                 {
1411                     fprintf(log, "%12s   ", egrp_nm[i]);
1412                 }
1413             }
1414             fprintf(log, "\n");
1415
1416             int n = 0;
1417             for (int i = 0; (i < nEg_); i++)
1418             {
1419                 int ni = groups->groups[SimulationAtomGroupType::EnergyOutput][i];
1420                 for (int j = i; (j < nEg_); j++)
1421                 {
1422                     int nj = groups->groups[SimulationAtomGroupType::EnergyOutput][j];
1423                     int padding =
1424                             14 - (strlen(*(groups->groupNames[ni])) + strlen(*(groups->groupNames[nj])));
1425                     fprintf(log, "%*s%s-%s", padding, "", *(groups->groupNames[ni]), *(groups->groupNames[nj]));
1426                     pr_ebin(log, ebin_, igrp_[n], nEc_, nEc_, eprAVER, false);
1427                     n++;
1428                 }
1429             }
1430             fprintf(log, "\n");
1431         }
1432         if (nTC_ > 1)
1433         {
1434             pr_ebin(log, ebin_, itemp_, nTC_, 4, eprAVER, true);
1435             fprintf(log, "\n");
1436         }
1437         if (nU_ > 1)
1438         {
1439             fprintf(log, "%15s   %12s   %12s   %12s\n", "Group", "Ux", "Uy", "Uz");
1440             for (int i = 0; (i < nU_); i++)
1441             {
1442                 int ni = groups->groups[SimulationAtomGroupType::Acceleration][i];
1443                 fprintf(log, "%15s", *groups->groupNames[ni]);
1444                 pr_ebin(log, ebin_, iu_ + 3 * i, 3, 3, eprAVER, false);
1445             }
1446             fprintf(log, "\n");
1447         }
1448     }
1449 }
1450
1451 void EnergyOutput::fillEnergyHistory(energyhistory_t* enerhist) const
1452 {
1453     const t_ebin* const ebin = ebin_;
1454
1455     enerhist->nsteps     = ebin->nsteps;
1456     enerhist->nsum       = ebin->nsum;
1457     enerhist->nsteps_sim = ebin->nsteps_sim;
1458     enerhist->nsum_sim   = ebin->nsum_sim;
1459
1460     if (ebin->nsum > 0)
1461     {
1462         /* This will only actually resize the first time */
1463         enerhist->ener_ave.resize(ebin->nener);
1464         enerhist->ener_sum.resize(ebin->nener);
1465
1466         for (int i = 0; i < ebin->nener; i++)
1467         {
1468             enerhist->ener_ave[i] = ebin->e[i].eav;
1469             enerhist->ener_sum[i] = ebin->e[i].esum;
1470         }
1471     }
1472
1473     if (ebin->nsum_sim > 0)
1474     {
1475         /* This will only actually resize the first time */
1476         enerhist->ener_sum_sim.resize(ebin->nener);
1477
1478         for (int i = 0; i < ebin->nener; i++)
1479         {
1480             enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1481         }
1482     }
1483     if (dhc_)
1484     {
1485         mde_delta_h_coll_update_energyhistory(dhc_, enerhist);
1486     }
1487 }
1488
1489 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t& enerhist)
1490 {
1491     unsigned int nener = static_cast<unsigned int>(ebin_->nener);
1492
1493     if ((enerhist.nsum > 0 && nener != enerhist.ener_sum.size())
1494         || (enerhist.nsum_sim > 0 && nener != enerhist.ener_sum_sim.size()))
1495     {
1496         gmx_fatal(FARGS,
1497                   "Mismatch between number of energies in run input (%u) and checkpoint file (%zu "
1498                   "or %zu).",
1499                   nener,
1500                   enerhist.ener_sum.size(),
1501                   enerhist.ener_sum_sim.size());
1502     }
1503
1504     ebin_->nsteps     = enerhist.nsteps;
1505     ebin_->nsum       = enerhist.nsum;
1506     ebin_->nsteps_sim = enerhist.nsteps_sim;
1507     ebin_->nsum_sim   = enerhist.nsum_sim;
1508
1509     for (int i = 0; i < ebin_->nener; i++)
1510     {
1511         ebin_->e[i].eav      = (enerhist.nsum > 0 ? enerhist.ener_ave[i] : 0);
1512         ebin_->e[i].esum     = (enerhist.nsum > 0 ? enerhist.ener_sum[i] : 0);
1513         ebin_->e_sim[i].esum = (enerhist.nsum_sim > 0 ? enerhist.ener_sum_sim[i] : 0);
1514     }
1515     if (dhc_)
1516     {
1517         mde_delta_h_coll_restore_energyhistory(dhc_, enerhist.deltaHForeignLambdas.get());
1518     }
1519 }
1520
1521 int EnergyOutput::numEnergyTerms() const
1522 {
1523     return ebin_->nener;
1524 }
1525
1526 void EnergyOutput::printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const
1527 {
1528     if (fplog == nullptr)
1529     {
1530         return;
1531     }
1532
1533     if (conservedEnergyTracker_)
1534     {
1535         std::string partName = formatString("simulation part #%d", simulationPart);
1536         fprintf(fplog, "\n%s\n", conservedEnergyTracker_->energyDriftString(partName).c_str());
1537     }
1538     else if (usingMdIntegrator)
1539     {
1540         fprintf(fplog,
1541                 "\nCannot report drift of the conserved energy quantity because simulations share "
1542                 "state\n\n");
1543     }
1544 }
1545
1546 } // namespace gmx