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