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