f90aabe1630fa82f45f46e729ae8187605d39da2
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "energyoutput.h"
40
41 #include <cfloat>
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <string>
46
47 #include "gromacs/awh/awh.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/listed_forces/disre.h"
53 #include "gromacs/listed_forces/orires.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/constr.h"
58 #include "gromacs/mdlib/ebin.h"
59 #include "gromacs/mdlib/mdebin_bar.h"
60 #include "gromacs/mdtypes/energyhistory.h"
61 #include "gromacs/mdtypes/fcdata.h"
62 #include "gromacs/mdtypes/group.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pulling/pull.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/trajectory/energyframe.h"
70 #include "gromacs/utility/arraysize.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
74
75 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
76
77 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
78
79 static const char *tricl_boxs_nm[] = {
80     "Box-XX", "Box-YY", "Box-ZZ",
81     "Box-YX", "Box-ZX", "Box-ZY"
82 };
83
84 static const char *vol_nm[] = { "Volume" };
85
86 static const char *dens_nm[] = {"Density" };
87
88 static const char *pv_nm[] = {"pV" };
89
90 static const char *enthalpy_nm[] = {"Enthalpy" };
91
92 static const char *boxvel_nm[] = {
93     "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
94     "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
95 };
96
97 #define NBOXS asize(boxs_nm)
98 #define NTRICLBOXS asize(tricl_boxs_nm)
99
100 const char *egrp_nm[egNR+1] = {
101     "Coul-SR", "LJ-SR", "Buck-SR",
102     "Coul-14", "LJ-14", nullptr
103 };
104
105 /* forward declaration */
106 typedef struct t_mde_delta_h_coll t_mde_delta_h_coll;
107
108 namespace gmx
109 {
110
111 namespace detail
112 {
113
114 /* This is the collection of energy averages collected during mdrun, and to
115    be written out to the .edr file. */
116 struct t_mdebin
117 {
118     double              delta_t;
119     t_ebin             *ebin;
120     int                 ie, iconrmsd, ib, ivol, idens, ipv, ienthalpy;
121     int                 isvir, ifvir, ipres, ivir, isurft, ipc, itemp, itc, itcb, iu, imu;
122     int                 ivcos, ivisc;
123     int                 nE, nEg, nEc, nTC, nTCP, nU, nNHC;
124     int                *igrp;
125     int                 mde_n, mdeb_n;
126     real               *tmp_r;
127     rvec               *tmp_v;
128     gmx_bool            bConstr;
129     gmx_bool            bConstrVir;
130     gmx_bool            bTricl;
131     gmx_bool            bDynBox;
132     gmx_bool            bNHC_trotter;
133     gmx_bool            bPrintNHChains;
134     gmx_bool            bMTTK;
135     gmx_bool            bMu; /* true if dipole is calculated */
136     gmx_bool            bDiagPres;
137     gmx_bool            bPres;
138     int                 f_nre;
139     int                 epc;
140     real                ref_p;
141     int                 etc;
142     int                 nCrmsd;
143     gmx_bool            bEner[F_NRE];
144     gmx_bool            bEInd[egNR];
145     char              **print_grpnms;
146
147     FILE               *fp_dhdl; /* the dhdl.xvg output file */
148     double             *dE;      /* energy components for dhdl.xvg output */
149     t_mde_delta_h_coll *dhc;     /* the delta U components (raw data + histogram) */
150     real               *temperatures;
151 };
152
153 }   // namespace detail
154
155 using detail::t_mdebin;
156
157 namespace
158 {
159
160 //! Legacy init function
161 t_mdebin *init_mdebin(ener_file_t       fp_ene,
162                       const gmx_mtop_t *mtop,
163                       const t_inputrec *ir,
164                       FILE             *fp_dhdl,
165                       bool              isRerun)
166 {
167     const char         *ener_nm[F_NRE];
168     static const char  *vir_nm[] = {
169         "Vir-XX", "Vir-XY", "Vir-XZ",
170         "Vir-YX", "Vir-YY", "Vir-YZ",
171         "Vir-ZX", "Vir-ZY", "Vir-ZZ"
172     };
173     static const char  *sv_nm[] = {
174         "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
175         "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
176         "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
177     };
178     static const char  *fv_nm[] = {
179         "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
180         "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
181         "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
182     };
183     static const char  *pres_nm[] = {
184         "Pres-XX", "Pres-XY", "Pres-XZ",
185         "Pres-YX", "Pres-YY", "Pres-YZ",
186         "Pres-ZX", "Pres-ZY", "Pres-ZZ"
187     };
188     static const char  *surft_nm[] = {
189         "#Surf*SurfTen"
190     };
191     static const char  *mu_nm[] = {
192         "Mu-X", "Mu-Y", "Mu-Z"
193     };
194     static const char  *vcos_nm[] = {
195         "2CosZ*Vel-X"
196     };
197     static const char  *visc_nm[] = {
198         "1/Viscosity"
199     };
200     static const char  *baro_nm[] = {
201         "Barostat"
202     };
203
204     const gmx_groups_t *groups;
205     char              **gnm;
206     char                buf[256];
207     const char         *bufi;
208     t_mdebin           *md;
209     int                 i, j, ni, nj, n, k, kk, ncon, nset;
210     gmx_bool            bBHAM, b14;
211
212     snew(md, 1);
213
214     if (EI_DYNAMICS(ir->eI))
215     {
216         md->delta_t = ir->delta_t;
217     }
218     else
219     {
220         md->delta_t = 0;
221     }
222
223     groups = &mtop->groups;
224
225     bBHAM = (mtop->ffparams.numTypes() > 0) && (mtop->ffparams.functype[0] == F_BHAM);
226     b14   = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
227              gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
228
229     ncon           = gmx_mtop_ftype_count(mtop, F_CONSTR);
230     nset           = gmx_mtop_ftype_count(mtop, F_SETTLE);
231     md->bConstr    = (ncon > 0 || nset > 0) && !isRerun;
232     md->bConstrVir = FALSE;
233     if (md->bConstr)
234     {
235         if (ncon > 0 && ir->eConstrAlg == econtLINCS)
236         {
237             md->nCrmsd = 1;
238         }
239         md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != nullptr);
240     }
241     else
242     {
243         md->nCrmsd = 0;
244     }
245
246     /* Energy monitoring */
247     for (i = 0; i < egNR; i++)
248     {
249         md->bEInd[i] = FALSE;
250     }
251
252     for (i = 0; i < F_NRE; i++)
253     {
254         md->bEner[i] = FALSE;
255         if (isRerun &&
256             (i == F_EKIN || i == F_ETOT || i == F_ECONSERVED ||
257              i == F_TEMP || i == F_PDISPCORR || i == F_PRES))
258         {
259             continue;
260         }
261         if (i == F_LJ)
262         {
263             md->bEner[i] = !bBHAM;
264         }
265         else if (i == F_BHAM)
266         {
267             md->bEner[i] = bBHAM;
268         }
269         else if (i == F_EQM)
270         {
271             md->bEner[i] = ir->bQMMM;
272         }
273         else if (i == F_RF_EXCL)
274         {
275             md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
276         }
277         else if (i == F_COUL_RECIP)
278         {
279             md->bEner[i] = EEL_FULL(ir->coulombtype);
280         }
281         else if (i == F_LJ_RECIP)
282         {
283             md->bEner[i] = EVDW_PME(ir->vdwtype);
284         }
285         else if (i == F_LJ14)
286         {
287             md->bEner[i] = b14;
288         }
289         else if (i == F_COUL14)
290         {
291             md->bEner[i] = b14;
292         }
293         else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
294         {
295             md->bEner[i] = FALSE;
296         }
297         else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
298                  (i == F_DVDL_VDW  && ir->fepvals->separate_dvdl[efptVDW]) ||
299                  (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
300                  (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
301                  (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
302                  (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
303         {
304             md->bEner[i] = (ir->efep != efepNO);
305         }
306         else if ((interaction_function[i].flags & IF_VSITE) ||
307                  (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
308         {
309             md->bEner[i] = FALSE;
310         }
311         else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES)  || (i == F_EQM))
312         {
313             md->bEner[i] = TRUE;
314         }
315         else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
316         {
317             md->bEner[i] = EI_DYNAMICS(ir->eI);
318         }
319         else if (i == F_DISPCORR || i == F_PDISPCORR)
320         {
321             md->bEner[i] = (ir->eDispCorr != edispcNO);
322         }
323         else if (i == F_DISRESVIOL)
324         {
325             md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
326         }
327         else if (i == F_ORIRESDEV)
328         {
329             md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
330         }
331         else if (i == F_CONNBONDS)
332         {
333             md->bEner[i] = FALSE;
334         }
335         else if (i == F_COM_PULL)
336         {
337             md->bEner[i] = ((ir->bPull && pull_have_potential(ir->pull_work)) ||
338                             ir->bRot);
339         }
340         else if (i == F_ECONSERVED)
341         {
342             md->bEner[i] = (integratorHasConservedEnergyQuantity(ir));
343         }
344         else
345         {
346             md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
347         }
348     }
349
350     md->f_nre = 0;
351     for (i = 0; i < F_NRE; i++)
352     {
353         if (md->bEner[i])
354         {
355             ener_nm[md->f_nre] = interaction_function[i].longname;
356             md->f_nre++;
357         }
358     }
359
360     md->epc            = isRerun ? epcNO : ir->epc;
361     md->bDiagPres      = !TRICLINIC(ir->ref_p) && !isRerun;
362     md->ref_p          = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
363     md->bTricl         = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
364     md->bDynBox        = inputrecDynamicBox(ir);
365     md->etc            = isRerun ? etcNO : ir->etc;
366     md->bNHC_trotter   = inputrecNvtTrotter(ir) && !isRerun;
367     md->bPrintNHChains = ir->bPrintNHChains && !isRerun;
368     md->bMTTK          = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)) && !isRerun;
369     md->bMu            = inputrecNeedMutot(ir);
370     md->bPres          = !isRerun;
371
372     md->ebin  = mk_ebin();
373     /* Pass NULL for unit to let get_ebin_space determine the units
374      * for interaction_function[i].longname
375      */
376     md->ie    = get_ebin_space(md->ebin, md->f_nre, ener_nm, nullptr);
377     if (md->nCrmsd)
378     {
379         /* This should be called directly after the call for md->ie,
380          * such that md->iconrmsd follows directly in the list.
381          */
382         md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
383     }
384     if (md->bDynBox)
385     {
386         md->ib    = get_ebin_space(md->ebin,
387                                    md->bTricl ? NTRICLBOXS : NBOXS,
388                                    md->bTricl ? tricl_boxs_nm : boxs_nm,
389                                    unit_length);
390         md->ivol  = get_ebin_space(md->ebin, 1, vol_nm,  unit_volume);
391         md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
392         if (md->bDiagPres)
393         {
394             md->ipv       = get_ebin_space(md->ebin, 1, pv_nm,   unit_energy);
395             md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm,   unit_energy);
396         }
397     }
398     if (md->bConstrVir)
399     {
400         md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
401         md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
402     }
403     if (md->bPres)
404     {
405         md->ivir   = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
406         md->ipres  = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
407         md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
408                                     unit_surft_bar);
409     }
410     if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
411     {
412         md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
413                                  boxvel_nm, unit_vel);
414     }
415     if (md->bMu)
416     {
417         md->imu    = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
418     }
419     if (ir->cos_accel != 0)
420     {
421         md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
422         md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
423                                    unit_invvisc_SI);
424     }
425
426     /* Energy monitoring */
427     for (i = 0; i < egNR; i++)
428     {
429         md->bEInd[i] = FALSE;
430     }
431     md->bEInd[egCOULSR] = TRUE;
432     md->bEInd[egLJSR  ] = TRUE;
433
434     if (bBHAM)
435     {
436         md->bEInd[egLJSR]   = FALSE;
437         md->bEInd[egBHAMSR] = TRUE;
438     }
439     if (b14)
440     {
441         md->bEInd[egLJ14]   = TRUE;
442         md->bEInd[egCOUL14] = TRUE;
443     }
444     md->nEc = 0;
445     for (i = 0; (i < egNR); i++)
446     {
447         if (md->bEInd[i])
448         {
449             md->nEc++;
450         }
451     }
452
453     n       = groups->grps[egcENER].nr;
454     md->nEg = n;
455     md->nE  = (n*(n+1))/2;
456
457     snew(md->igrp, md->nE);
458     if (md->nE > 1)
459     {
460         n = 0;
461         snew(gnm, md->nEc);
462         for (k = 0; (k < md->nEc); k++)
463         {
464             snew(gnm[k], STRLEN);
465         }
466         for (i = 0; (i < groups->grps[egcENER].nr); i++)
467         {
468             ni = groups->grps[egcENER].nm_ind[i];
469             for (j = i; (j < groups->grps[egcENER].nr); j++)
470             {
471                 nj = groups->grps[egcENER].nm_ind[j];
472                 for (k = kk = 0; (k < egNR); k++)
473                 {
474                     if (md->bEInd[k])
475                     {
476                         sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
477                                 *(groups->grpname[ni]), *(groups->grpname[nj]));
478                         kk++;
479                     }
480                 }
481                 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
482                                              gnm, unit_energy);
483                 n++;
484             }
485         }
486         for (k = 0; (k < md->nEc); k++)
487         {
488             sfree(gnm[k]);
489         }
490         sfree(gnm);
491
492         if (n != md->nE)
493         {
494             gmx_incons("Number of energy terms wrong");
495         }
496     }
497
498     md->nTC  = isRerun ? 0 : groups->grps[egcTC].nr;
499     md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
500     if (md->bMTTK)
501     {
502         md->nTCP = 1;  /* assume only one possible coupling system for barostat
503                           for now */
504     }
505     else
506     {
507         md->nTCP = 0;
508     }
509     if (md->etc == etcNOSEHOOVER)
510     {
511         if (md->bNHC_trotter)
512         {
513             md->mde_n = 2*md->nNHC*md->nTC;
514         }
515         else
516         {
517             md->mde_n = 2*md->nTC;
518         }
519         if (md->epc == epcMTTK)
520         {
521             md->mdeb_n = 2*md->nNHC*md->nTCP;
522         }
523     }
524     else
525     {
526         md->mde_n  = md->nTC;
527         md->mdeb_n = 0;
528     }
529
530     snew(md->tmp_r, md->mde_n);
531     snew(md->tmp_v, md->mde_n);
532     char **grpnms;
533     snew(grpnms, md->mde_n);
534
535     for (i = 0; (i < md->nTC); i++)
536     {
537         ni = groups->grps[egcTC].nm_ind[i];
538         sprintf(buf, "T-%s", *(groups->grpname[ni]));
539         grpnms[i] = gmx_strdup(buf);
540     }
541     md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
542                                unit_temp_K);
543
544     if (md->etc == etcNOSEHOOVER)
545     {
546         if (md->bPrintNHChains)
547         {
548             if (md->bNHC_trotter)
549             {
550                 for (i = 0; (i < md->nTC); i++)
551                 {
552                     ni   = groups->grps[egcTC].nm_ind[i];
553                     bufi = *(groups->grpname[ni]);
554                     for (j = 0; (j < md->nNHC); j++)
555                     {
556                         sprintf(buf, "Xi-%d-%s", j, bufi);
557                         grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
558                         sprintf(buf, "vXi-%d-%s", j, bufi);
559                         grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
560                     }
561                 }
562                 md->itc = get_ebin_space(md->ebin, md->mde_n,
563                                          grpnms, unit_invtime);
564                 if (md->bMTTK)
565                 {
566                     for (i = 0; (i < md->nTCP); i++)
567                     {
568                         bufi = baro_nm[0];  /* All barostat DOF's together for now. */
569                         for (j = 0; (j < md->nNHC); j++)
570                         {
571                             sprintf(buf, "Xi-%d-%s", j, bufi);
572                             grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
573                             sprintf(buf, "vXi-%d-%s", j, bufi);
574                             grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
575                         }
576                     }
577                     md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
578                                               grpnms, unit_invtime);
579                 }
580             }
581             else
582             {
583                 for (i = 0; (i < md->nTC); i++)
584                 {
585                     ni   = groups->grps[egcTC].nm_ind[i];
586                     bufi = *(groups->grpname[ni]);
587                     sprintf(buf, "Xi-%s", bufi);
588                     grpnms[2*i] = gmx_strdup(buf);
589                     sprintf(buf, "vXi-%s", bufi);
590                     grpnms[2*i+1] = gmx_strdup(buf);
591                 }
592                 md->itc = get_ebin_space(md->ebin, md->mde_n,
593                                          grpnms, unit_invtime);
594             }
595         }
596     }
597     else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
598              md->etc == etcVRESCALE)
599     {
600         for (i = 0; (i < md->nTC); i++)
601         {
602             ni = groups->grps[egcTC].nm_ind[i];
603             sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
604             grpnms[i] = gmx_strdup(buf);
605         }
606         md->itc = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
607     }
608
609     for (i = 0; i < md->mde_n; i++)
610     {
611         sfree(grpnms[i]);
612     }
613     sfree(grpnms);
614
615     md->nU = groups->grps[egcACC].nr;
616     if (md->nU > 1)
617     {
618         snew(grpnms, 3*md->nU);
619         for (i = 0; (i < md->nU); i++)
620         {
621             ni = groups->grps[egcACC].nm_ind[i];
622             sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
623             grpnms[3*i+XX] = gmx_strdup(buf);
624             sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
625             grpnms[3*i+YY] = gmx_strdup(buf);
626             sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
627             grpnms[3*i+ZZ] = gmx_strdup(buf);
628         }
629         md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
630         sfree(grpnms);
631     }
632
633     if (fp_ene)
634     {
635         do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
636     }
637
638     md->print_grpnms = nullptr;
639
640     /* check whether we're going to write dh histograms */
641     md->dhc = nullptr;
642     if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
643     {
644         /* Currently dh histograms are only written with dynamics */
645         if (EI_DYNAMICS(ir->eI))
646         {
647             snew(md->dhc, 1);
648
649             mde_delta_h_coll_init(md->dhc, ir);
650         }
651         md->fp_dhdl = nullptr;
652         snew(md->dE, ir->fepvals->n_lambda);
653     }
654     else
655     {
656         md->fp_dhdl = fp_dhdl;
657         snew(md->dE, ir->fepvals->n_lambda);
658     }
659     if (ir->bSimTemp)
660     {
661         int i;
662         snew(md->temperatures, ir->fepvals->n_lambda);
663         for (i = 0; i < ir->fepvals->n_lambda; i++)
664         {
665             md->temperatures[i] = ir->simtempvals->temperatures[i];
666         }
667     }
668     return md;
669 }
670
671 //! Legacy cleanup function
672 void done_mdebin(t_mdebin *mdebin)
673 {
674     sfree(mdebin->igrp);
675     sfree(mdebin->tmp_r);
676     sfree(mdebin->tmp_v);
677     done_ebin(mdebin->ebin);
678     done_mde_delta_h_coll(mdebin->dhc);
679     sfree(mdebin->dE);
680     sfree(mdebin->temperatures);
681     sfree(mdebin);
682 }
683
684 } // namespace
685 } // namespace gmx
686
687 /* print a lambda vector to a string
688    fep = the inputrec's FEP input data
689    i = the index of the lambda vector
690    get_native_lambda = whether to print the native lambda
691    get_names = whether to print the names rather than the values
692    str = the pre-allocated string buffer to print to. */
693 static void print_lambda_vector(t_lambda *fep, int i,
694                                 gmx_bool get_native_lambda, gmx_bool get_names,
695                                 char *str)
696 {
697     int    j, k = 0;
698     int    Nsep = 0;
699
700     for (j = 0; j < efptNR; j++)
701     {
702         if (fep->separate_dvdl[j])
703         {
704             Nsep++;
705         }
706     }
707     str[0] = 0; /* reset the string */
708     if (Nsep > 1)
709     {
710         str += sprintf(str, "("); /* set the opening parenthesis*/
711     }
712     for (j = 0; j < efptNR; j++)
713     {
714         if (fep->separate_dvdl[j])
715         {
716             if (!get_names)
717             {
718                 if (get_native_lambda && fep->init_lambda >= 0)
719                 {
720                     str += sprintf(str, "%.4f", fep->init_lambda);
721                 }
722                 else
723                 {
724                     str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
725                 }
726             }
727             else
728             {
729                 str += sprintf(str, "%s", efpt_singular_names[j]);
730             }
731             /* print comma for the next item */
732             if (k < Nsep-1)
733             {
734                 str += sprintf(str, ", ");
735             }
736             k++;
737         }
738     }
739     if (Nsep > 1)
740     {
741         /* and add the closing parenthesis */
742         sprintf(str, ")");
743     }
744 }
745
746 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
747                 const gmx_output_env_t *oenv)
748 {
749     FILE       *fp;
750     const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
751     *lambdastate     = "\\lambda state";
752     int         i, nsets, nsets_de, nsetsbegin;
753     int         n_lambda_terms = 0;
754     t_lambda   *fep            = ir->fepvals; /* for simplicity */
755     t_expanded *expand         = ir->expandedvals;
756     char        lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
757
758     int         nsets_dhdl = 0;
759     int         s          = 0;
760     int         nsetsextend;
761     gmx_bool    write_pV = FALSE;
762
763     /* count the number of different lambda terms */
764     for (i = 0; i < efptNR; i++)
765     {
766         if (fep->separate_dvdl[i])
767         {
768             n_lambda_terms++;
769         }
770     }
771
772     std::string title, label_x, label_y;
773     if (fep->n_lambda == 0)
774     {
775         title   = gmx::formatString("%s", dhdl);
776         label_x = gmx::formatString("Time (ps)");
777         label_y = gmx::formatString("%s (%s %s)",
778                                     dhdl, unit_energy, "[\\lambda]\\S-1\\N");
779     }
780     else
781     {
782         title   = gmx::formatString("%s and %s", dhdl, deltag);
783         label_x = gmx::formatString("Time (ps)");
784         label_y = gmx::formatString("%s and %s (%s %s)",
785                                     dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
786     }
787     fp = gmx_fio_fopen(filename, "w+");
788     xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
789
790     std::string buf;
791     if (!(ir->bSimTemp))
792     {
793         buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
794     }
795     if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
796     {
797         if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
798         {
799             /* compatibility output */
800             buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
801         }
802         else
803         {
804             print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
805                                 lambda_vec_str);
806             print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
807                                 lambda_name_str);
808             buf += gmx::formatString("%s %d: %s = %s",
809                                      lambdastate, fep->init_fep_state,
810                                      lambda_name_str, lambda_vec_str);
811         }
812     }
813     xvgr_subtitle(fp, buf.c_str(), oenv);
814
815
816     nsets_dhdl = 0;
817     if (fep->dhdl_derivatives == edhdlderivativesYES)
818     {
819         nsets_dhdl = n_lambda_terms;
820     }
821     /* count the number of delta_g states */
822     nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
823
824     nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
825
826     if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
827     {
828         nsets += 1;   /*add fep state for expanded ensemble */
829     }
830
831     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
832     {
833         nsets += 1;  /* add energy to the dhdl as well */
834     }
835
836     nsetsextend = nsets;
837     if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
838     {
839         nsetsextend += 1; /* for PV term, other terms possible if required for
840                              the reduced potential (only needed with foreign
841                              lambda, and only output when init_lambda is not
842                              set in order to maintain compatibility of the
843                              dhdl.xvg file) */
844         write_pV     = TRUE;
845     }
846     std::vector<std::string> setname(nsetsextend);
847
848     if (expand->elmcmove > elmcmoveNO)
849     {
850         /* state for the fep_vals, if we have alchemical sampling */
851         setname[s++] = "Thermodynamic state";
852     }
853
854     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
855     {
856         std::string energy;
857         switch (fep->edHdLPrintEnergy)
858         {
859             case edHdLPrintEnergyPOTENTIAL:
860                 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
861                 break;
862             case edHdLPrintEnergyTOTAL:
863             case edHdLPrintEnergyYES:
864             default:
865                 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
866         }
867         setname[s++] = energy;
868     }
869
870     if (fep->dhdl_derivatives == edhdlderivativesYES)
871     {
872         for (i = 0; i < efptNR; i++)
873         {
874             if (fep->separate_dvdl[i])
875             {
876                 std::string derivative;
877                 if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
878                 {
879                     /* compatibility output */
880                     derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
881                 }
882                 else
883                 {
884                     double lam = fep->init_lambda;
885                     if (fep->init_lambda < 0)
886                     {
887                         lam = fep->all_lambda[i][fep->init_fep_state];
888                     }
889                     derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
890                                                    lam);
891                 }
892                 setname[s++] = derivative;
893             }
894         }
895     }
896
897     if (fep->n_lambda > 0)
898     {
899         /* g_bar has to determine the lambda values used in this simulation
900          * from this xvg legend.
901          */
902
903         if (expand->elmcmove > elmcmoveNO)
904         {
905             nsetsbegin = 1;  /* for including the expanded ensemble */
906         }
907         else
908         {
909             nsetsbegin = 0;
910         }
911
912         if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
913         {
914             nsetsbegin += 1;
915         }
916         nsetsbegin += nsets_dhdl;
917
918         for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
919         {
920             print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
921             std::string buf;
922             if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
923             {
924                 /* for compatible dhdl.xvg files */
925                 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
926             }
927             else
928             {
929                 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
930             }
931
932             if (ir->bSimTemp)
933             {
934                 /* print the temperature for this state if doing simulated annealing */
935                 buf += gmx::formatString("T = %g (%s)",
936                                          ir->simtempvals->temperatures[s-(nsetsbegin)],
937                                          unit_temp_K);
938             }
939             setname[s++] = buf;
940         }
941         if (write_pV)
942         {
943             setname[s++] = gmx::formatString("pV (%s)", unit_energy);
944         }
945
946         xvgrLegend(fp, setname, oenv);
947     }
948
949     return fp;
950 }
951
952 namespace gmx
953 {
954 namespace
955 {
956
957 //! Legacy update function
958 void upd_mdebin(t_mdebin               *md,
959                 gmx_bool                bDoDHDL,
960                 gmx_bool                bSum,
961                 double                  time,
962                 real                    tmass,
963                 gmx_enerdata_t         *enerd,
964                 t_state                *state,
965                 t_lambda               *fep,
966                 t_expanded             *expand,
967                 matrix                  box,
968                 tensor                  svir,
969                 tensor                  fvir,
970                 tensor                  vir,
971                 tensor                  pres,
972                 gmx_ekindata_t         *ekind,
973                 rvec                    mu_tot,
974                 const gmx::Constraints *constr)
975 {
976     int    i, j, k, kk, n, gid;
977     real   crmsd[2], tmp6[6];
978     real   bs[NTRICLBOXS], vol, dens, pv, enthalpy;
979     real   eee[egNR];
980     double store_dhdl[efptNR];
981     real   store_energy = 0;
982     real   tmp;
983
984     /* Do NOT use the box in the state variable, but the separate box provided
985      * as an argument. This is because we sometimes need to write the box from
986      * the last timestep to match the trajectory frames.
987      */
988     add_ebin_indexed(md->ebin, md->ie, gmx::ArrayRef<gmx_bool>(md->bEner), enerd->term, bSum);
989     if (md->nCrmsd)
990     {
991         crmsd[0] = constr->rmsd();
992         add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
993     }
994     if (md->bDynBox)
995     {
996         int nboxs;
997         if (md->bTricl)
998         {
999             bs[0] = box[XX][XX];
1000             bs[1] = box[YY][YY];
1001             bs[2] = box[ZZ][ZZ];
1002             bs[3] = box[YY][XX];
1003             bs[4] = box[ZZ][XX];
1004             bs[5] = box[ZZ][YY];
1005             nboxs = NTRICLBOXS;
1006         }
1007         else
1008         {
1009             bs[0] = box[XX][XX];
1010             bs[1] = box[YY][YY];
1011             bs[2] = box[ZZ][ZZ];
1012             nboxs = NBOXS;
1013         }
1014         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1015         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1016         add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1017         add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1018         add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1019
1020         if (md->bDiagPres)
1021         {
1022             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
1023                not the instantaneous pressure */
1024             pv = vol*md->ref_p/PRESFAC;
1025
1026             add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1027             enthalpy = pv + enerd->term[F_ETOT];
1028             add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1029         }
1030     }
1031     if (md->bConstrVir)
1032     {
1033         add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1034         add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1035     }
1036     if (md->bPres)
1037     {
1038         add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1039         add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1040         tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1041         add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1042     }
1043     if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1044     {
1045         tmp6[0] = state->boxv[XX][XX];
1046         tmp6[1] = state->boxv[YY][YY];
1047         tmp6[2] = state->boxv[ZZ][ZZ];
1048         tmp6[3] = state->boxv[YY][XX];
1049         tmp6[4] = state->boxv[ZZ][XX];
1050         tmp6[5] = state->boxv[ZZ][YY];
1051         add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1052     }
1053     if (md->bMu)
1054     {
1055         add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1056     }
1057     if (ekind && ekind->cosacc.cos_accel != 0)
1058     {
1059         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1060         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1061         add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1062         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1063         tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1064                  *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1065         add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1066     }
1067     if (md->nE > 1)
1068     {
1069         n = 0;
1070         for (i = 0; (i < md->nEg); i++)
1071         {
1072             for (j = i; (j < md->nEg); j++)
1073             {
1074                 gid = GID(i, j, md->nEg);
1075                 for (k = kk = 0; (k < egNR); k++)
1076                 {
1077                     if (md->bEInd[k])
1078                     {
1079                         eee[kk++] = enerd->grpp.ener[k][gid];
1080                     }
1081                 }
1082                 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1083                 n++;
1084             }
1085         }
1086     }
1087
1088     if (ekind)
1089     {
1090         for (i = 0; (i < md->nTC); i++)
1091         {
1092             md->tmp_r[i] = ekind->tcstat[i].T;
1093         }
1094         add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1095
1096         if (md->etc == etcNOSEHOOVER)
1097         {
1098             /* whether to print Nose-Hoover chains: */
1099             if (md->bPrintNHChains)
1100             {
1101                 if (md->bNHC_trotter)
1102                 {
1103                     for (i = 0; (i < md->nTC); i++)
1104                     {
1105                         for (j = 0; j < md->nNHC; j++)
1106                         {
1107                             k                = i*md->nNHC+j;
1108                             md->tmp_r[2*k]   = state->nosehoover_xi[k];
1109                             md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1110                         }
1111                     }
1112                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1113
1114                     if (md->bMTTK)
1115                     {
1116                         for (i = 0; (i < md->nTCP); i++)
1117                         {
1118                             for (j = 0; j < md->nNHC; j++)
1119                             {
1120                                 k                = i*md->nNHC+j;
1121                                 md->tmp_r[2*k]   = state->nhpres_xi[k];
1122                                 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1123                             }
1124                         }
1125                         add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1126                     }
1127                 }
1128                 else
1129                 {
1130                     for (i = 0; (i < md->nTC); i++)
1131                     {
1132                         md->tmp_r[2*i]   = state->nosehoover_xi[i];
1133                         md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1134                     }
1135                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1136                 }
1137             }
1138         }
1139         else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1140                  md->etc == etcVRESCALE)
1141         {
1142             for (i = 0; (i < md->nTC); i++)
1143             {
1144                 md->tmp_r[i] = ekind->tcstat[i].lambda;
1145             }
1146             add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1147         }
1148     }
1149
1150     if (ekind && md->nU > 1)
1151     {
1152         for (i = 0; (i < md->nU); i++)
1153         {
1154             copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1155         }
1156         add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1157     }
1158
1159     ebin_increase_count(md->ebin, bSum);
1160
1161     /* BAR + thermodynamic integration values */
1162     if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1163     {
1164         for (i = 0; i < enerd->n_lambda-1; i++)
1165         {
1166             /* zero for simulated tempering */
1167             md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1168             if (md->temperatures != nullptr)
1169             {
1170                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1171                 /* is this even useful to have at all? */
1172                 md->dE[i] += (md->temperatures[i]/
1173                               md->temperatures[state->fep_state]-1.0)*
1174                     enerd->term[F_EKIN];
1175             }
1176         }
1177
1178         if (md->fp_dhdl)
1179         {
1180             fprintf(md->fp_dhdl, "%.4f", time);
1181             /* the current free energy state */
1182
1183             /* print the current state if we are doing expanded ensemble */
1184             if (expand->elmcmove > elmcmoveNO)
1185             {
1186                 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1187             }
1188             /* total energy (for if the temperature changes */
1189
1190             if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1191             {
1192                 switch (fep->edHdLPrintEnergy)
1193                 {
1194                     case edHdLPrintEnergyPOTENTIAL:
1195                         store_energy = enerd->term[F_EPOT];
1196                         break;
1197                     case edHdLPrintEnergyTOTAL:
1198                     case edHdLPrintEnergyYES:
1199                     default:
1200                         store_energy = enerd->term[F_ETOT];
1201                 }
1202                 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1203             }
1204
1205             if (fep->dhdl_derivatives == edhdlderivativesYES)
1206             {
1207                 for (i = 0; i < efptNR; i++)
1208                 {
1209                     if (fep->separate_dvdl[i])
1210                     {
1211                         /* assumes F_DVDL is first */
1212                         fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1213                     }
1214                 }
1215             }
1216             for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1217             {
1218                 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1219             }
1220             if (md->bDynBox &&
1221                 md->bDiagPres &&
1222                 (md->epc != epcNO) &&
1223                 (enerd->n_lambda > 0) &&
1224                 (fep->init_lambda < 0))
1225             {
1226                 fprintf(md->fp_dhdl, " %#.8g", pv);  /* PV term only needed when
1227                                                         there are alternate state
1228                                                         lambda and we're not in
1229                                                         compatibility mode */
1230             }
1231             fprintf(md->fp_dhdl, "\n");
1232             /* and the binary free energy output */
1233         }
1234         if (md->dhc && bDoDHDL)
1235         {
1236             int idhdl = 0;
1237             for (i = 0; i < efptNR; i++)
1238             {
1239                 if (fep->separate_dvdl[i])
1240                 {
1241                     /* assumes F_DVDL is first */
1242                     store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1243                     idhdl            += 1;
1244                 }
1245             }
1246             store_energy = enerd->term[F_ETOT];
1247             /* store_dh is dE */
1248             mde_delta_h_coll_add_dh(md->dhc,
1249                                     static_cast<double>(state->fep_state),
1250                                     store_energy,
1251                                     pv,
1252                                     store_dhdl,
1253                                     md->dE + fep->lambda_start_n,
1254                                     time);
1255         }
1256     }
1257 }
1258
1259 }   // namespace
1260
1261 void EnergyOutput::recordNonEnergyStep()
1262 {
1263     ebin_increase_count(mdebin->ebin, false);
1264 }
1265
1266 namespace
1267 {
1268
1269 //! Legacy output function
1270 void npr(FILE *log, int n, char c)
1271 {
1272     for (; (n > 0); n--)
1273     {
1274         fprintf(log, "%c", c);
1275     }
1276 }
1277
1278 //! Legacy output function
1279 void pprint(FILE *log, const char *s, t_mdebin *md)
1280 {
1281     char CHAR = '#';
1282     int  slen;
1283     char buf1[22], buf2[22];
1284
1285     slen = strlen(s);
1286     fprintf(log, "\t<======  ");
1287     npr(log, slen, CHAR);
1288     fprintf(log, "  ==>\n");
1289     fprintf(log, "\t<====  %s  ====>\n", s);
1290     fprintf(log, "\t<==  ");
1291     npr(log, slen, CHAR);
1292     fprintf(log, "  ======>\n\n");
1293
1294     fprintf(log, "\tStatistics over %s steps using %s frames\n",
1295             gmx_step_str(md->ebin->nsteps_sim, buf1),
1296             gmx_step_str(md->ebin->nsum_sim, buf2));
1297     fprintf(log, "\n");
1298 }
1299
1300 }   // namespace
1301
1302 void print_ebin_header(FILE *log, int64_t steps, double time)
1303 {
1304     char buf[22];
1305
1306     fprintf(log, "   %12s   %12s\n"
1307             "   %12s   %12.5f\n\n",
1308             "Step", "Time", gmx_step_str(steps, buf), time);
1309 }
1310
1311 namespace
1312 {
1313
1314 // TODO It is too many responsibilities for this function to handle
1315 // both .edr and .log output for both per-time and time-average data.
1316 //! Legacy ebin output function
1317 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1318                 FILE *log,
1319                 int64_t step, double time,
1320                 int mode,
1321                 t_mdebin *md, t_fcdata *fcd,
1322                 gmx_groups_t *groups, t_grpopts *opts,
1323                 gmx::Awh *awh)
1324 {
1325     /*static char **grpnms=NULL;*/
1326     char         buf[246];
1327     int          i, j, n, ni, nj, b;
1328     int          ndisre = 0;
1329
1330     /* these are for the old-style blocks (1 subblock, only reals), because
1331        there can be only one per ID for these */
1332     int          nr[enxNR];
1333     int          id[enxNR];
1334     real        *block[enxNR];
1335
1336     t_enxframe   fr;
1337
1338     if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1339     {
1340         if (log)
1341         {
1342             fprintf(log, "Not enough data recorded to report energy averages\n");
1343         }
1344         return;
1345     }
1346
1347     switch (mode)
1348     {
1349         case eprNORMAL:
1350             init_enxframe(&fr);
1351             fr.t            = time;
1352             fr.step         = step;
1353             fr.nsteps       = md->ebin->nsteps;
1354             fr.dt           = md->delta_t;
1355             fr.nsum         = md->ebin->nsum;
1356             fr.nre          = (bEne) ? md->ebin->nener : 0;
1357             fr.ener         = md->ebin->e;
1358             ndisre          = bDR ? fcd->disres.npair : 0;
1359             /* Optional additional old-style (real-only) blocks. */
1360             for (i = 0; i < enxNR; i++)
1361             {
1362                 nr[i] = 0;
1363             }
1364             if (bOR && fcd->orires.nr > 0)
1365             {
1366                 diagonalize_orires_tensors(&(fcd->orires));
1367                 nr[enxOR]     = fcd->orires.nr;
1368                 block[enxOR]  = fcd->orires.otav;
1369                 id[enxOR]     = enxOR;
1370                 nr[enxORI]    = (fcd->orires.oinsl != fcd->orires.otav) ?
1371                     fcd->orires.nr : 0;
1372                 block[enxORI] = fcd->orires.oinsl;
1373                 id[enxORI]    = enxORI;
1374                 nr[enxORT]    = fcd->orires.nex*12;
1375                 block[enxORT] = fcd->orires.eig;
1376                 id[enxORT]    = enxORT;
1377             }
1378
1379             /* whether we are going to wrte anything out: */
1380             if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1381             {
1382
1383                 /* the old-style blocks go first */
1384                 fr.nblock = 0;
1385                 for (i = 0; i < enxNR; i++)
1386                 {
1387                     if (nr[i] > 0)
1388                     {
1389                         fr.nblock = i + 1;
1390                     }
1391                 }
1392                 add_blocks_enxframe(&fr, fr.nblock);
1393                 for (b = 0; b < fr.nblock; b++)
1394                 {
1395                     add_subblocks_enxblock(&(fr.block[b]), 1);
1396                     fr.block[b].id        = id[b];
1397                     fr.block[b].sub[0].nr = nr[b];
1398 #if !GMX_DOUBLE
1399                     fr.block[b].sub[0].type = xdr_datatype_float;
1400                     fr.block[b].sub[0].fval = block[b];
1401 #else
1402                     fr.block[b].sub[0].type = xdr_datatype_double;
1403                     fr.block[b].sub[0].dval = block[b];
1404 #endif
1405                 }
1406
1407                 /* check for disre block & fill it. */
1408                 if (ndisre > 0)
1409                 {
1410                     int db = fr.nblock;
1411                     fr.nblock += 1;
1412                     add_blocks_enxframe(&fr, fr.nblock);
1413
1414                     add_subblocks_enxblock(&(fr.block[db]), 2);
1415                     fr.block[db].id        = enxDISRE;
1416                     fr.block[db].sub[0].nr = ndisre;
1417                     fr.block[db].sub[1].nr = ndisre;
1418 #if !GMX_DOUBLE
1419                     fr.block[db].sub[0].type = xdr_datatype_float;
1420                     fr.block[db].sub[1].type = xdr_datatype_float;
1421                     fr.block[db].sub[0].fval = fcd->disres.rt;
1422                     fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1423 #else
1424                     fr.block[db].sub[0].type = xdr_datatype_double;
1425                     fr.block[db].sub[1].type = xdr_datatype_double;
1426                     fr.block[db].sub[0].dval = fcd->disres.rt;
1427                     fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1428 #endif
1429                 }
1430                 /* here we can put new-style blocks */
1431
1432                 /* Free energy perturbation blocks */
1433                 if (md->dhc)
1434                 {
1435                     mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1436                 }
1437
1438                 /* we can now free & reset the data in the blocks */
1439                 if (md->dhc)
1440                 {
1441                     mde_delta_h_coll_reset(md->dhc);
1442                 }
1443
1444                 /* AWH bias blocks. */
1445                 if (awh != nullptr)  // TODO: add boolean in t_mdebin. See in mdebin.h.
1446                 {
1447                     awh->writeToEnergyFrame(step, &fr);
1448                 }
1449
1450                 /* do the actual I/O */
1451                 do_enx(fp_ene, &fr);
1452                 if (fr.nre)
1453                 {
1454                     /* We have stored the sums, so reset the sum history */
1455                     reset_ebin_sums(md->ebin);
1456                 }
1457             }
1458             free_enxframe(&fr);
1459             break;
1460         case eprAVER:
1461             if (log)
1462             {
1463                 pprint(log, "A V E R A G E S", md);
1464             }
1465             break;
1466         case eprRMS:
1467             if (log)
1468             {
1469                 pprint(log, "R M S - F L U C T U A T I O N S", md);
1470             }
1471             break;
1472         default:
1473             gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1474     }
1475
1476     if (log)
1477     {
1478         if (opts)
1479         {
1480             for (i = 0; i < opts->ngtc; i++)
1481             {
1482                 if (opts->annealing[i] != eannNO)
1483                 {
1484                     fprintf(log, "Current ref_t for group %s: %8.1f\n",
1485                             *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1486                             opts->ref_t[i]);
1487                 }
1488             }
1489         }
1490         if (mode == eprNORMAL && bOR && fcd->orires.nr > 0)
1491         {
1492             print_orires_log(log, &(fcd->orires));
1493         }
1494         fprintf(log, "   Energies (%s)\n", unit_energy);
1495         pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1496         fprintf(log, "\n");
1497
1498         if (mode == eprAVER)
1499         {
1500             if (md->bDynBox)
1501             {
1502                 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1503                         mode, TRUE);
1504                 fprintf(log, "\n");
1505             }
1506             if (md->bConstrVir)
1507             {
1508                 fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1509                 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1510                 fprintf(log, "\n");
1511                 fprintf(log, "   Force Virial (%s)\n", unit_energy);
1512                 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1513                 fprintf(log, "\n");
1514             }
1515             if (md->bPres)
1516             {
1517                 fprintf(log, "   Total Virial (%s)\n", unit_energy);
1518                 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1519                 fprintf(log, "\n");
1520                 fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1521                 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1522                 fprintf(log, "\n");
1523             }
1524             if (md->bMu)
1525             {
1526                 fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1527                 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1528                 fprintf(log, "\n");
1529             }
1530
1531             if (md->nE > 1)
1532             {
1533                 if (md->print_grpnms == nullptr)
1534                 {
1535                     snew(md->print_grpnms, md->nE);
1536                     n = 0;
1537                     for (i = 0; (i < md->nEg); i++)
1538                     {
1539                         ni = groups->grps[egcENER].nm_ind[i];
1540                         for (j = i; (j < md->nEg); j++)
1541                         {
1542                             nj = groups->grps[egcENER].nm_ind[j];
1543                             sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1544                                     *(groups->grpname[nj]));
1545                             md->print_grpnms[n++] = gmx_strdup(buf);
1546                         }
1547                     }
1548                 }
1549                 sprintf(buf, "Epot (%s)", unit_energy);
1550                 fprintf(log, "%15s   ", buf);
1551                 for (i = 0; (i < egNR); i++)
1552                 {
1553                     if (md->bEInd[i])
1554                     {
1555                         fprintf(log, "%12s   ", egrp_nm[i]);
1556                     }
1557                 }
1558                 fprintf(log, "\n");
1559                 for (i = 0; (i < md->nE); i++)
1560                 {
1561                     fprintf(log, "%15s", md->print_grpnms[i]);
1562                     pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1563                             FALSE);
1564                 }
1565                 fprintf(log, "\n");
1566             }
1567             if (md->nTC > 1)
1568             {
1569                 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1570                 fprintf(log, "\n");
1571             }
1572             if (md->nU > 1)
1573             {
1574                 fprintf(log, "%15s   %12s   %12s   %12s\n",
1575                         "Group", "Ux", "Uy", "Uz");
1576                 for (i = 0; (i < md->nU); i++)
1577                 {
1578                     ni = groups->grps[egcACC].nm_ind[i];
1579                     fprintf(log, "%15s", *groups->grpname[ni]);
1580                     pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1581                 }
1582                 fprintf(log, "\n");
1583             }
1584         }
1585     }
1586
1587 }
1588
1589 //! Legacy update function
1590 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1591 {
1592     const t_ebin * const ebin = mdebin->ebin;
1593
1594     enerhist->nsteps     = ebin->nsteps;
1595     enerhist->nsum       = ebin->nsum;
1596     enerhist->nsteps_sim = ebin->nsteps_sim;
1597     enerhist->nsum_sim   = ebin->nsum_sim;
1598
1599     if (ebin->nsum > 0)
1600     {
1601         /* This will only actually resize the first time */
1602         enerhist->ener_ave.resize(ebin->nener);
1603         enerhist->ener_sum.resize(ebin->nener);
1604
1605         for (int i = 0; i < ebin->nener; i++)
1606         {
1607             enerhist->ener_ave[i] = ebin->e[i].eav;
1608             enerhist->ener_sum[i] = ebin->e[i].esum;
1609         }
1610     }
1611
1612     if (ebin->nsum_sim > 0)
1613     {
1614         /* This will only actually resize the first time */
1615         enerhist->ener_sum_sim.resize(ebin->nener);
1616
1617         for (int i = 0; i < ebin->nener; i++)
1618         {
1619             enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1620         }
1621     }
1622     if (mdebin->dhc)
1623     {
1624         mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1625     }
1626 }
1627
1628 //! Legacy restore function
1629 void restore_energyhistory_from_state(t_mdebin              * mdebin,
1630                                       const energyhistory_t * enerhist)
1631 {
1632     unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1633
1634     GMX_RELEASE_ASSERT(enerhist, "Need valid history to restore");
1635
1636     if ((enerhist->nsum     > 0 && nener != enerhist->ener_sum.size()) ||
1637         (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1638     {
1639         gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1640                   nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1641     }
1642
1643     mdebin->ebin->nsteps     = enerhist->nsteps;
1644     mdebin->ebin->nsum       = enerhist->nsum;
1645     mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1646     mdebin->ebin->nsum_sim   = enerhist->nsum_sim;
1647
1648     for (int i = 0; i < mdebin->ebin->nener; i++)
1649     {
1650         mdebin->ebin->e[i].eav  =
1651             (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1652         mdebin->ebin->e[i].esum =
1653             (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1654         mdebin->ebin->e_sim[i].esum =
1655             (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1656     }
1657     if (mdebin->dhc)
1658     {
1659         mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());
1660     }
1661 }
1662
1663 }   // namespace
1664
1665 EnergyOutput::EnergyOutput()
1666     : mdebin(nullptr)
1667 {
1668 }
1669
1670 void EnergyOutput::prepare(ener_file        *fp_ene,
1671                            const gmx_mtop_t *mtop,
1672                            const t_inputrec *ir,
1673                            FILE             *fp_dhdl,
1674                            bool              isRerun)
1675 {
1676     mdebin = init_mdebin(fp_ene, mtop, ir, fp_dhdl, isRerun);
1677 }
1678
1679 EnergyOutput::~EnergyOutput()
1680 {
1681     done_mdebin(mdebin);
1682 }
1683
1684 t_ebin *EnergyOutput::getEbin()
1685 {
1686     return mdebin->ebin;
1687 }
1688
1689 void EnergyOutput::addDataAtEnergyStep(bool                    bDoDHDL,
1690                                        bool                    bSum,
1691                                        double                  time,
1692                                        real                    tmass,
1693                                        gmx_enerdata_t         *enerd,
1694                                        t_state                *state,
1695                                        t_lambda               *fep,
1696                                        t_expanded             *expand,
1697                                        matrix                  box,
1698                                        tensor                  svir,
1699                                        tensor                  fvir,
1700                                        tensor                  vir,
1701                                        tensor                  pres,
1702                                        gmx_ekindata_t         *ekind,
1703                                        rvec                    mu_tot,
1704                                        const gmx::Constraints *constr)
1705 {
1706     upd_mdebin(mdebin, bDoDHDL, bSum, time, tmass, enerd, state, fep,
1707                expand, box, svir, fvir, vir, pres, ekind, mu_tot, constr);
1708 }
1709
1710 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1711                                          FILE *log,
1712                                          int64_t step, double time,
1713                                          int mode,
1714                                          t_fcdata *fcd,
1715                                          gmx_groups_t *groups, t_grpopts *opts,
1716                                          gmx::Awh *awh)
1717 {
1718     print_ebin(fp_ene, bEne, bDR, bOR, log, step, time, mode,
1719                mdebin, fcd, groups, opts, awh);
1720 }
1721
1722 int EnergyOutput::numEnergyTerms() const
1723 {
1724     return mdebin->ebin->nener;
1725 }
1726
1727 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1728 {
1729     update_energyhistory(enerhist, mdebin);
1730 }
1731
1732 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1733 {
1734     restore_energyhistory_from_state(mdebin, &enerhist);
1735 }
1736
1737 } // namespace gmx