f6905f20c747a31fa50ab2005c3170594d4c24b7
[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,2020,2021, 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 \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  * This header is only used within the modular simulator module
42  */
43
44 #ifndef GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H
45 #define GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H
46
47 #include "gromacs/mdlib/simulationsignal.h"
48 #include "gromacs/mdlib/vcm.h"
49
50 #include "energydata.h"
51 #include "modularsimulatorinterfaces.h"
52 #include "statepropagatordata.h"
53 #include "topologyholder.h"
54
55 struct gmx_global_stat;
56 struct gmx_wallcycle;
57 struct t_nrnb;
58
59 namespace gmx
60 {
61 class FreeEnergyPerturbationData;
62 class LegacySimulatorData;
63 class MDAtoms;
64 class MDLogger;
65 class ObservablesReducer;
66
67 //! \addtogroup module_modularsimulator
68 //! \{
69
70 //! The different global reduction schemes we know about
71 enum class ComputeGlobalsAlgorithm
72 {
73     LeapFrog,
74     VelocityVerlet
75 };
76
77 //! The function type allowing to request a check of the number of bonded interactions
78 typedef std::function<void()> CheckBondedInteractionsCallback;
79
80 /*! \internal
81  * \brief Encapsulate the calls to `compute_globals`
82  *
83  * This element aims at offering an interface to the legacy
84  * implementation which is compatible with the new simulator approach.
85  *
86  * The element comes in 3 (templated) flavors: the leap-frog case, the first
87  * call during a velocity-verlet integrator, and the second call during a
88  * velocity-verlet integrator. In velocity verlet, the state at the beginning
89  * of the step corresponds to
90  *     positions at time t
91  *     velocities at time t - dt/2
92  * The first velocity propagation (+dt/2) therefore actually corresponds to the
93  * previous step, bringing the state to the full timestep at time t. Most global
94  * reductions are made at this point. The second call is needed to correct the
95  * constraint virial after the second propagation of velocities (+dt/2) and of
96  * the positions (+dt).
97  *
98  * \tparam algorithm  The global reduction scheme
99  */
100 template<ComputeGlobalsAlgorithm algorithm>
101 class ComputeGlobalsElement final :
102     public ISimulatorElement,
103     public IEnergySignallerClient,
104     public ITrajectorySignallerClient,
105     public ITopologyHolderClient
106 {
107 public:
108     //! Constructor
109     ComputeGlobalsElement(StatePropagatorData*        statePropagatorData,
110                           EnergyData*                 energyData,
111                           FreeEnergyPerturbationData* freeEnergyPerturbationData,
112                           SimulationSignals*          signals,
113                           int                         nstglobalcomm,
114                           FILE*                       fplog,
115                           const MDLogger&             mdlog,
116                           t_commrec*                  cr,
117                           const t_inputrec*           inputrec,
118                           const MDAtoms*              mdAtoms,
119                           t_nrnb*                     nrnb,
120                           gmx_wallcycle*              wcycle,
121                           t_forcerec*                 fr,
122                           const gmx_mtop_t&           global_top,
123                           Constraints*                constr,
124                           ObservablesReducer*         observablesReducer);
125
126     //! Destructor
127     ~ComputeGlobalsElement() override;
128
129     /*! \brief Element setup - first call to compute_globals
130      *
131      */
132     void elementSetup() override;
133
134     /*! \brief Register run function for step / time
135      *
136      * This registers the call to compute_globals when needed.
137      *
138      * \param step                 The step number
139      * \param time                 The time
140      * \param registerRunFunction  Function allowing to register a run function
141      */
142     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
143
144     //! Get callback to request checking of bonded interactions
145     CheckBondedInteractionsCallback getCheckNumberOfBondedInteractionsCallback();
146
147     //! No element teardown needed
148     void elementTeardown() override {}
149
150     /*! \brief Factory method implementation
151      *
152      * \param legacySimulatorData  Pointer allowing access to simulator level data
153      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
154      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
155      * \param energyData  Pointer to the \c EnergyData object
156      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
157      * \param globalCommunicationHelper   Pointer to the \c GlobalCommunicationHelper object
158      * \param observablesReducer          Pointer to the \c ObservablesReducer object
159      *
160      * \throws std::bad_any_cast  on internal error in VelocityVerlet algorithm builder.
161      * \throws std::bad_alloc  when out of memory.
162      *
163      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
164      */
165     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
166                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
167                                                     StatePropagatorData*        statePropagatorData,
168                                                     EnergyData*                 energyData,
169                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
170                                                     GlobalCommunicationHelper* globalCommunicationHelper,
171                                                     ObservablesReducer*        observablesReducer);
172
173 private:
174     //! ITopologyClient implementation
175     void setTopology(const gmx_localtop_t* top) override;
176     //! IEnergySignallerClient implementation
177     std::optional<SignallerCallback> registerEnergyCallback(EnergySignallerEvent event) override;
178     //! ITrajectorySignallerClient implementation
179     std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
180     //! The compute_globals call
181     void compute(Step step, unsigned int flags, SimulationSignaller* signaller, bool useLastBox, bool isInit = false);
182
183     //! Next step at which energy needs to be reduced
184     Step energyReductionStep_;
185     //! Next step at which virial needs to be reduced
186     Step virialReductionStep_;
187
188     //! For VV only, we need to schedule twice per step. This keeps track of the scheduling stage.
189     Step vvSchedulingStep_;
190
191     //! Whether center of mass motion stopping is enabled
192     const bool doStopCM_;
193     //! Number of steps after which center of mass motion is removed
194     int nstcomm_;
195     //! Compute globals communication period
196     int nstglobalcomm_;
197     //! The last (planned) step (only used for LF)
198     const Step lastStep_;
199     //! The initial step (only used for VV)
200     const Step initStep_;
201     //! A dummy signaller (used for setup and VV)
202     std::unique_ptr<SimulationSignaller> nullSignaller_;
203
204     /*! \brief Check that DD doesn't miss bonded interactions
205      *
206      * Domain decomposition could incorrectly miss a bonded
207      * interaction, but checking for that requires a global
208      * communication stage, which does not otherwise happen in DD
209      * code. So we do that alongside the first global energy reduction
210      * after a new DD is made. These variables handle whether the
211      * check happens, and the result it returns.
212      */
213     //! \{
214     int  totalNumberOfBondedInteractions_;
215     bool shouldCheckNumberOfBondedInteractions_;
216     //! \}
217
218     /*! \brief Signal to ComputeGlobalsElement that it should check for DD errors
219      *
220      * Note that this should really be the responsibility of the DD element.
221      * MDLogger, global and local topology are only needed due to the call to
222      * checkNumberOfBondedInteractions(...).
223      *
224      * The DD element should have a single variable which gets reduced, and then
225      * be responsible for the checking after a global reduction has happened.
226      * This would, however, require a new approach for the compute_globals calls,
227      * which is not yet implemented. So for now, we're leaving this here.
228      */
229     void needToCheckNumberOfBondedInteractions();
230
231     //! Global reduction struct
232     gmx_global_stat* gstat_;
233
234     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
235     //! Pointer to the microstate
236     StatePropagatorData* statePropagatorData_;
237     //! Pointer to the energy data (needed for the tensors and mu_tot)
238     EnergyData* energyData_;
239     //! Pointer to the local topology (only needed for checkNumberOfBondedInteractions)
240     const gmx_localtop_t* localTopology_;
241     //! Pointer to the free energy perturbation data
242     FreeEnergyPerturbationData* freeEnergyPerturbationData_;
243
244     //! Center of mass motion removal
245     t_vcm vcm_;
246     //! Signals
247     SimulationSignals* signals_;
248
249     // Access to ISimulator data
250     //! Handles logging.
251     FILE* fplog_;
252     //! Handles logging.
253     const MDLogger& mdlog_;
254     //! Handles communication.
255     t_commrec* cr_;
256     //! Contains user input mdp options.
257     const t_inputrec* inputrec_;
258     //! Full system topology - only needed for checkNumberOfBondedInteractions.
259     const gmx_mtop_t& top_global_;
260     //! Atom parameters for this domain.
261     const MDAtoms* mdAtoms_;
262     //! Handles constraints.
263     Constraints* constr_;
264     //! Manages flop accounting.
265     t_nrnb* nrnb_;
266     //! Manages wall cycle accounting.
267     gmx_wallcycle* wcycle_;
268     //! Parameters for force calculations.
269     t_forcerec* fr_;
270     //! Coordinates reduction for observables
271     ObservablesReducer* observablesReducer_;
272 };
273
274 //! \}
275 } // namespace gmx
276
277 #endif // GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H