Add class ListedForces
[alexxy/gromacs.git] / src / gromacs / mdrun / isimulator.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 /*! \internal
36  * \brief Declares the general simulator interface
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_mdrun
40  */
41 #ifndef GMX_MDRUN_ISIMULATOR_H
42 #define GMX_MDRUN_ISIMULATOR_H
43
44 #include "gromacs/mdlib/stophandler.h"
45
46 class energyhistory_t;
47 struct gmx_ekindata_t;
48 struct gmx_enerdata_t;
49 struct gmx_enfrot;
50 struct gmx_mtop_t;
51 struct gmx_membed_t;
52 struct gmx_multisim_t;
53 struct gmx_output_env_t;
54 struct gmx_wallcycle;
55 struct gmx_walltime_accounting;
56 struct ObservablesHistory;
57 struct pull_t;
58 struct ReplicaExchangeParameters;
59 struct t_commrec;
60 struct t_forcerec;
61 struct t_filenm;
62 struct t_inputrec;
63 struct t_nrnb;
64 struct t_swap;
65 class t_state;
66
67 namespace gmx
68 {
69 enum class StartingBehavior;
70 class BoxDeformation;
71 class Constraints;
72 class MdrunScheduleWorkload;
73 class IMDOutputProvider;
74 struct MdModulesNotifier;
75 class ImdSession;
76 class MDLogger;
77 class MDAtoms;
78 class StopHandlerBuilder;
79 struct MdrunOptions;
80 class VirtualSitesHandler;
81
82 /*! \internal
83  * \brief The Simulator interface
84  *
85  * This is the general interface for any type of simulation type
86  * ran with GROMACS. This allows having a builder return different
87  * Simulator objects based on user input.
88  */
89 class ISimulator
90 {
91 public:
92     /*! \brief The simulation run
93      *
94      * This will be called by the owner of the simulator object. To be redefined
95      * by the child classes. This function is expected to run the simulation.
96      */
97     virtual void run() = 0;
98     //! Standard destructor
99     virtual ~ISimulator() = default;
100     //! The constructor
101     ISimulator(FILE*                               fplog,
102                t_commrec*                          cr,
103                const gmx_multisim_t*               ms,
104                const MDLogger&                     mdlog,
105                int                                 nfile,
106                const t_filenm*                     fnm,
107                const gmx_output_env_t*             oenv,
108                const MdrunOptions&                 mdrunOptions,
109                StartingBehavior                    startingBehavior,
110                VirtualSitesHandler*                vsite,
111                Constraints*                        constr,
112                gmx_enfrot*                         enforcedRotation,
113                BoxDeformation*                     deform,
114                IMDOutputProvider*                  outputProvider,
115                const MdModulesNotifier&            mdModulesNotifier,
116                t_inputrec*                         inputrec,
117                ImdSession*                         imdSession,
118                pull_t*                             pull_work,
119                t_swap*                             swap,
120                gmx_mtop_t*                         top_global,
121                t_state*                            state_global,
122                ObservablesHistory*                 observablesHistory,
123                MDAtoms*                            mdAtoms,
124                t_nrnb*                             nrnb,
125                gmx_wallcycle*                      wcycle,
126                t_forcerec*                         fr,
127                gmx_enerdata_t*                     enerd,
128                gmx_ekindata_t*                     ekind,
129                MdrunScheduleWorkload*              runScheduleWork,
130                const ReplicaExchangeParameters&    replExParams,
131                gmx_membed_t*                       membed,
132                gmx_walltime_accounting*            walltime_accounting,
133                std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder,
134                bool                                doRerun) :
135         fplog(fplog),
136         cr(cr),
137         ms(ms),
138         mdlog(mdlog),
139         nfile(nfile),
140         fnm(fnm),
141         oenv(oenv),
142         mdrunOptions(mdrunOptions),
143         startingBehavior(startingBehavior),
144         vsite(vsite),
145         constr(constr),
146         enforcedRotation(enforcedRotation),
147         deform(deform),
148         outputProvider(outputProvider),
149         mdModulesNotifier(mdModulesNotifier),
150         inputrec(inputrec),
151         imdSession(imdSession),
152         pull_work(pull_work),
153         swap(swap),
154         top_global(top_global),
155         state_global(state_global),
156         observablesHistory(observablesHistory),
157         mdAtoms(mdAtoms),
158         nrnb(nrnb),
159         wcycle(wcycle),
160         fr(fr),
161         enerd(enerd),
162         ekind(ekind),
163         runScheduleWork(runScheduleWork),
164         replExParams(replExParams),
165         membed(membed),
166         walltime_accounting(walltime_accounting),
167         stopHandlerBuilder(std::move(stopHandlerBuilder)),
168         doRerun(doRerun)
169     {
170     }
171
172 protected:
173     //! Handles logging.
174     FILE* fplog;
175     //! Handles communication.
176     t_commrec* cr;
177     //! Coordinates multi-simulations.
178     const gmx_multisim_t* ms;
179     //! Handles logging.
180     const MDLogger& mdlog;
181     //! Count of input file options.
182     int nfile;
183     //! Content of input file options.
184     const t_filenm* fnm;
185     //! Handles writing text output.
186     const gmx_output_env_t* oenv;
187     //! Contains command-line options to mdrun.
188     const MdrunOptions& mdrunOptions;
189     //! Whether the simulation will start afresh, or restart with/without appending.
190     const StartingBehavior startingBehavior;
191     //! Handles virtual sites.
192     VirtualSitesHandler* vsite;
193     //! Handles constraints.
194     Constraints* constr;
195     //! Handles enforced rotation.
196     gmx_enfrot* enforcedRotation;
197     //! Handles box deformation.
198     BoxDeformation* deform;
199     //! Handles writing output files.
200     IMDOutputProvider* outputProvider;
201     //! Handles notifications to MdModules for checkpoint writing
202     const MdModulesNotifier& mdModulesNotifier;
203     //! Contains user input mdp options.
204     t_inputrec* inputrec;
205     //! The Interactive Molecular Dynamics session.
206     ImdSession* imdSession;
207     //! The pull work object.
208     pull_t* pull_work;
209     //! The coordinate-swapping session.
210     t_swap* swap;
211     //! Full system topology.
212     const gmx_mtop_t* top_global;
213     //! Full simulation state (only non-nullptr on master rank).
214     t_state* state_global;
215     //! History of simulation observables.
216     ObservablesHistory* observablesHistory;
217     //! Atom parameters for this domain.
218     MDAtoms* mdAtoms;
219     //! Manages flop accounting.
220     t_nrnb* nrnb;
221     //! Manages wall cycle accounting.
222     gmx_wallcycle* wcycle;
223     //! Parameters for force calculations.
224     t_forcerec* fr;
225     //! Data for energy output.
226     gmx_enerdata_t* enerd;
227     //! Kinetic energy data.
228     gmx_ekindata_t* ekind;
229     //! Schedule of work for each MD step for this task.
230     MdrunScheduleWorkload* runScheduleWork;
231     //! Parameters for replica exchange algorihtms.
232     const ReplicaExchangeParameters& replExParams;
233     //! Parameters for membrane embedding.
234     gmx_membed_t* membed;
235     //! Manages wall time accounting.
236     gmx_walltime_accounting* walltime_accounting;
237     //! Registers stop conditions
238     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder;
239     //! Whether we're doing a rerun.
240     bool doRerun;
241 };
242
243 } // namespace gmx
244
245 #endif // GMX_MDRUN_ISIMULATOR_H