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