Use test with parameters for EnergyOutput
[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                       const pull_t     *pull_work,
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 SimulationGroups    *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(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     n       = groups->groups[SimulationAtomGroupType::EnergyOutput].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->groups[SimulationAtomGroupType::EnergyOutput].nr); i++)
467         {
468             ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
469             for (j = i; (j < groups->groups[SimulationAtomGroupType::EnergyOutput].nr); j++)
470             {
471                 nj = groups->groups[SimulationAtomGroupType::EnergyOutput].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->groupNames[ni]), *(groups->groupNames[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->groups[SimulationAtomGroupType::TemperatureCoupling].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
532     // TODO redo the group name memory management to make it more clear
533     char **grpnms;
534     snew(grpnms, std::max(md->mde_n, md->mdeb_n)); // Just in case md->mdeb_n > md->mde_n
535
536     for (i = 0; (i < md->nTC); i++)
537     {
538         ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
539         sprintf(buf, "T-%s", *(groups->groupNames[ni]));
540         grpnms[i] = gmx_strdup(buf);
541     }
542     md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
543                                unit_temp_K);
544     for (i = 0; i < md->nTC; i++)
545     {
546         sfree(grpnms[i]);
547     }
548
549     int allocated = 0;
550     if (md->etc == etcNOSEHOOVER)
551     {
552         if (md->bPrintNHChains)
553         {
554             if (md->bNHC_trotter)
555             {
556                 for (i = 0; (i < md->nTC); i++)
557                 {
558                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
559                     bufi = *(groups->groupNames[ni]);
560                     for (j = 0; (j < md->nNHC); j++)
561                     {
562                         sprintf(buf, "Xi-%d-%s", j, bufi);
563                         grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
564                         sprintf(buf, "vXi-%d-%s", j, bufi);
565                         grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
566                     }
567                 }
568                 md->itc = get_ebin_space(md->ebin, md->mde_n,
569                                          grpnms, unit_invtime);
570                 allocated = md->mde_n;
571                 if (md->bMTTK)
572                 {
573                     for (i = 0; (i < md->nTCP); i++)
574                     {
575                         bufi = baro_nm[0];  /* All barostat DOF's together for now. */
576                         for (j = 0; (j < md->nNHC); j++)
577                         {
578                             sprintf(buf, "Xi-%d-%s", j, bufi);
579                             grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
580                             sprintf(buf, "vXi-%d-%s", j, bufi);
581                             grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
582                         }
583                     }
584                     md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
585                                               grpnms, unit_invtime);
586                     allocated = md->mdeb_n;
587                 }
588             }
589             else
590             {
591                 for (i = 0; (i < md->nTC); i++)
592                 {
593                     ni   = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
594                     bufi = *(groups->groupNames[ni]);
595                     sprintf(buf, "Xi-%s", bufi);
596                     grpnms[2*i] = gmx_strdup(buf);
597                     sprintf(buf, "vXi-%s", bufi);
598                     grpnms[2*i+1] = gmx_strdup(buf);
599                 }
600                 md->itc = get_ebin_space(md->ebin, md->mde_n,
601                                          grpnms, unit_invtime);
602                 allocated = md->mde_n;
603             }
604         }
605     }
606     else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
607              md->etc == etcVRESCALE)
608     {
609         for (i = 0; (i < md->nTC); i++)
610         {
611             ni = groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i];
612             sprintf(buf, "Lamb-%s", *(groups->groupNames[ni]));
613             grpnms[i] = gmx_strdup(buf);
614         }
615         md->itc   = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
616         allocated = md->mde_n;
617     }
618
619     for (i = 0; i < allocated; i++)
620     {
621         sfree(grpnms[i]);
622     }
623     sfree(grpnms);
624
625     md->nU = groups->groups[SimulationAtomGroupType::Acceleration].nr;
626     snew(md->tmp_v, md->nU);
627     if (md->nU > 1)
628     {
629         snew(grpnms, 3*md->nU);
630         for (i = 0; (i < md->nU); i++)
631         {
632             ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
633             sprintf(buf, "Ux-%s", *(groups->groupNames[ni]));
634             grpnms[3*i+XX] = gmx_strdup(buf);
635             sprintf(buf, "Uy-%s", *(groups->groupNames[ni]));
636             grpnms[3*i+YY] = gmx_strdup(buf);
637             sprintf(buf, "Uz-%s", *(groups->groupNames[ni]));
638             grpnms[3*i+ZZ] = gmx_strdup(buf);
639         }
640         md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
641         for (i = 0; i < 3*md->nU; i++)
642         {
643             sfree(grpnms[i]);
644         }
645         sfree(grpnms);
646     }
647
648     if (fp_ene)
649     {
650         do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
651     }
652
653     md->print_grpnms = nullptr;
654
655     /* check whether we're going to write dh histograms */
656     md->dhc = nullptr;
657     if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
658     {
659         /* Currently dh histograms are only written with dynamics */
660         if (EI_DYNAMICS(ir->eI))
661         {
662             snew(md->dhc, 1);
663
664             mde_delta_h_coll_init(md->dhc, ir);
665         }
666         md->fp_dhdl = nullptr;
667         snew(md->dE, ir->fepvals->n_lambda);
668     }
669     else
670     {
671         md->fp_dhdl = fp_dhdl;
672         snew(md->dE, ir->fepvals->n_lambda);
673     }
674     if (ir->bSimTemp)
675     {
676         int i;
677         snew(md->temperatures, ir->fepvals->n_lambda);
678         for (i = 0; i < ir->fepvals->n_lambda; i++)
679         {
680             md->temperatures[i] = ir->simtempvals->temperatures[i];
681         }
682     }
683     return md;
684 }
685
686 //! Legacy cleanup function
687 void done_mdebin(t_mdebin *mdebin)
688 {
689     sfree(mdebin->igrp);
690     sfree(mdebin->tmp_r);
691     sfree(mdebin->tmp_v);
692     done_ebin(mdebin->ebin);
693     done_mde_delta_h_coll(mdebin->dhc);
694     sfree(mdebin->dE);
695     sfree(mdebin->temperatures);
696     if (mdebin->nE > 1 && mdebin->print_grpnms != nullptr)
697     {
698         for (int n = 0; n < mdebin->nE; n++)
699         {
700             sfree(mdebin->print_grpnms[n]);
701         }
702         sfree(mdebin->print_grpnms);
703     }
704     sfree(mdebin);
705 }
706
707 } // namespace
708 } // namespace gmx
709
710 /* print a lambda vector to a string
711    fep = the inputrec's FEP input data
712    i = the index of the lambda vector
713    get_native_lambda = whether to print the native lambda
714    get_names = whether to print the names rather than the values
715    str = the pre-allocated string buffer to print to. */
716 static void print_lambda_vector(t_lambda *fep, int i,
717                                 gmx_bool get_native_lambda, gmx_bool get_names,
718                                 char *str)
719 {
720     int    j, k = 0;
721     int    Nsep = 0;
722
723     for (j = 0; j < efptNR; j++)
724     {
725         if (fep->separate_dvdl[j])
726         {
727             Nsep++;
728         }
729     }
730     str[0] = 0; /* reset the string */
731     if (Nsep > 1)
732     {
733         str += sprintf(str, "("); /* set the opening parenthesis*/
734     }
735     for (j = 0; j < efptNR; j++)
736     {
737         if (fep->separate_dvdl[j])
738         {
739             if (!get_names)
740             {
741                 if (get_native_lambda && fep->init_lambda >= 0)
742                 {
743                     str += sprintf(str, "%.4f", fep->init_lambda);
744                 }
745                 else
746                 {
747                     str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
748                 }
749             }
750             else
751             {
752                 str += sprintf(str, "%s", efpt_singular_names[j]);
753             }
754             /* print comma for the next item */
755             if (k < Nsep-1)
756             {
757                 str += sprintf(str, ", ");
758             }
759             k++;
760         }
761     }
762     if (Nsep > 1)
763     {
764         /* and add the closing parenthesis */
765         sprintf(str, ")");
766     }
767 }
768
769 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
770                 const gmx_output_env_t *oenv)
771 {
772     FILE       *fp;
773     const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
774     *lambdastate     = "\\lambda state";
775     int         i, nsets, nsets_de, nsetsbegin;
776     int         n_lambda_terms = 0;
777     t_lambda   *fep            = ir->fepvals; /* for simplicity */
778     t_expanded *expand         = ir->expandedvals;
779     char        lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
780
781     int         nsets_dhdl = 0;
782     int         s          = 0;
783     int         nsetsextend;
784     gmx_bool    write_pV = FALSE;
785
786     /* count the number of different lambda terms */
787     for (i = 0; i < efptNR; i++)
788     {
789         if (fep->separate_dvdl[i])
790         {
791             n_lambda_terms++;
792         }
793     }
794
795     std::string title, label_x, label_y;
796     if (fep->n_lambda == 0)
797     {
798         title   = gmx::formatString("%s", dhdl);
799         label_x = gmx::formatString("Time (ps)");
800         label_y = gmx::formatString("%s (%s %s)",
801                                     dhdl, unit_energy, "[\\lambda]\\S-1\\N");
802     }
803     else
804     {
805         title   = gmx::formatString("%s and %s", dhdl, deltag);
806         label_x = gmx::formatString("Time (ps)");
807         label_y = gmx::formatString("%s and %s (%s %s)",
808                                     dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
809     }
810     fp = gmx_fio_fopen(filename, "w+");
811     xvgr_header(fp, title.c_str(), label_x, label_y, exvggtXNY, oenv);
812
813     std::string buf;
814     if (!(ir->bSimTemp))
815     {
816         buf = gmx::formatString("T = %g (K) ", ir->opts.ref_t[0]);
817     }
818     if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
819     {
820         if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
821         {
822             /* compatibility output */
823             buf += gmx::formatString("%s = %.4f", lambda, fep->init_lambda);
824         }
825         else
826         {
827             print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
828                                 lambda_vec_str);
829             print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
830                                 lambda_name_str);
831             buf += gmx::formatString("%s %d: %s = %s",
832                                      lambdastate, fep->init_fep_state,
833                                      lambda_name_str, lambda_vec_str);
834         }
835     }
836     xvgr_subtitle(fp, buf.c_str(), oenv);
837
838
839     nsets_dhdl = 0;
840     if (fep->dhdl_derivatives == edhdlderivativesYES)
841     {
842         nsets_dhdl = n_lambda_terms;
843     }
844     /* count the number of delta_g states */
845     nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
846
847     nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
848
849     if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
850     {
851         nsets += 1;   /*add fep state for expanded ensemble */
852     }
853
854     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
855     {
856         nsets += 1;  /* add energy to the dhdl as well */
857     }
858
859     nsetsextend = nsets;
860     if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
861     {
862         nsetsextend += 1; /* for PV term, other terms possible if required for
863                              the reduced potential (only needed with foreign
864                              lambda, and only output when init_lambda is not
865                              set in order to maintain compatibility of the
866                              dhdl.xvg file) */
867         write_pV     = TRUE;
868     }
869     std::vector<std::string> setname(nsetsextend);
870
871     if (expand->elmcmove > elmcmoveNO)
872     {
873         /* state for the fep_vals, if we have alchemical sampling */
874         setname[s++] = "Thermodynamic state";
875     }
876
877     if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
878     {
879         std::string energy;
880         switch (fep->edHdLPrintEnergy)
881         {
882             case edHdLPrintEnergyPOTENTIAL:
883                 energy = gmx::formatString("%s (%s)", "Potential Energy", unit_energy);
884                 break;
885             case edHdLPrintEnergyTOTAL:
886             case edHdLPrintEnergyYES:
887             default:
888                 energy = gmx::formatString("%s (%s)", "Total Energy", unit_energy);
889         }
890         setname[s++] = energy;
891     }
892
893     if (fep->dhdl_derivatives == edhdlderivativesYES)
894     {
895         for (i = 0; i < efptNR; i++)
896         {
897             if (fep->separate_dvdl[i])
898             {
899                 std::string derivative;
900                 if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
901                 {
902                     /* compatibility output */
903                     derivative = gmx::formatString("%s %s %.4f", dhdl, lambda, fep->init_lambda);
904                 }
905                 else
906                 {
907                     double lam = fep->init_lambda;
908                     if (fep->init_lambda < 0)
909                     {
910                         lam = fep->all_lambda[i][fep->init_fep_state];
911                     }
912                     derivative = gmx::formatString("%s %s = %.4f", dhdl, efpt_singular_names[i],
913                                                    lam);
914                 }
915                 setname[s++] = derivative;
916             }
917         }
918     }
919
920     if (fep->n_lambda > 0)
921     {
922         /* g_bar has to determine the lambda values used in this simulation
923          * from this xvg legend.
924          */
925
926         if (expand->elmcmove > elmcmoveNO)
927         {
928             nsetsbegin = 1;  /* for including the expanded ensemble */
929         }
930         else
931         {
932             nsetsbegin = 0;
933         }
934
935         if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
936         {
937             nsetsbegin += 1;
938         }
939         nsetsbegin += nsets_dhdl;
940
941         for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
942         {
943             print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
944             std::string buf;
945             if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
946             {
947                 /* for compatible dhdl.xvg files */
948                 buf = gmx::formatString("%s %s %s", deltag, lambda, lambda_vec_str);
949             }
950             else
951             {
952                 buf = gmx::formatString("%s %s to %s", deltag, lambda, lambda_vec_str);
953             }
954
955             if (ir->bSimTemp)
956             {
957                 /* print the temperature for this state if doing simulated annealing */
958                 buf += gmx::formatString("T = %g (%s)",
959                                          ir->simtempvals->temperatures[s-(nsetsbegin)],
960                                          unit_temp_K);
961             }
962             setname[s++] = buf;
963         }
964         if (write_pV)
965         {
966             setname[s++] = gmx::formatString("pV (%s)", unit_energy);
967         }
968
969         xvgrLegend(fp, setname, oenv);
970     }
971
972     return fp;
973 }
974
975 namespace gmx
976 {
977 namespace
978 {
979
980 //! Legacy update function
981 void upd_mdebin(t_mdebin               *md,
982                 gmx_bool                bDoDHDL,
983                 gmx_bool                bSum,
984                 double                  time,
985                 real                    tmass,
986                 gmx_enerdata_t         *enerd,
987                 t_state                *state,
988                 t_lambda               *fep,
989                 t_expanded             *expand,
990                 matrix                  box,
991                 tensor                  svir,
992                 tensor                  fvir,
993                 tensor                  vir,
994                 tensor                  pres,
995                 gmx_ekindata_t         *ekind,
996                 rvec                    mu_tot,
997                 const gmx::Constraints *constr)
998 {
999     int    i, j, k, kk, n, gid;
1000     real   crmsd[2], tmp6[6];
1001     real   bs[NTRICLBOXS], vol, dens, pv, enthalpy;
1002     real   eee[egNR];
1003     double store_dhdl[efptNR];
1004     real   store_energy = 0;
1005     real   tmp;
1006
1007     /* Do NOT use the box in the state variable, but the separate box provided
1008      * as an argument. This is because we sometimes need to write the box from
1009      * the last timestep to match the trajectory frames.
1010      */
1011     add_ebin_indexed(md->ebin, md->ie, gmx::ArrayRef<gmx_bool>(md->bEner), enerd->term, bSum);
1012     if (md->nCrmsd)
1013     {
1014         crmsd[0] = constr->rmsd();
1015         add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1016     }
1017     if (md->bDynBox)
1018     {
1019         int nboxs;
1020         if (md->bTricl)
1021         {
1022             bs[0] = box[XX][XX];
1023             bs[1] = box[YY][YY];
1024             bs[2] = box[ZZ][ZZ];
1025             bs[3] = box[YY][XX];
1026             bs[4] = box[ZZ][XX];
1027             bs[5] = box[ZZ][YY];
1028             nboxs = NTRICLBOXS;
1029         }
1030         else
1031         {
1032             bs[0] = box[XX][XX];
1033             bs[1] = box[YY][YY];
1034             bs[2] = box[ZZ][ZZ];
1035             nboxs = NBOXS;
1036         }
1037         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1038         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1039         add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1040         add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1041         add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1042
1043         if (md->bDiagPres)
1044         {
1045             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
1046                not the instantaneous pressure */
1047             pv = vol*md->ref_p/PRESFAC;
1048
1049             add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1050             enthalpy = pv + enerd->term[F_ETOT];
1051             add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1052         }
1053     }
1054     if (md->bConstrVir)
1055     {
1056         add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1057         add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1058     }
1059     if (md->bPres)
1060     {
1061         add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1062         add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1063         tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1064         add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1065     }
1066     if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1067     {
1068         tmp6[0] = state->boxv[XX][XX];
1069         tmp6[1] = state->boxv[YY][YY];
1070         tmp6[2] = state->boxv[ZZ][ZZ];
1071         tmp6[3] = state->boxv[YY][XX];
1072         tmp6[4] = state->boxv[ZZ][XX];
1073         tmp6[5] = state->boxv[ZZ][YY];
1074         add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1075     }
1076     if (md->bMu)
1077     {
1078         add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1079     }
1080     if (ekind && ekind->cosacc.cos_accel != 0)
1081     {
1082         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1083         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1084         add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1085         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1086         tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1087                  *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1088         add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1089     }
1090     if (md->nE > 1)
1091     {
1092         n = 0;
1093         for (i = 0; (i < md->nEg); i++)
1094         {
1095             for (j = i; (j < md->nEg); j++)
1096             {
1097                 gid = GID(i, j, md->nEg);
1098                 for (k = kk = 0; (k < egNR); k++)
1099                 {
1100                     if (md->bEInd[k])
1101                     {
1102                         eee[kk++] = enerd->grpp.ener[k][gid];
1103                     }
1104                 }
1105                 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1106                 n++;
1107             }
1108         }
1109     }
1110
1111     if (ekind)
1112     {
1113         for (i = 0; (i < md->nTC); i++)
1114         {
1115             md->tmp_r[i] = ekind->tcstat[i].T;
1116         }
1117         add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1118
1119         if (md->etc == etcNOSEHOOVER)
1120         {
1121             /* whether to print Nose-Hoover chains: */
1122             if (md->bPrintNHChains)
1123             {
1124                 if (md->bNHC_trotter)
1125                 {
1126                     for (i = 0; (i < md->nTC); i++)
1127                     {
1128                         for (j = 0; j < md->nNHC; j++)
1129                         {
1130                             k                = i*md->nNHC+j;
1131                             md->tmp_r[2*k]   = state->nosehoover_xi[k];
1132                             md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1133                         }
1134                     }
1135                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1136
1137                     if (md->bMTTK)
1138                     {
1139                         for (i = 0; (i < md->nTCP); i++)
1140                         {
1141                             for (j = 0; j < md->nNHC; j++)
1142                             {
1143                                 k                = i*md->nNHC+j;
1144                                 md->tmp_r[2*k]   = state->nhpres_xi[k];
1145                                 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1146                             }
1147                         }
1148                         add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1149                     }
1150                 }
1151                 else
1152                 {
1153                     for (i = 0; (i < md->nTC); i++)
1154                     {
1155                         md->tmp_r[2*i]   = state->nosehoover_xi[i];
1156                         md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1157                     }
1158                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1159                 }
1160             }
1161         }
1162         else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1163                  md->etc == etcVRESCALE)
1164         {
1165             for (i = 0; (i < md->nTC); i++)
1166             {
1167                 md->tmp_r[i] = ekind->tcstat[i].lambda;
1168             }
1169             add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1170         }
1171     }
1172
1173     if (ekind && md->nU > 1)
1174     {
1175         for (i = 0; (i < md->nU); i++)
1176         {
1177             copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1178         }
1179         add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1180     }
1181
1182     ebin_increase_count(md->ebin, bSum);
1183
1184     /* BAR + thermodynamic integration values */
1185     if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1186     {
1187         for (i = 0; i < enerd->n_lambda-1; i++)
1188         {
1189             /* zero for simulated tempering */
1190             md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1191             if (md->temperatures != nullptr)
1192             {
1193                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1194                 /* is this even useful to have at all? */
1195                 md->dE[i] += (md->temperatures[i]/
1196                               md->temperatures[state->fep_state]-1.0)*
1197                     enerd->term[F_EKIN];
1198             }
1199         }
1200
1201         if (md->fp_dhdl)
1202         {
1203             fprintf(md->fp_dhdl, "%.4f", time);
1204             /* the current free energy state */
1205
1206             /* print the current state if we are doing expanded ensemble */
1207             if (expand->elmcmove > elmcmoveNO)
1208             {
1209                 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1210             }
1211             /* total energy (for if the temperature changes */
1212
1213             if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1214             {
1215                 switch (fep->edHdLPrintEnergy)
1216                 {
1217                     case edHdLPrintEnergyPOTENTIAL:
1218                         store_energy = enerd->term[F_EPOT];
1219                         break;
1220                     case edHdLPrintEnergyTOTAL:
1221                     case edHdLPrintEnergyYES:
1222                     default:
1223                         store_energy = enerd->term[F_ETOT];
1224                 }
1225                 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1226             }
1227
1228             if (fep->dhdl_derivatives == edhdlderivativesYES)
1229             {
1230                 for (i = 0; i < efptNR; i++)
1231                 {
1232                     if (fep->separate_dvdl[i])
1233                     {
1234                         /* assumes F_DVDL is first */
1235                         fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1236                     }
1237                 }
1238             }
1239             for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1240             {
1241                 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1242             }
1243             if (md->bDynBox &&
1244                 md->bDiagPres &&
1245                 (md->epc != epcNO) &&
1246                 (enerd->n_lambda > 0) &&
1247                 (fep->init_lambda < 0))
1248             {
1249                 fprintf(md->fp_dhdl, " %#.8g", pv);  /* PV term only needed when
1250                                                         there are alternate state
1251                                                         lambda and we're not in
1252                                                         compatibility mode */
1253             }
1254             fprintf(md->fp_dhdl, "\n");
1255             /* and the binary free energy output */
1256         }
1257         if (md->dhc && bDoDHDL)
1258         {
1259             int idhdl = 0;
1260             for (i = 0; i < efptNR; i++)
1261             {
1262                 if (fep->separate_dvdl[i])
1263                 {
1264                     /* assumes F_DVDL is first */
1265                     store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1266                     idhdl            += 1;
1267                 }
1268             }
1269             store_energy = enerd->term[F_ETOT];
1270             /* store_dh is dE */
1271             mde_delta_h_coll_add_dh(md->dhc,
1272                                     static_cast<double>(state->fep_state),
1273                                     store_energy,
1274                                     pv,
1275                                     store_dhdl,
1276                                     md->dE + fep->lambda_start_n,
1277                                     time);
1278         }
1279     }
1280 }
1281
1282 }   // namespace
1283
1284 void EnergyOutput::recordNonEnergyStep()
1285 {
1286     ebin_increase_count(mdebin->ebin, false);
1287 }
1288
1289 namespace
1290 {
1291
1292 //! Legacy output function
1293 void npr(FILE *log, int n, char c)
1294 {
1295     for (; (n > 0); n--)
1296     {
1297         fprintf(log, "%c", c);
1298     }
1299 }
1300
1301 //! Legacy output function
1302 void pprint(FILE *log, const char *s, t_mdebin *md)
1303 {
1304     char CHAR = '#';
1305     int  slen;
1306     char buf1[22], buf2[22];
1307
1308     slen = strlen(s);
1309     fprintf(log, "\t<======  ");
1310     npr(log, slen, CHAR);
1311     fprintf(log, "  ==>\n");
1312     fprintf(log, "\t<====  %s  ====>\n", s);
1313     fprintf(log, "\t<==  ");
1314     npr(log, slen, CHAR);
1315     fprintf(log, "  ======>\n\n");
1316
1317     fprintf(log, "\tStatistics over %s steps using %s frames\n",
1318             gmx_step_str(md->ebin->nsteps_sim, buf1),
1319             gmx_step_str(md->ebin->nsum_sim, buf2));
1320     fprintf(log, "\n");
1321 }
1322
1323 }   // namespace
1324
1325 void print_ebin_header(FILE *log, int64_t steps, double time)
1326 {
1327     char buf[22];
1328
1329     fprintf(log, "   %12s   %12s\n"
1330             "   %12s   %12.5f\n\n",
1331             "Step", "Time", gmx_step_str(steps, buf), time);
1332 }
1333
1334 namespace
1335 {
1336
1337 // TODO It is too many responsibilities for this function to handle
1338 // both .edr and .log output for both per-time and time-average data.
1339 //! Legacy ebin output function
1340 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1341                 FILE *log,
1342                 int64_t step, double time,
1343                 int mode,
1344                 t_mdebin *md, t_fcdata *fcd,
1345                 SimulationGroups *groups, t_grpopts *opts,
1346                 gmx::Awh *awh)
1347 {
1348     /*static char **grpnms=NULL;*/
1349     char         buf[246];
1350     int          i, j, n, ni, nj, b;
1351     int          ndisre = 0;
1352
1353     /* these are for the old-style blocks (1 subblock, only reals), because
1354        there can be only one per ID for these */
1355     int          nr[enxNR];
1356     int          id[enxNR];
1357     real        *block[enxNR];
1358
1359     t_enxframe   fr;
1360
1361     if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1362     {
1363         if (log)
1364         {
1365             fprintf(log, "Not enough data recorded to report energy averages\n");
1366         }
1367         return;
1368     }
1369
1370     switch (mode)
1371     {
1372         case eprNORMAL:
1373             init_enxframe(&fr);
1374             fr.t            = time;
1375             fr.step         = step;
1376             fr.nsteps       = md->ebin->nsteps;
1377             fr.dt           = md->delta_t;
1378             fr.nsum         = md->ebin->nsum;
1379             fr.nre          = (bEne) ? md->ebin->nener : 0;
1380             fr.ener         = md->ebin->e;
1381             ndisre          = bDR ? fcd->disres.npair : 0;
1382             /* Optional additional old-style (real-only) blocks. */
1383             for (i = 0; i < enxNR; i++)
1384             {
1385                 nr[i] = 0;
1386             }
1387             if (bOR && fcd->orires.nr > 0)
1388             {
1389                 diagonalize_orires_tensors(&(fcd->orires));
1390                 nr[enxOR]     = fcd->orires.nr;
1391                 block[enxOR]  = fcd->orires.otav;
1392                 id[enxOR]     = enxOR;
1393                 nr[enxORI]    = (fcd->orires.oinsl != fcd->orires.otav) ?
1394                     fcd->orires.nr : 0;
1395                 block[enxORI] = fcd->orires.oinsl;
1396                 id[enxORI]    = enxORI;
1397                 nr[enxORT]    = fcd->orires.nex*12;
1398                 block[enxORT] = fcd->orires.eig;
1399                 id[enxORT]    = enxORT;
1400             }
1401
1402             /* whether we are going to wrte anything out: */
1403             if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1404             {
1405
1406                 /* the old-style blocks go first */
1407                 fr.nblock = 0;
1408                 for (i = 0; i < enxNR; i++)
1409                 {
1410                     if (nr[i] > 0)
1411                     {
1412                         fr.nblock = i + 1;
1413                     }
1414                 }
1415                 add_blocks_enxframe(&fr, fr.nblock);
1416                 for (b = 0; b < fr.nblock; b++)
1417                 {
1418                     add_subblocks_enxblock(&(fr.block[b]), 1);
1419                     fr.block[b].id        = id[b];
1420                     fr.block[b].sub[0].nr = nr[b];
1421 #if !GMX_DOUBLE
1422                     fr.block[b].sub[0].type = xdr_datatype_float;
1423                     fr.block[b].sub[0].fval = block[b];
1424 #else
1425                     fr.block[b].sub[0].type = xdr_datatype_double;
1426                     fr.block[b].sub[0].dval = block[b];
1427 #endif
1428                 }
1429
1430                 /* check for disre block & fill it. */
1431                 if (ndisre > 0)
1432                 {
1433                     int db = fr.nblock;
1434                     fr.nblock += 1;
1435                     add_blocks_enxframe(&fr, fr.nblock);
1436
1437                     add_subblocks_enxblock(&(fr.block[db]), 2);
1438                     fr.block[db].id        = enxDISRE;
1439                     fr.block[db].sub[0].nr = ndisre;
1440                     fr.block[db].sub[1].nr = ndisre;
1441 #if !GMX_DOUBLE
1442                     fr.block[db].sub[0].type = xdr_datatype_float;
1443                     fr.block[db].sub[1].type = xdr_datatype_float;
1444                     fr.block[db].sub[0].fval = fcd->disres.rt;
1445                     fr.block[db].sub[1].fval = fcd->disres.rm3tav;
1446 #else
1447                     fr.block[db].sub[0].type = xdr_datatype_double;
1448                     fr.block[db].sub[1].type = xdr_datatype_double;
1449                     fr.block[db].sub[0].dval = fcd->disres.rt;
1450                     fr.block[db].sub[1].dval = fcd->disres.rm3tav;
1451 #endif
1452                 }
1453                 /* here we can put new-style blocks */
1454
1455                 /* Free energy perturbation blocks */
1456                 if (md->dhc)
1457                 {
1458                     mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1459                 }
1460
1461                 /* we can now free & reset the data in the blocks */
1462                 if (md->dhc)
1463                 {
1464                     mde_delta_h_coll_reset(md->dhc);
1465                 }
1466
1467                 /* AWH bias blocks. */
1468                 if (awh != nullptr)  // TODO: add boolean in t_mdebin. See in mdebin.h.
1469                 {
1470                     awh->writeToEnergyFrame(step, &fr);
1471                 }
1472
1473                 /* do the actual I/O */
1474                 do_enx(fp_ene, &fr);
1475                 if (fr.nre)
1476                 {
1477                     /* We have stored the sums, so reset the sum history */
1478                     reset_ebin_sums(md->ebin);
1479                 }
1480             }
1481             free_enxframe(&fr);
1482             break;
1483         case eprAVER:
1484             if (log)
1485             {
1486                 pprint(log, "A V E R A G E S", md);
1487             }
1488             break;
1489         case eprRMS:
1490             if (log)
1491             {
1492                 pprint(log, "R M S - F L U C T U A T I O N S", md);
1493             }
1494             break;
1495         default:
1496             gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1497     }
1498
1499     if (log)
1500     {
1501         if (opts)
1502         {
1503             for (i = 0; i < opts->ngtc; i++)
1504             {
1505                 if (opts->annealing[i] != eannNO)
1506                 {
1507                     fprintf(log, "Current ref_t for group %s: %8.1f\n",
1508                             *(groups->groupNames[groups->groups[SimulationAtomGroupType::TemperatureCoupling].nm_ind[i]]),
1509                             opts->ref_t[i]);
1510                 }
1511             }
1512         }
1513         if (mode == eprNORMAL && bOR && fcd->orires.nr > 0)
1514         {
1515             print_orires_log(log, &(fcd->orires));
1516         }
1517         fprintf(log, "   Energies (%s)\n", unit_energy);
1518         pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1519         fprintf(log, "\n");
1520
1521         if (mode == eprAVER)
1522         {
1523             if (md->bDynBox)
1524             {
1525                 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1526                         mode, TRUE);
1527                 fprintf(log, "\n");
1528             }
1529             if (md->bConstrVir)
1530             {
1531                 fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1532                 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1533                 fprintf(log, "\n");
1534                 fprintf(log, "   Force Virial (%s)\n", unit_energy);
1535                 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1536                 fprintf(log, "\n");
1537             }
1538             if (md->bPres)
1539             {
1540                 fprintf(log, "   Total Virial (%s)\n", unit_energy);
1541                 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1542                 fprintf(log, "\n");
1543                 fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1544                 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1545                 fprintf(log, "\n");
1546             }
1547             if (md->bMu)
1548             {
1549                 fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1550                 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1551                 fprintf(log, "\n");
1552             }
1553
1554             if (md->nE > 1)
1555             {
1556                 if (md->print_grpnms == nullptr)
1557                 {
1558                     snew(md->print_grpnms, md->nE);
1559                     n = 0;
1560                     for (i = 0; (i < md->nEg); i++)
1561                     {
1562                         ni = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i];
1563                         for (j = i; (j < md->nEg); j++)
1564                         {
1565                             nj = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[j];
1566                             sprintf(buf, "%s-%s", *(groups->groupNames[ni]),
1567                                     *(groups->groupNames[nj]));
1568                             md->print_grpnms[n++] = gmx_strdup(buf);
1569                         }
1570                     }
1571                 }
1572                 sprintf(buf, "Epot (%s)", unit_energy);
1573                 fprintf(log, "%15s   ", buf);
1574                 for (i = 0; (i < egNR); i++)
1575                 {
1576                     if (md->bEInd[i])
1577                     {
1578                         fprintf(log, "%12s   ", egrp_nm[i]);
1579                     }
1580                 }
1581                 fprintf(log, "\n");
1582                 for (i = 0; (i < md->nE); i++)
1583                 {
1584                     fprintf(log, "%15s", md->print_grpnms[i]);
1585                     pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1586                             FALSE);
1587                 }
1588                 fprintf(log, "\n");
1589             }
1590             if (md->nTC > 1)
1591             {
1592                 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1593                 fprintf(log, "\n");
1594             }
1595             if (md->nU > 1)
1596             {
1597                 fprintf(log, "%15s   %12s   %12s   %12s\n",
1598                         "Group", "Ux", "Uy", "Uz");
1599                 for (i = 0; (i < md->nU); i++)
1600                 {
1601                     ni = groups->groups[SimulationAtomGroupType::Acceleration].nm_ind[i];
1602                     fprintf(log, "%15s", *groups->groupNames[ni]);
1603                     pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1604                 }
1605                 fprintf(log, "\n");
1606             }
1607         }
1608     }
1609
1610 }
1611
1612 //! Legacy update function
1613 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1614 {
1615     const t_ebin * const ebin = mdebin->ebin;
1616
1617     enerhist->nsteps     = ebin->nsteps;
1618     enerhist->nsum       = ebin->nsum;
1619     enerhist->nsteps_sim = ebin->nsteps_sim;
1620     enerhist->nsum_sim   = ebin->nsum_sim;
1621
1622     if (ebin->nsum > 0)
1623     {
1624         /* This will only actually resize the first time */
1625         enerhist->ener_ave.resize(ebin->nener);
1626         enerhist->ener_sum.resize(ebin->nener);
1627
1628         for (int i = 0; i < ebin->nener; i++)
1629         {
1630             enerhist->ener_ave[i] = ebin->e[i].eav;
1631             enerhist->ener_sum[i] = ebin->e[i].esum;
1632         }
1633     }
1634
1635     if (ebin->nsum_sim > 0)
1636     {
1637         /* This will only actually resize the first time */
1638         enerhist->ener_sum_sim.resize(ebin->nener);
1639
1640         for (int i = 0; i < ebin->nener; i++)
1641         {
1642             enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1643         }
1644     }
1645     if (mdebin->dhc)
1646     {
1647         mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1648     }
1649 }
1650
1651 //! Legacy restore function
1652 void restore_energyhistory_from_state(t_mdebin              * mdebin,
1653                                       const energyhistory_t * enerhist)
1654 {
1655     unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1656
1657     GMX_RELEASE_ASSERT(enerhist, "Need valid history to restore");
1658
1659     if ((enerhist->nsum     > 0 && nener != enerhist->ener_sum.size()) ||
1660         (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1661     {
1662         gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
1663                   nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1664     }
1665
1666     mdebin->ebin->nsteps     = enerhist->nsteps;
1667     mdebin->ebin->nsum       = enerhist->nsum;
1668     mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1669     mdebin->ebin->nsum_sim   = enerhist->nsum_sim;
1670
1671     for (int i = 0; i < mdebin->ebin->nener; i++)
1672     {
1673         mdebin->ebin->e[i].eav  =
1674             (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1675         mdebin->ebin->e[i].esum =
1676             (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1677         mdebin->ebin->e_sim[i].esum =
1678             (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1679     }
1680     if (mdebin->dhc)
1681     {
1682         mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());
1683     }
1684 }
1685
1686 }   // namespace
1687
1688 EnergyOutput::EnergyOutput()
1689     : mdebin(nullptr)
1690 {
1691 }
1692
1693 void EnergyOutput::prepare(ener_file        *fp_ene,
1694                            const gmx_mtop_t *mtop,
1695                            const t_inputrec *ir,
1696                            const pull_t     *pull_work,
1697                            FILE             *fp_dhdl,
1698                            bool              isRerun)
1699 {
1700     mdebin = init_mdebin(fp_ene, mtop, ir, pull_work, fp_dhdl, isRerun);
1701 }
1702
1703 EnergyOutput::~EnergyOutput()
1704 {
1705     done_mdebin(mdebin);
1706 }
1707
1708 t_ebin *EnergyOutput::getEbin()
1709 {
1710     return mdebin->ebin;
1711 }
1712
1713 void EnergyOutput::addDataAtEnergyStep(bool                    bDoDHDL,
1714                                        bool                    bSum,
1715                                        double                  time,
1716                                        real                    tmass,
1717                                        gmx_enerdata_t         *enerd,
1718                                        t_state                *state,
1719                                        t_lambda               *fep,
1720                                        t_expanded             *expand,
1721                                        matrix                  box,
1722                                        tensor                  svir,
1723                                        tensor                  fvir,
1724                                        tensor                  vir,
1725                                        tensor                  pres,
1726                                        gmx_ekindata_t         *ekind,
1727                                        rvec                    mu_tot,
1728                                        const gmx::Constraints *constr)
1729 {
1730     upd_mdebin(mdebin, bDoDHDL, bSum, time, tmass, enerd, state, fep,
1731                expand, box, svir, fvir, vir, pres, ekind, mu_tot, constr);
1732 }
1733
1734 void EnergyOutput::printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
1735                                          FILE *log,
1736                                          int64_t step, double time,
1737                                          int mode,
1738                                          t_fcdata *fcd,
1739                                          SimulationGroups *groups, t_grpopts *opts,
1740                                          gmx::Awh *awh)
1741 {
1742     print_ebin(fp_ene, bEne, bDR, bOR, log, step, time, mode,
1743                mdebin, fcd, groups, opts, awh);
1744 }
1745
1746 int EnergyOutput::numEnergyTerms() const
1747 {
1748     return mdebin->ebin->nener;
1749 }
1750
1751 void EnergyOutput::fillEnergyHistory(energyhistory_t *enerhist) const
1752 {
1753     update_energyhistory(enerhist, mdebin);
1754 }
1755
1756 void EnergyOutput::restoreFromEnergyHistory(const energyhistory_t &enerhist)
1757 {
1758     restore_energyhistory_from_state(mdebin, &enerhist);
1759 }
1760
1761 } // namespace gmx