Apply clang-format to source tree
[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, 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 "gromacs/mdtypes/enerdata.h"
51
52 class energyhistory_t;
53 struct ener_file;
54 struct gmx_ekindata_t;
55 struct gmx_enerdata_t;
56 struct SimulationGroups;
57 struct gmx_mtop_t;
58 struct gmx_output_env_t;
59 struct pull_t;
60 struct t_ebin;
61 struct t_expanded;
62 struct t_fcdata;
63 struct t_grpopts;
64 struct t_inputrec;
65 struct t_lambda;
66 class t_state;
67
68 struct t_mde_delta_h_coll;
69
70 namespace gmx
71 {
72 class Awh;
73 class Constraints;
74 struct MdModulesNotifier;
75 } // namespace gmx
76
77 //! \brief Printed names for intergroup energies
78 extern const char* egrp_nm[egNR + 1];
79
80 /* \brief delta_h block type enum: the kinds of energies written out. */
81 enum
82 {
83     //! Delta H BAR energy difference
84     dhbtDH = 0,
85     //! dH/dlambda derivative
86     dhbtDHDL = 1,
87     //! System energy
88     dhbtEN,
89     //! pV term
90     dhbtPV,
91     //! Expanded ensemble statistics
92     dhbtEXPANDED,
93     //! Total number of energy types in this enum
94     dhbtNR
95 };
96
97 namespace gmx
98 {
99
100 /* Energy output class
101  *
102  * This is the collection of energy averages collected during mdrun, and to
103  * be written out to the .edr file.
104  *
105  * \todo Use more std containers.
106  * \todo Remove GMX_CONSTRAINTVIR
107  * \todo Write free-energy output also to energy file (after adding more tests)
108  */
109 class EnergyOutput
110 {
111 public:
112     /*! \brief Initiate MD energy bin
113      *
114      * \param[in] fp_ene     Energy output file.
115      * \param[in] mtop       Topology.
116      * \param[in] ir         Input parameters.
117      * \param[in] pull_work  Pulling simulations data
118      * \param[in] fp_dhdl    FEP file.
119      * \param[in] isRerun    Is this is a rerun instead of the simulations.
120      * \param[in] mdModulesNotifier Notifications to MD modules.
121      */
122     EnergyOutput(ener_file*               fp_ene,
123                  const gmx_mtop_t*        mtop,
124                  const t_inputrec*        ir,
125                  const pull_t*            pull_work,
126                  FILE*                    fp_dhdl,
127                  bool                     isRerun,
128                  const MdModulesNotifier& mdModulesNotifier);
129
130     ~EnergyOutput();
131
132     /*! \brief Update the averaging structures.
133      *
134      * Called every step on which the thermodynamic values are evaluated.
135      *
136      * \param[in] bDoDHDL  Whether the FEP is enabled.
137      * \param[in] bSum     If this stepshould be recorded to compute sums and averaes.
138      * \param[in] time     Current simulation time.
139      * \param[in] tmass    Total mass
140      * \param[in] enerd    Energy data object.
141      * \param[in] state    System state.
142      * \param[in] fep      FEP data.
143      * \param[in] expand   Expanded ensemble (for FEP).
144      * \param[in] lastbox  PBC data.
145      * \param[in] svir     Constraint virial.
146      * \param[in] fvir     Force virial.
147      * \param[in] vir      Total virial.
148      * \param[in] pres     Pressure.
149      * \param[in] ekind    Kinetic energy data.
150      * \param[in] mu_tot   Total dipole.
151      * \param[in] constr   Constraints object to get RMSD from (for LINCS only).
152      */
153     void addDataAtEnergyStep(bool                    bDoDHDL,
154                              bool                    bSum,
155                              double                  time,
156                              real                    tmass,
157                              const gmx_enerdata_t*   enerd,
158                              const t_state*          state,
159                              const t_lambda*         fep,
160                              const t_expanded*       expand,
161                              const matrix            lastbox,
162                              const tensor            svir,
163                              const tensor            fvir,
164                              const tensor            vir,
165                              const tensor            pres,
166                              const gmx_ekindata_t*   ekind,
167                              const rvec              mu_tot,
168                              const gmx::Constraints* constr);
169
170     /*! \brief Update the data averaging structure counts.
171      *
172      * Updates the number of steps, the values have not being computed.
173      */
174     void recordNonEnergyStep();
175
176     /*! \brief Writes current quantites to log and energy files.
177      *
178      * Prints current values of energies, pressure, temperature, restraint
179      * data, etc. to energy output file and to the log file (if not nullptr).
180      *
181      * This function only does something useful when bEne || bDR || bOR || log.
182      *
183      * \todo Perhaps this responsibility should involve some other
184      *       object visiting all the contributing objects.
185      *
186      * \param[in] fp_ene   Energy file for the output.
187      * \param[in] bEne     If it is a step for energy output or last step.
188      * \param[in] bDR      If it is a step of writing distance restraints.
189      * \param[in] bOR      If it is a step of writing orientation restraints.
190      * \param[in] log      Pointer to the log file.
191      * \param[in] step     Current step.
192      * \param[in] time     Current simulation time.
193      * \param[in] fcd      Bonded force computation data,
194      *                     including orientation and distance restraints.
195      * \param[in] awh      AWH data.
196      */
197     void printStepToEnergyFile(ener_file* fp_ene,
198                                bool       bEne,
199                                bool       bDR,
200                                bool       bOR,
201                                FILE*      log,
202                                int64_t    step,
203                                double     time,
204                                t_fcdata*  fcd,
205                                gmx::Awh*  awh);
206
207     /*! \brief Print reference temperatures for annealing groups.
208      *
209      * Nothing is done if log is nullptr.
210      *
211      * \param[in] log     Log file to print to.
212      * \param[in] groups  Information on atom groups.
213      * \param[in] opts    Atom temperature coupling groups options
214      *                    (annealing is done by groups).
215      */
216     void printAnnealingTemperatures(FILE* log, SimulationGroups* groups, t_grpopts* opts);
217
218     /*! \brief Prints average values to log file.
219      *
220      * This is called at the end of the simulation run to print accumulated average values.
221      * Nothing it done if log is nullptr.
222      *
223      * \param[in]   log      Where to print.
224      * \param[in]   groups   Atom groups.
225      */
226     void printAverages(FILE* log, const SimulationGroups* groups);
227
228     /*! \brief Get the number of thermodynamic terms recorded.
229      *
230      * \todo Refactor this to return the expected output size,
231      *       rather than exposing the implementation details about
232      *       thermodynamic terms.
233      *
234      * \returns Number of thermodynamic terms.
235      */
236     int numEnergyTerms() const;
237
238     /*! \brief Fill the energyhistory_t data.
239      *
240      * Between .edr writes, the averages are history dependent,
241      * and that history needs to be retained in checkpoints.
242      * These functions set/read the energyhistory_t class
243      * that is written to checkpoints.
244      *
245      * \param[out] enerhist  Energy history data structure.
246      */
247     void fillEnergyHistory(energyhistory_t* enerhist) const;
248
249     /*! \brief Restore from energyhistory_t data.
250      *
251      * \param[in] enerhist  Energy history data structure.
252      */
253     void restoreFromEnergyHistory(const energyhistory_t& enerhist);
254
255     //! Print an output header to the log file.
256     void printHeader(FILE* log, int64_t steps, double time);
257
258 private:
259     //! Timestep
260     double delta_t_ = 0;
261
262     //! Structure to store energy components and their running averages
263     t_ebin* ebin_ = nullptr;
264
265     //! Is the periodic box triclinic
266     bool bTricl_ = false;
267     //! NHC trotter is used
268     bool bNHC_trotter_ = false;
269     //! If Nose-Hoover chains should be printed
270     bool bPrintNHChains_ = false;
271     //! If MTTK integrator was used
272     bool bMTTK_ = false;
273
274     //! Temperature control scheme
275     int etc_ = 0;
276
277     //! Which of the main energy terms should be printed
278     bool bEner_[F_NRE] = { false };
279     //! Index for main energy terms
280     int ie_ = 0;
281     //! Number of energy terms from F_NRE list to be saved (i.e. number of 'true' in bEner)
282     int f_nre_ = 0;
283
284     //! Index for constraints RMSD
285     int iconrmsd_ = 0;
286     /* !\brief How many constraints RMSD terms there are.
287      * Usually 1 if LINCS is used and 0 otherwise)
288      * nCrmsd > 0 indicates when constraints RMSD is saves, hence no boolean
289      */
290     int nCrmsd_ = 0;
291
292     //! Is the periodic box dynamic
293     bool bDynBox_ = false;
294     //! Index for box dimensions
295     int ib_ = 0;
296     //! Index for box volume
297     int ivol_ = 0;
298     //! Index for density
299     int idens_ = 0;
300     //! Triclinic box and not a rerun
301     bool bDiagPres_ = false;
302     //! Reference pressure, averaged over dimensions
303     real ref_p_ = 0.0;
304     //! Index for thermodynamic work (pV)
305     int ipv_ = 0;
306     //! Index for entalpy (pV + total energy)
307     int ienthalpy_ = 0;
308
309     /*! \brief If the constraints virial should be printed.
310      * Can only be true if "GMX_CONSTRAINTVIR" environmental variable is set */
311     bool bConstrVir_ = false;
312     //! Index for constrains virial
313     int isvir_ = 0;
314     //! Index for force virial
315     int ifvir_ = 0;
316
317     //! If we have pressure computed
318     bool bPres_ = false;
319     //! Index for total virial
320     int ivir_ = 0;
321     //! Index for pressure
322     int ipres_ = 0;
323     /*! \brief Index for surface tension
324      * [(pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ]]*/
325     int isurft_ = 0;
326
327     //! Pressure control scheme
328     int epc_ = 0;
329     //! Index for velocity of the box borders
330     int ipc_ = 0;
331
332     //! If dipole was calculated
333     bool bMu_ = false;
334     //! Index for the dipole
335     int imu_ = 0;
336
337     //! Index for coseine acceleration used for viscocity calculation
338     int ivcos_ = 0;
339     //! Index for viscocity
340     int ivisc_ = 0;
341
342     //! Which energy terms from egNR list should be printed in group-to-group block
343     bool bEInd_[egNR] = { false };
344     //! Number of energy terms to be printed (these, for which bEInd[] == true)
345     int nEc_ = 0;
346     //! Number of energy output groups
347     int nEg_ = 0;
348     //! Number of intergroup energy sets to be printed for each energy term (nE = (nEg*(nEg+1))/2)
349     int nE_ = 0;
350     //! Indexes for integroup energy sets (each set with nEc energies)
351     int* igrp_ = nullptr;
352
353     //! Number of temperature coupling groups
354     int nTC_ = 0;
355     //! Index for temperature
356     int itemp_ = 0;
357
358     //! Number of Nose-Hoover chains
359     int nNHC_ = 0;
360     //! Number of temperature coupling coefficients in case of NH Chains
361     int mde_n_ = 0;
362     //! Index for temperature coupling coefficient in case of NH chains
363     int itc_ = 0;
364
365     //! Number of temperature coupling terms if the temperature control is dealt by barostat (MTTK case)
366     int nTCP_ = 0;
367     //! Scalling factors for temperaturs control in MTTK
368     int mdeb_n_ = 0;
369     //! Index for scalling factor of MTTK
370     int itcb_ = 0;
371
372     //! Number of acceleration groups
373     int nU_ = 0;
374     //! Index for group velocities
375     int iu_ = 0;
376
377     //! Array to accumulate values during update
378     real* tmp_r_ = nullptr;
379     //! Array to accumulate values during update
380     rvec* tmp_v_ = nullptr;
381
382     //! The dhdl.xvg output file
383     FILE* fp_dhdl_ = nullptr;
384     //! Energy components for dhdl.xvg output
385     double* dE_ = nullptr;
386     //! The delta U components (raw data + histogram)
387     t_mde_delta_h_coll* dhc_ = nullptr;
388     //! Temperatures for simulated tempering groups
389     real* temperatures_ = nullptr;
390     //! Number of temperatures actually saved
391     int numTemperatures_ = 0;
392 };
393
394 } // namespace gmx
395
396 //! Open the dhdl file for output
397 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv);
398
399 #endif