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