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