Add lambda state to dhdl.xvg with AWH fep
[alexxy/gromacs.git] / src / gromacs / mdlib / energyoutput.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \libinternal \file
36  * \brief Header for the code that writes energy-like quantities.
37  *
38  * \author Mark Abraham <mark.j.abraham@gmail.com>
39  * \author Paul Bauer <paul.bauer.q@gmail.com>
40  * \author Artem Zhmurov <zhmurov@gmail.com>
41  *
42  * \ingroup module_mdlib
43  * \inlibraryapi
44  */
45 #ifndef GMX_MDLIB_ENERGYOUTPUT_H
46 #define GMX_MDLIB_ENERGYOUTPUT_H
47
48 #include <cstdio>
49
50 #include <memory>
51
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/mdtypes/enerdata.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/real.h"
58
59 class energyhistory_t;
60 struct ener_file;
61 class gmx_ekindata_t;
62 struct gmx_enerdata_t;
63 struct SimulationGroups;
64 struct gmx_mtop_t;
65 struct gmx_output_env_t;
66 struct pull_t;
67 struct t_ebin;
68 struct t_fcdata;
69 struct t_grpopts;
70 struct t_inputrec;
71 struct t_lambda;
72 class t_state;
73
74 struct t_mde_delta_h_coll;
75
76 namespace gmx
77 {
78 class Awh;
79 class Constraints;
80 struct MDModulesNotifiers;
81 enum class StartingBehavior;
82 } // namespace gmx
83
84 //! \brief Printed names for intergroup energies
85 const char* enumValueToString(NonBondedEnergyTerms enumValue);
86
87 /* \brief delta_h block type enum: the kinds of energies written out. */
88 enum
89 {
90     //! Delta H BAR energy difference
91     dhbtDH = 0,
92     //! dH/dlambda derivative
93     dhbtDHDL = 1,
94     //! System energy
95     dhbtEN,
96     //! pV term
97     dhbtPV,
98     //! Expanded ensemble statistics
99     dhbtEXPANDED,
100     //! Total number of energy types in this enum
101     dhbtNR
102 };
103
104 namespace gmx
105 {
106 class EnergyDriftTracker;
107
108 /*! \internal
109  * \brief Arrays connected to Pressure and Temperature coupling
110  */
111 struct PTCouplingArrays
112 {
113     //! Box velocities for Parrinello-Rahman P-coupling.
114     const rvec* boxv;
115     //! Nose-Hoover coordinates.
116     ArrayRef<const double> nosehoover_xi;
117     //! Nose-Hoover velocities.
118     ArrayRef<const double> nosehoover_vxi;
119     //! Pressure Nose-Hoover coordinates.
120     ArrayRef<const double> nhpres_xi;
121     //! Pressure Nose-Hoover velocities.
122     ArrayRef<const double> nhpres_vxi;
123 };
124
125 /* Energy output class
126  *
127  * This is the collection of energy averages collected during mdrun, and to
128  * be written out to the .edr file.
129  *
130  * \todo Use more std containers.
131  * \todo Write free-energy output also to energy file (after adding more tests)
132  */
133 class EnergyOutput
134 {
135 public:
136     /*! \brief Initiate MD energy bin
137      *
138      * \param[in] fp_ene     Energy output file.
139      * \param[in] mtop       Topology.
140      * \param[in] inputrec   Input parameters.
141      * \param[in] pull_work  Pulling simulations data
142      * \param[in] fp_dhdl    FEP file.
143      * \param[in] isRerun    Is this is a rerun instead of the simulations.
144      * \param[in] startingBehavior  Run starting behavior.
145      * \param[in] simulationsShareState  Tells whether the physical state is shared over simulations
146      * \param[in] mdModulesNotifiers Notifications to MD modules.
147      */
148     EnergyOutput(ener_file*                fp_ene,
149                  const gmx_mtop_t&         mtop,
150                  const t_inputrec&         inputrec,
151                  const pull_t*             pull_work,
152                  FILE*                     fp_dhdl,
153                  bool                      isRerun,
154                  StartingBehavior          startingBehavior,
155                  bool                      simulationsShareState,
156                  const MDModulesNotifiers& mdModulesNotifiers);
157
158     ~EnergyOutput();
159
160     /*! \brief Update the averaging structures.
161      *
162      * Called every step on which the thermodynamic values are evaluated.
163      *
164      * \param[in] bDoDHDL           Whether the FEP is enabled.
165      * \param[in] bSum              If this stepshould be recorded to compute sums and averages.
166      * \param[in] time              Current simulation time.
167      * \param[in] tmass             Total mass
168      * \param[in] enerd             Energy data object.
169      * \param[in] fep               FEP data.
170      * \param[in] lastbox           PBC data.
171      * \param[in] ptCouplingArrays  Arrays connected to pressure and temperature coupling.
172      * \param[in] fep_state         The current alchemical state we are in.
173      * \param[in] vir               Total virial.
174      * \param[in] pres              Pressure.
175      * \param[in] ekind             Kinetic energy data.
176      * \param[in] mu_tot            Total dipole.
177      * \param[in] constr            Constraints object to get RMSD from (for LINCS only).
178      */
179     void addDataAtEnergyStep(bool                    bDoDHDL,
180                              bool                    bSum,
181                              double                  time,
182                              real                    tmass,
183                              const gmx_enerdata_t*   enerd,
184                              const t_lambda*         fep,
185                              const matrix            lastbox,
186                              PTCouplingArrays        ptCouplingArrays,
187                              int                     fep_state,
188                              const tensor            vir,
189                              const tensor            pres,
190                              const gmx_ekindata_t*   ekind,
191                              const rvec              mu_tot,
192                              const gmx::Constraints* constr);
193
194     /*! \brief Update the data averaging structure counts.
195      *
196      * Updates the number of steps, the values have not being computed.
197      */
198     void recordNonEnergyStep();
199
200     /*! \brief Writes current quantites to log and energy files.
201      *
202      * Prints current values of energies, pressure, temperature, restraint
203      * data, etc. to energy output file and to the log file (if not nullptr).
204      *
205      * This function only does something useful when bEne || bDR || bOR || log.
206      *
207      * \todo Perhaps this responsibility should involve some other
208      *       object visiting all the contributing objects.
209      *
210      * \param[in] fp_ene   Energy file for the output.
211      * \param[in] bEne     If it is a step for energy output or last step.
212      * \param[in] bDR      If it is a step of writing distance restraints.
213      * \param[in] bOR      If it is a step of writing orientation restraints.
214      * \param[in] log      Pointer to the log file.
215      * \param[in] step     Current step.
216      * \param[in] time     Current simulation time.
217      * \param[in] fcd      Bonded force computation data,
218      *                     including orientation and distance restraints.
219      * \param[in] awh      AWH data.
220      */
221     void printStepToEnergyFile(ener_file* fp_ene,
222                                bool       bEne,
223                                bool       bDR,
224                                bool       bOR,
225                                FILE*      log,
226                                int64_t    step,
227                                double     time,
228                                t_fcdata*  fcd,
229                                gmx::Awh*  awh);
230
231     /*! \brief Print reference temperatures for annealing groups.
232      *
233      * Nothing is done if log is nullptr.
234      *
235      * \param[in] log     Log file to print to.
236      * \param[in] groups  Information on atom groups.
237      * \param[in] opts    Atom temperature coupling groups options
238      *                    (annealing is done by groups).
239      */
240     static void printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts);
241
242     /*! \brief Prints average values to log file.
243      *
244      * This is called at the end of the simulation run to print accumulated average values.
245      * Nothing it done if log is nullptr.
246      *
247      * \param[in]   log      Where to print.
248      * \param[in]   groups   Atom groups.
249      */
250     void printAverages(FILE* log, const SimulationGroups* groups);
251
252     /*! \brief Get the number of thermodynamic terms recorded.
253      *
254      * \todo Refactor this to return the expected output size,
255      *       rather than exposing the implementation details about
256      *       thermodynamic terms.
257      *
258      * \returns Number of thermodynamic terms.
259      */
260     int numEnergyTerms() const;
261
262     /*! \brief Fill the energyhistory_t data.
263      *
264      * Between .edr writes, the averages are history dependent,
265      * and that history needs to be retained in checkpoints.
266      * These functions set/read the energyhistory_t class
267      * that is written to checkpoints.
268      *
269      * \param[out] enerhist  Energy history data structure.
270      */
271     void fillEnergyHistory(energyhistory_t* enerhist) const;
272
273     /*! \brief Restore from energyhistory_t data.
274      *
275      * \param[in] enerhist  Energy history data structure.
276      */
277     void restoreFromEnergyHistory(const energyhistory_t& enerhist);
278
279     //! Print an output header to the log file.
280     static void printHeader(FILE* log, int64_t steps, double time);
281
282     /*! \brief Print conserved energy drift message to \p fplog
283      *
284      * Note that this is only over the current run (times are printed),
285      * this is not from the original start time for runs with continuation.
286      * This has the advantage that one can find if conservation issues are
287      * from the current run with the current settings on the current hardware.
288      */
289     void printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const;
290
291 private:
292     //! Timestep
293     double delta_t_ = 0;
294
295     //! Structure to store energy components and their running averages
296     t_ebin* ebin_ = nullptr;
297
298     //! Is the periodic box triclinic
299     bool bTricl_ = false;
300     //! NHC trotter is used
301     bool bNHC_trotter_ = false;
302     //! If Nose-Hoover chains should be printed
303     bool bPrintNHChains_ = false;
304     //! If MTTK integrator was used
305     bool bMTTK_ = false;
306
307     //! Temperature control scheme
308     TemperatureCoupling etc_ = TemperatureCoupling::No;
309
310     //! Which of the main energy terms should be printed
311     bool bEner_[F_NRE] = { false };
312     //! Index for main energy terms
313     int ie_ = 0;
314     //! Number of energy terms from F_NRE list to be saved (i.e. number of 'true' in bEner)
315     int f_nre_ = 0;
316
317     //! Index for constraints RMSD
318     int iconrmsd_ = 0;
319     /* !\brief How many constraints RMSD terms there are.
320      * Usually 1 if LINCS is used and 0 otherwise)
321      * nCrmsd > 0 indicates when constraints RMSD is saves, hence no boolean
322      */
323     int nCrmsd_ = 0;
324
325     //! Is the periodic box dynamic
326     bool bDynBox_ = false;
327     //! Index for box dimensions
328     int ib_ = 0;
329     //! Index for box volume
330     int ivol_ = 0;
331     //! Index for density
332     int idens_ = 0;
333     //! Triclinic box and not a rerun
334     bool bDiagPres_ = false;
335     //! Reference pressure, averaged over dimensions
336     real ref_p_ = 0.0;
337     //! Index for thermodynamic work (pV)
338     int ipv_ = 0;
339     //! Index for entalpy (pV + total energy)
340     int ienthalpy_ = 0;
341
342     //! If we have pressure computed
343     bool bPres_ = false;
344     //! Index for total virial
345     int ivir_ = 0;
346     //! Index for pressure
347     int ipres_ = 0;
348     /*! \brief Index for surface tension
349      * [(pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ]]*/
350     int isurft_ = 0;
351
352     //! Pressure control scheme
353     PressureCoupling epc_ = PressureCoupling::No;
354     //! Index for velocity of the box borders
355     int ipc_ = 0;
356
357     //! If dipole was calculated
358     bool bMu_ = false;
359     //! Index for the dipole
360     int imu_ = 0;
361
362     //! Index for coseine acceleration used for viscocity calculation
363     int ivcos_ = 0;
364     //! Index for viscocity
365     int ivisc_ = 0;
366
367     //! Which energy terms from NonBondedEnergyTerms list should be printed in group-to-group block
368     gmx::EnumerationArray<NonBondedEnergyTerms, bool> bEInd_;
369     //! Number of energy terms to be printed (these, for which bEInd[] == true)
370     int nEc_ = 0;
371     //! Number of energy output groups
372     int nEg_ = 0;
373     //! Number of intergroup energy sets to be printed for each energy term (nE = (nEg*(nEg+1))/2)
374     int nE_ = 0;
375     //! Indexes for integroup energy sets (each set with nEc energies)
376     std::vector<int> igrp_;
377
378     //! Number of temperature coupling groups
379     int nTC_ = 0;
380     //! Index for temperature
381     int itemp_ = 0;
382
383     //! Number of Nose-Hoover chains
384     int nNHC_ = 0;
385     //! Number of temperature coupling coefficients in case of NH Chains
386     int mde_n_ = 0;
387     //! Index for temperature coupling coefficient in case of NH chains
388     int itc_ = 0;
389
390     //! Number of temperature coupling terms if the temperature control is dealt by barostat (MTTK case)
391     int nTCP_ = 0;
392     //! Scalling factors for temperaturs control in MTTK
393     int mdeb_n_ = 0;
394     //! Index for scalling factor of MTTK
395     int itcb_ = 0;
396
397     //! Array to accumulate values during update
398     std::vector<real> tmp_r_;
399
400     //! The dhdl.xvg output file
401     FILE* fp_dhdl_ = nullptr;
402     //! Whether the free-energy lambda moves dynamically between lambda states
403     bool haveFepLambdaMoves_;
404     //! Energy components for dhdl.xvg output
405     std::vector<double> dE_;
406     //! The delta U components (raw data + histogram)
407     std::unique_ptr<t_mde_delta_h_coll> dhc_;
408     //! Temperatures for simulated tempering groups
409     std::vector<real> temperatures_;
410
411     //! For tracking the conserved or total energy
412     std::unique_ptr<EnergyDriftTracker> conservedEnergyTracker_;
413 };
414
415 } // namespace gmx
416
417 //! Open the dhdl file for output
418 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv);
419
420 #endif