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