Eliminate mdlib/mdrun.h
[alexxy/gromacs.git] / src / gromacs / mdrun / integrator.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,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 integrator interface for mdrun
37  *
38  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \ingroup module_mdrun
41  */
42 #ifndef GMX_MDRUN_INTEGRATOR_H
43 #define GMX_MDRUN_INTEGRATOR_H
44
45 #include <cstdio>
46
47 #include <memory>
48
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/real.h"
51
52 class energyhistory_t;
53 struct gmx_enfrot;
54 struct gmx_mtop_t;
55 struct gmx_membed_t;
56 struct gmx_multisim_t;
57 struct gmx_output_env_t;
58 struct gmx_vsite_t;
59 struct gmx_wallcycle;
60 struct gmx_walltime_accounting;
61 struct ObservablesHistory;
62 struct ReplicaExchangeParameters;
63 struct t_commrec;
64 struct t_fcdata;
65 struct t_forcerec;
66 struct t_filenm;
67 struct t_inputrec;
68 struct t_nrnb;
69 class t_state;
70
71 namespace gmx
72 {
73
74 class BoxDeformation;
75 class Constraints;
76 class PpForceWorkload;
77 class IMDOutputProvider;
78 class MDLogger;
79 class MDAtoms;
80 class StopHandlerBuilder;
81 struct MdrunOptions;
82
83 //! Function type for integrator code.
84 using IntegratorFunctionType = void();
85
86 /*! \internal
87  * \brief Struct to handle setting up and running the different "integrators".
88  *
89  * This struct is a mere aggregate of parameters to pass to evaluate an
90  * energy, so that future changes to names and types of them consume
91  * less time when refactoring other code.
92  *
93  * Aggregate initialization is used, for which the chief risk is that
94  * if a member is added at the end and not all initializer lists are
95  * updated, then the member will be value initialized, which will
96  * typically mean initialization to zero.
97  *
98  * Having multiple integrators as member functions isn't a good
99  * design, and we definitely only intend one to be called, but the
100  * goal is to make it easy to change the names and types of members
101  * without having to make identical changes in several places in the
102  * code. Once many of them have become modules, we should change this
103  * approach.
104  *
105  * Note that the presence of const reference members means that the
106  * default constructor would be implicitly deleted. But we only want
107  * to make one of these when we know how to initialize these members,
108  * so that is perfect. To ensure this remains true even if we would
109  * remove those members, we explicitly delete this constructor.
110  * Other constructors, copies and moves are OK. */
111 struct Integrator
112 {
113     //! Handles logging.
114     FILE                               *fplog;
115     //! Handles communication.
116     t_commrec                          *cr;
117     //! Coordinates multi-simulations.
118     const gmx_multisim_t               *ms;
119     //! Handles logging.
120     const MDLogger                     &mdlog;
121     //! Count of input file options.
122     int                                 nfile;
123     //! Content of input file options.
124     const t_filenm                     *fnm;
125     //! Handles writing text output.
126     const gmx_output_env_t             *oenv;
127     //! Contains command-line options to mdrun.
128     const MdrunOptions                 &mdrunOptions;
129     //! Handles virtual sites.
130     gmx_vsite_t                        *vsite;
131     //! Handles constraints.
132     Constraints                        *constr;
133     //! Handles enforced rotation.
134     gmx_enfrot                         *enforcedRotation;
135     //! Handles box deformation.
136     BoxDeformation                     *deform;
137     //! Handles writing output files.
138     IMDOutputProvider                  *outputProvider;
139     //! Contains user input mdp options.
140     t_inputrec                         *inputrec;
141     //! Full system topology.
142     gmx_mtop_t                         *top_global;
143     //! Helper struct for force calculations.
144     t_fcdata                           *fcd;
145     //! Full simulation state (only non-nullptr on master rank).
146     t_state                            *state_global;
147     //! History of simulation observables.
148     ObservablesHistory                 *observablesHistory;
149     //! Atom parameters for this domain.
150     MDAtoms                            *mdAtoms;
151     //! Manages flop accounting.
152     t_nrnb                             *nrnb;
153     //! Manages wall cycle accounting.
154     gmx_wallcycle                      *wcycle;
155     //! Parameters for force calculations.
156     t_forcerec                         *fr;
157     //! Schedule of force-calculation work each step for this task.
158     PpForceWorkload                    *ppForceWorkload;
159     //! Parameters for replica exchange algorihtms.
160     const ReplicaExchangeParameters    &replExParams;
161     //! Parameters for membrane embedding.
162     gmx_membed_t                       *membed;
163     //! Manages wall time accounting.
164     gmx_walltime_accounting            *walltime_accounting;
165     //! Registers stop conditions
166     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder;
167     //! Implements the normal MD integrators.
168     IntegratorFunctionType              do_md;
169     //! Implements the rerun functionality.
170     IntegratorFunctionType              do_rerun;
171     //! Implements steepest descent EM.
172     IntegratorFunctionType              do_steep;
173     //! Implements conjugate gradient energy minimization
174     IntegratorFunctionType              do_cg;
175     //! Implements onjugate gradient energy minimization using the L-BFGS algorithm
176     IntegratorFunctionType              do_lbfgs;
177     //! Implements normal mode analysis
178     IntegratorFunctionType              do_nm;
179     //! Implements test particle insertion
180     IntegratorFunctionType              do_tpi;
181     //! Implements MiMiC QM/MM workflow
182     IntegratorFunctionType              do_mimic;
183     /*! \brief Function to run the correct IntegratorFunctionType,
184      * based on the .mdp integrator field. */
185     void run(unsigned int ei, bool doRerun);
186     //! We only intend to construct such objects with an initializer list.
187     Integrator() = delete;
188 };
189
190 }      // namespace gmx
191
192 #endif // GMX_MDRUN_INTEGRATOR_H