Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / modularsimulator / computeglobalselement.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 /*! \libinternal \file
36  * \brief Declares the global reduction element for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  */
41
42 #ifndef GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H
43 #define GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H
44
45 #include "gromacs/mdlib/simulationsignal.h"
46 #include "gromacs/mdlib/vcm.h"
47
48 #include "energyelement.h"
49 #include "modularsimulatorinterfaces.h"
50 #include "statepropagatordata.h"
51 #include "topologyholder.h"
52
53 struct gmx_global_stat;
54 struct gmx_wallcycle;
55 struct t_nrnb;
56
57 namespace gmx
58 {
59 class FreeEnergyPerturbationElement;
60 class MDAtoms;
61 class MDLogger;
62
63 //! \addtogroup module_modularsimulator
64 //! \{
65
66 //! The different global reduction schemes we know about
67 enum class ComputeGlobalsAlgorithm
68 {
69     LeapFrog,
70     VelocityVerletAtFullTimeStep,
71     VelocityVerletAfterCoordinateUpdate
72 };
73
74 //! The function type allowing to request a check of the number of bonded interactions
75 typedef std::function<void()> CheckBondedInteractionsCallback;
76 //! Pointer to the function type allowing to request a check of the number of bonded interactions
77 typedef std::unique_ptr<CheckBondedInteractionsCallback> CheckBondedInteractionsCallbackPtr;
78
79 /*! \libinternal
80  * \brief Encapsulate the calls to `compute_globals`
81  *
82  * This element aims at offering an interface to the legacy
83  * implementation which is compatible with the new simulator approach.
84  *
85  * The element comes in 3 (templated) flavors: the leap-frog case, the first
86  * call during a velocity-verlet integrator, and the second call during a
87  * velocity-verlet integrator. In velocity verlet, the state at the beginning
88  * of the step corresponds to
89  *     positions at time t
90  *     velocities at time t - dt/2
91  * The first velocity propagation (+dt/2) therefore actually corresponds to the
92  * previous step, bringing the state to the full timestep at time t. Most global
93  * reductions are made at this point. The second call is needed to correct the
94  * constraint virial after the second propagation of velocities (+dt/2) and of
95  * the positions (+dt).
96  *
97  * @tparam algorithm  The global reduction scheme
98  */
99 template<ComputeGlobalsAlgorithm algorithm>
100 class ComputeGlobalsElement final :
101     public ISimulatorElement,
102     public IEnergySignallerClient,
103     public ITrajectorySignallerClient,
104     public ITopologyHolderClient
105 {
106 public:
107     //! Constructor
108     ComputeGlobalsElement(StatePropagatorData*           statePropagatorData,
109                           EnergyElement*                 energyElement,
110                           FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
111                           SimulationSignals*             signals,
112                           int                            nstglobalcomm,
113                           FILE*                          fplog,
114                           const MDLogger&                mdlog,
115                           t_commrec*                     cr,
116                           t_inputrec*                    inputrec,
117                           const MDAtoms*                 mdAtoms,
118                           t_nrnb*                        nrnb,
119                           gmx_wallcycle*                 wcycle,
120                           t_forcerec*                    fr,
121                           const gmx_mtop_t*              global_top,
122                           Constraints*                   constr,
123                           bool                           hasReadEkinState);
124
125     //! Destructor
126     ~ComputeGlobalsElement() override;
127
128     /*! \brief Element setup - first call to compute_globals
129      *
130      */
131     void elementSetup() override;
132
133     /*! \brief Register run function for step / time
134      *
135      * This registers the call to compute_globals when needed.
136      *
137      * @param step                 The step number
138      * @param time                 The time
139      * @param registerRunFunction  Function allowing to register a run function
140      */
141     void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
142
143     //! Get callback to request checking of bonded interactions
144     CheckBondedInteractionsCallbackPtr getCheckNumberOfBondedInteractionsCallback();
145
146     //! No element teardown needed
147     void elementTeardown() override {}
148
149 private:
150     //! ITopologyClient implementation
151     void setTopology(const gmx_localtop_t* top) override;
152     //! IEnergySignallerClient implementation
153     SignallerCallbackPtr registerEnergyCallback(EnergySignallerEvent event) override;
154     //! ITrajectorySignallerClient implementation
155     SignallerCallbackPtr registerTrajectorySignallerCallback(TrajectoryEvent event) override;
156     //! The compute_globals call
157     void compute(Step step, unsigned int flags, SimulationSignaller* signaller, bool useLastBox, bool isInit = false);
158
159     //! Next step at which energy needs to be reduced
160     Step energyReductionStep_;
161     //! Next step at which virial needs to be reduced
162     Step virialReductionStep_;
163
164     //! Whether center of mass motion stopping is enabled
165     const bool doStopCM_;
166     //! Number of steps after which center of mass motion is removed
167     int nstcomm_;
168     //! Compute globals communication period
169     int nstglobalcomm_;
170     //! The last (planned) step (only used for LF)
171     const Step lastStep_;
172     //! The initial step (only used for VV)
173     const Step initStep_;
174     //! A dummy signaller (used for setup and VV)
175     std::unique_ptr<SimulationSignaller> nullSignaller_;
176     //! Whether we read kinetic energy from checkpoint
177     const bool hasReadEkinState_;
178
179     /*! \brief Check that DD doesn't miss bonded interactions
180      *
181      * Domain decomposition could incorrectly miss a bonded
182      * interaction, but checking for that requires a global
183      * communication stage, which does not otherwise happen in DD
184      * code. So we do that alongside the first global energy reduction
185      * after a new DD is made. These variables handle whether the
186      * check happens, and the result it returns.
187      */
188     //! @{
189     int  totalNumberOfBondedInteractions_;
190     bool shouldCheckNumberOfBondedInteractions_;
191     //! @}
192
193     /*! \brief Signal to ComputeGlobalsElement that it should check for DD errors
194      *
195      * Note that this should really be the responsibility of the DD element.
196      * MDLogger, global and local topology are only needed due to the call to
197      * checkNumberOfBondedInteractions(...).
198      *
199      * The DD element should have a single variable which gets reduced, and then
200      * be responsible for the checking after a global reduction has happened.
201      * This would, however, require a new approach for the compute_globals calls,
202      * which is not yet implemented. So for now, we're leaving this here.
203      */
204     void needToCheckNumberOfBondedInteractions();
205
206     //! Global reduction struct
207     gmx_global_stat* gstat_;
208
209     //! Pointer to the microstate
210     StatePropagatorData* statePropagatorData_;
211     //! Pointer to the energy element (needed for the tensors and mu_tot)
212     EnergyElement* energyElement_;
213     //! Pointer to the local topology (only needed for checkNumberOfBondedInteractions)
214     const gmx_localtop_t* localTopology_;
215     //! Pointer to the free energy perturbation element
216     FreeEnergyPerturbationElement* freeEnergyPerturbationElement_;
217
218     //! Center of mass motion removal
219     t_vcm vcm_;
220     //! Signals
221     SimulationSignals* signals_;
222
223     // Access to ISimulator data
224     //! Handles logging.
225     FILE* fplog_;
226     //! Handles logging.
227     const MDLogger& mdlog_;
228     //! Handles communication.
229     t_commrec* cr_;
230     //! Contains user input mdp options.
231     t_inputrec* inputrec_;
232     //! Full system topology - only needed for checkNumberOfBondedInteractions.
233     const gmx_mtop_t* top_global_;
234     //! Atom parameters for this domain.
235     const MDAtoms* mdAtoms_;
236     //! Handles constraints.
237     Constraints* constr_;
238     //! Manages flop accounting.
239     t_nrnb* nrnb_;
240     //! Manages wall cycle accounting.
241     gmx_wallcycle* wcycle_;
242     //! Parameters for force calculations.
243     t_forcerec* fr_;
244 };
245
246 //! \}
247 } // namespace gmx
248
249 #endif // GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H