Remove unused declaration of function
[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 : public ISimulatorElement, public IEnergySignallerClient, public ITrajectorySignallerClient
102 {
103 public:
104     //! Constructor
105     ComputeGlobalsElement(StatePropagatorData*        statePropagatorData,
106                           EnergyData*                 energyData,
107                           FreeEnergyPerturbationData* freeEnergyPerturbationData,
108                           SimulationSignals*          signals,
109                           int                         nstglobalcomm,
110                           FILE*                       fplog,
111                           const MDLogger&             mdlog,
112                           t_commrec*                  cr,
113                           const t_inputrec*           inputrec,
114                           const MDAtoms*              mdAtoms,
115                           t_nrnb*                     nrnb,
116                           gmx_wallcycle*              wcycle,
117                           t_forcerec*                 fr,
118                           const gmx_mtop_t&           global_top,
119                           Constraints*                constr,
120                           ObservablesReducer*         observablesReducer);
121
122     //! Destructor
123     ~ComputeGlobalsElement() override;
124
125     /*! \brief Element setup - first call to compute_globals
126      *
127      */
128     void elementSetup() override;
129
130     /*! \brief Register run function for step / time
131      *
132      * This registers the call to compute_globals when needed.
133      *
134      * \param step                 The step number
135      * \param time                 The time
136      * \param registerRunFunction  Function allowing to register a run function
137      */
138     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
139
140     //! No element teardown needed
141     void elementTeardown() override {}
142
143     /*! \brief Factory method implementation
144      *
145      * \param legacySimulatorData  Pointer allowing access to simulator level data
146      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
147      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
148      * \param energyData  Pointer to the \c EnergyData object
149      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
150      * \param globalCommunicationHelper   Pointer to the \c GlobalCommunicationHelper object
151      * \param observablesReducer          Pointer to the \c ObservablesReducer object
152      *
153      * \throws std::bad_any_cast  on internal error in VelocityVerlet algorithm builder.
154      * \throws std::bad_alloc  when out of memory.
155      *
156      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
157      */
158     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
159                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
160                                                     StatePropagatorData*        statePropagatorData,
161                                                     EnergyData*                 energyData,
162                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
163                                                     GlobalCommunicationHelper* globalCommunicationHelper,
164                                                     ObservablesReducer*        observablesReducer);
165
166 private:
167     //! IEnergySignallerClient implementation
168     std::optional<SignallerCallback> registerEnergyCallback(EnergySignallerEvent event) override;
169     //! ITrajectorySignallerClient implementation
170     std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
171     //! The compute_globals call
172     void compute(Step step, unsigned int flags, SimulationSignaller* signaller, bool useLastBox, bool isInit = false);
173
174     //! Next step at which energy needs to be reduced
175     Step energyReductionStep_;
176     //! Next step at which virial needs to be reduced
177     Step virialReductionStep_;
178
179     //! For VV only, we need to schedule twice per step. This keeps track of the scheduling stage.
180     Step vvSchedulingStep_;
181
182     //! Whether center of mass motion stopping is enabled
183     const bool doStopCM_;
184     //! Number of steps after which center of mass motion is removed
185     int nstcomm_;
186     //! Compute globals communication period
187     int nstglobalcomm_;
188     //! The last (planned) step (only used for LF)
189     const Step lastStep_;
190     //! The initial step (only used for VV)
191     const Step initStep_;
192     //! A dummy signaller (used for setup and VV)
193     std::unique_ptr<SimulationSignaller> nullSignaller_;
194
195     /*! \brief Check that DD doesn't miss bonded interactions
196      *
197      * Domain decomposition could incorrectly miss a bonded
198      * interaction, but checking for that requires a global
199      * communication stage, which does not otherwise happen in DD
200      * code. So we do that alongside the first global energy reduction
201      * after a new DD is made. These variables handle whether the
202      * check happens, and the result it returns.
203      */
204     //! \{
205     int  totalNumberOfBondedInteractions_;
206     bool shouldCheckNumberOfBondedInteractions_;
207     //! \}
208
209     /*! \brief Signal to ComputeGlobalsElement that it should check for DD errors
210      *
211      * Note that this should really be the responsibility of the DD element.
212      * MDLogger, global and local topology are only needed due to the call to
213      * checkNumberOfBondedInteractions(...).
214      *
215      * The DD element should have a single variable which gets reduced, and then
216      * be responsible for the checking after a global reduction has happened.
217      * This would, however, require a new approach for the compute_globals calls,
218      * which is not yet implemented. So for now, we're leaving this here.
219      */
220     void needToCheckNumberOfBondedInteractions();
221
222     //! Global reduction struct
223     gmx_global_stat* gstat_;
224
225     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
226     //! Pointer to the microstate
227     StatePropagatorData* statePropagatorData_;
228     //! Pointer to the energy data (needed for the tensors and mu_tot)
229     EnergyData* energyData_;
230     //! Pointer to the free energy perturbation data
231     FreeEnergyPerturbationData* freeEnergyPerturbationData_;
232
233     //! Center of mass motion removal
234     t_vcm vcm_;
235     //! Signals
236     SimulationSignals* signals_;
237
238     // Access to ISimulator data
239     //! Handles logging.
240     FILE* fplog_;
241     //! Handles logging.
242     const MDLogger& mdlog_;
243     //! Handles communication.
244     t_commrec* cr_;
245     //! Contains user input mdp options.
246     const t_inputrec* inputrec_;
247     //! Full system topology - only needed for checkNumberOfBondedInteractions.
248     const gmx_mtop_t& top_global_;
249     //! Atom parameters for this domain.
250     const MDAtoms* mdAtoms_;
251     //! Handles constraints.
252     Constraints* constr_;
253     //! Manages flop accounting.
254     t_nrnb* nrnb_;
255     //! Manages wall cycle accounting.
256     gmx_wallcycle* wcycle_;
257     //! Parameters for force calculations.
258     t_forcerec* fr_;
259     //! Coordinates reduction for observables
260     ObservablesReducer* observablesReducer_;
261 };
262
263 //! \}
264 } // namespace gmx
265
266 #endif // GMX_MODULARSIMULATOR_COMPUTEGLOBALSELEMENT_H