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