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