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