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