Move pull work struct out of inputrec
[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 #ifndef GMX_MDLIB_ENERGYOUTPUT_H
36 #define GMX_MDLIB_ENERGYOUTPUT_H
37
38 #include <cstdio>
39
40 #include "gromacs/mdtypes/enerdata.h"
41
42 class energyhistory_t;
43 struct ener_file;
44 struct gmx_ekindata_t;
45 struct gmx_enerdata_t;
46 struct SimulationGroups;
47 struct gmx_mtop_t;
48 struct gmx_output_env_t;
49 struct pull_t;
50 struct t_ebin;
51 struct t_expanded;
52 struct t_fcdata;
53 struct t_grpopts;
54 struct t_inputrec;
55 struct t_lambda;
56 class t_state;
57
58 namespace gmx
59 {
60 class Awh;
61 class Constraints;
62 }
63
64 extern const char *egrp_nm[egNR+1];
65
66 /* delta_h block type enum: the kinds of energies written out. */
67 enum
68 {
69     dhbtDH   = 0, /* delta H BAR energy difference*/
70     dhbtDHDL = 1, /* dH/dlambda derivative */
71     dhbtEN,       /* System energy */
72     dhbtPV,       /* pV term */
73     dhbtEXPANDED, /* expanded ensemble statistics */
74     dhbtNR
75 };
76
77 namespace gmx
78 {
79
80 // TODO remove use of detail namespace when removing t_mdebin in
81 // favour of an Impl class.
82 namespace detail
83 {
84 struct t_mdebin;
85 }
86
87 /* The functions & data structures here determine the content for outputting
88    the .edr file; the file format and actual writing is done with functions
89    defined in enxio.h */
90
91 class EnergyOutput
92 {
93     public:
94         EnergyOutput();
95         /*! \brief Initiate MD energy bin
96          *
97          * This second phase of construction is needed until we have
98          * modules that understand how to request output from
99          * EnergyOutput.
100          *
101          * \todo Refactor to separate a function to write the energy
102          * file header. Perhaps transform the remainder into a factory
103          * function.
104          */
105         void prepare(ener_file        *fp_ene,
106                      const gmx_mtop_t *mtop,
107                      const t_inputrec *ir,
108                      const pull_t     *pull_work,
109                      FILE             *fp_dhdl,
110                      bool              isRerun = false);
111         ~EnergyOutput();
112         /*! \brief Update the averaging structures.
113          *
114          * Called every step on which the energies are evaluated. */
115         void addDataAtEnergyStep(bool                    bDoDHDL,
116                                  bool                    bSum,
117                                  double                  time,
118                                  real                    tmass,
119                                  gmx_enerdata_t         *enerd,
120                                  t_state                *state,
121                                  t_lambda               *fep,
122                                  t_expanded             *expand,
123                                  matrix                  lastbox,
124                                  tensor                  svir,
125                                  tensor                  fvir,
126                                  tensor                  vir,
127                                  tensor                  pres,
128                                  gmx_ekindata_t         *ekind,
129                                  rvec                    mu_tot,
130                                  const gmx::Constraints *constr);
131         /*! \brief Updated the averaging structures
132          *
133          * Called every step on which the energies are not evaluated.
134          *
135          * \todo This schedule is known in advance, and should be made
136          * an intrinsic behaviour of EnergyOutput, rather than being
137          * wastefully called every step. */
138         void recordNonEnergyStep();
139
140         /*! \brief Help write quantites to the energy file
141          *
142          * \todo Perhaps this responsibility should involve some other
143          * object visiting all the contributing objects. */
144         void printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
145                                    FILE *log,
146                                    int64_t step, double time,
147                                    int mode,
148                                    t_fcdata *fcd,
149                                    SimulationGroups *groups, t_grpopts *opts,
150                                    gmx::Awh *awh);
151         /*! \brief Get the number of energy terms recorded.
152          *
153          * \todo Refactor this to return the expected output size,
154          * rather than exposing the implementation details about
155          * energy terms. */
156         int numEnergyTerms() const;
157         /*! \brief Getter used for testing t_ebin
158          *
159          * \todo Find a better approach for this. */
160         t_ebin *getEbin();
161
162         /* Between .edr writes, the averages are history dependent,
163            and that history needs to be retained in checkpoints.
164            These functions set/read the energyhistory_t class
165            that is written to checkpoints in checkpoint.c */
166
167         //! Fill the energyhistory_t data.
168         void fillEnergyHistory(energyhistory_t * enerhist) const;
169         //! Restore from energyhistory_t data.
170         void restoreFromEnergyHistory(const energyhistory_t &enerhist);
171
172     private:
173         // TODO transform this into an impl class.
174         detail::t_mdebin *mdebin = nullptr;
175 };
176
177 } // namespace gmx
178
179 //! Open the dhdl file for output
180 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
181                 const gmx_output_env_t *oenv);
182
183 namespace gmx
184 {
185
186 //! Print an energy-output header to the log file
187 void print_ebin_header(FILE *log, int64_t steps, double time);
188
189 } // namespace gmx
190
191 #endif