Use proper doxygen tags in modular simulator
[gromacs.git] / src / gromacs / modularsimulator / signallers.h
blob2c74ab9454d1b92a3f8bf69fef552a2b70acd5a4
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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 /*! \internal \file
36 * \brief Declares the signallers for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
41 * This header is only used within the modular simulator module
44 #ifndef GMX_MODULARSIMULATOR_SIGNALLERS_H
45 #define GMX_MODULARSIMULATOR_SIGNALLERS_H
47 #include <vector>
49 #include "gromacs/compat/pointers.h"
51 #include "modularsimulatorinterfaces.h"
53 namespace gmx
55 class StopHandler;
56 class TrajectoryElement;
58 /*! \internal
59 * \ingroup module_modularsimulator
60 * \brief Builder for signallers
62 * This builder allows clients to register, and then builds the signaller
63 * passing on the list of clients.
65 * \tparam Signaller The signaller to be built
67 template<typename Signaller>
68 class SignallerBuilder final
70 public:
71 //! Allows clients to register to the signaller
72 void registerSignallerClient(typename Signaller::Client* client);
74 //! Build the signaller
75 template<typename... Args>
76 std::unique_ptr<Signaller> build(Args&&... args);
78 private:
79 //! List of signaller clients
80 std::vector<typename Signaller::Client*> signallerClients_;
81 //! The state of the builder
82 ModularSimulatorBuilderState state_ = ModularSimulatorBuilderState::AcceptingClientRegistrations;
84 //! Helper function to get the callbacks from the clients
85 template<typename... Args>
86 std::vector<SignallerCallback> buildCallbackVector(Args&&... args);
88 /*! \brief Get a callback from a single client
90 * This is in a separate function, as the exact call depends on the
91 * specific signaller / client.
93 template<typename... Args>
94 std::optional<SignallerCallback> getSignallerCallback(typename Signaller::Client* client,
95 Args&&... args);
98 /*! \internal
99 * \ingroup module_modularsimulator
100 * \brief Element signalling a neighbor search step
102 * This element informs its clients via callbacks
103 * when a neighbor-searching step is happening.
105 class NeighborSearchSignaller final : public ISignaller
107 public:
108 /*! \brief Run the signaller at a specific step / time
110 * Informs callbacks if step % nstlist_ == 0
112 * \param step The current time step
113 * \param time The current time
115 void signal(Step step, Time time) override;
117 //! Do nothing at setup time
118 void setup() override{};
120 //! Allow builder to construct
121 friend class SignallerBuilder<NeighborSearchSignaller>;
122 //! Define client type
123 typedef INeighborSearchSignallerClient Client;
125 private:
126 /*! \brief Constructor
128 * \param callbacks A vector of pointers to callbacks
129 * \param nstlist The frequency at which neighbor search is performed
130 * \param initStep The first step of the simulation
131 * \param initTime The start time of the simulation
133 NeighborSearchSignaller(std::vector<SignallerCallback> callbacks, Step nstlist, Step initStep, Time initTime);
135 //! Client callbacks
136 std::vector<SignallerCallback> callbacks_;
138 //! The NS frequency
139 const Step nstlist_;
140 //! The initial step of the simulation
141 const Step initStep_;
142 //! The initial time of the simulation
143 const Time initTime_;
146 /*! \internal
147 * \ingroup module_modularsimulator
148 * \brief Element signalling the last step
150 * This element informs its clients via callbacks
151 * when the last step is happening.
153 class LastStepSignaller final : public ISignaller, public INeighborSearchSignallerClient
155 public:
156 /*! \brief Run the signaller at a specific step / time
158 * Informs callbacks if this is the last step
160 * \param step The current time step
161 * \param time The current time
163 void signal(Step step, Time time) override;
165 //! Check that necessary registration was done
166 void setup() override;
168 //! Allow builder to construct
169 friend class SignallerBuilder<LastStepSignaller>;
170 //! Define client type
171 typedef ILastStepSignallerClient Client;
173 private:
174 /*! \brief Constructor
176 * \param callbacks A vector of pointers to callbacks
177 * \param nsteps The total number of steps for the simulation
178 * \param initStep The first step of the simulation
179 * \param stopHandler A pointer to the stop handler (LastStepSignaller takes ownership)
181 LastStepSignaller(std::vector<SignallerCallback> callbacks, Step nsteps, Step initStep, StopHandler* stopHandler);
183 //! Client callbacks
184 std::vector<SignallerCallback> callbacks_;
186 //! The last step of the simulation
187 const Step stopStep_;
188 //! Whether we signalled last step due to stop condition
189 bool signalledStopCondition_;
190 //! A pointer to the stop handler communicating signal and time-related stops
191 StopHandler* stopHandler_;
193 //! INeighborSearchSignallerClient implementation
194 std::optional<SignallerCallback> registerNSCallback() override;
195 //! The next NS step (notified by NS signaller)
196 Step nextNSStep_;
197 //! Whether we registered to the NS signaller
198 bool nsStepRegistrationDone_;
201 /*! \internal
202 * \ingroup module_modularsimulator
203 * \brief Element signalling a logging step
205 * This element informs its clients via callbacks
206 * when a logging step is happening.
208 class LoggingSignaller final : public ISignaller, public ILastStepSignallerClient
210 public:
211 /*! \brief Run the signaller at a specific step / time
213 * Informs callbacks if step % nstlog_ == 0
215 * \param step The current time step
216 * \param time The current time
218 void signal(Step step, Time time) override;
220 //! Check that necessary registration was done
221 void setup() override;
223 //! Allow builder to construct
224 friend class SignallerBuilder<LoggingSignaller>;
225 //! Define client type
226 typedef ILoggingSignallerClient Client;
228 private:
229 /*! \brief Constructor
231 * \param callbacks A vector of pointers to callbacks
232 * \param nstlog The logging frequency
233 * \param initStep The first step of the simulation
234 * \param initTime The start time of the simulation
236 LoggingSignaller(std::vector<SignallerCallback> callbacks, Step nstlog, Step initStep, Time initTime);
238 //! Client callbacks
239 std::vector<SignallerCallback> callbacks_;
241 //! The logging frequency
242 const Step nstlog_;
243 //! The initial step of the simulation
244 const Step initStep_;
245 //! The initial time of the simulation
246 const Time initTime_;
248 //! ILastStepSignallerClient implementation
249 std::optional<SignallerCallback> registerLastStepCallback() override;
250 //! The last step (notified by signaller)
251 Step lastStep_;
252 //! Whether we registered to the last step signaller
253 bool lastStepRegistrationDone_;
256 /*! \internal
257 * \ingroup module_modularsimulator
258 * \brief Element signalling trajectory writing
260 * During signalling phase, it checks whether the current step is a writing
261 * step for either the energy or the state (position, velocity, forces)
262 * trajectory. It then notifies the signaller clients of the upcoming step.
264 * The TrajectorySignaller works in close collaboration with the TrajectoryElement
265 * which does the actual trajectory writing during the simulation step.
267 class TrajectorySignaller final : public ISignaller, public ILastStepSignallerClient
269 public:
270 /*! \brief Prepare signaller
272 * Check that necessary registration was done
274 void setup() override;
276 /*! \brief Run the signaller at a specific step / time
278 * Informs clients when energy or state will be written.
280 * \param step The current time step
281 * \param time The current time
283 void signal(Step step, Time time) override;
285 //! Allow builder to construct
286 friend class SignallerBuilder<TrajectorySignaller>;
287 //! Define client type
288 typedef ITrajectorySignallerClient Client;
290 private:
291 //! Constructor
292 TrajectorySignaller(std::vector<SignallerCallback> signalEnergyCallbacks,
293 std::vector<SignallerCallback> signalStateCallbacks,
294 int nstxout,
295 int nstvout,
296 int nstfout,
297 int nstxoutCompressed,
298 int tngBoxOut,
299 int tngLambdaOut,
300 int tngBoxOutCompressed,
301 int tngLambdaOutCompressed,
302 int nstenergy);
304 //! Output frequencies
305 //! {
306 const int nstxout_;
307 const int nstvout_;
308 const int nstfout_;
309 const int nstxoutCompressed_;
310 const int tngBoxOut_;
311 const int tngLambdaOut_;
312 const int tngBoxOutCompressed_;
313 const int tngLambdaOutCompressed_;
314 const int nstenergy_;
315 //! }
317 //! Callbacks to signal events
318 //! {
319 std::vector<SignallerCallback> signalEnergyCallbacks_;
320 std::vector<SignallerCallback> signalStateCallbacks_;
321 //! }
324 * Last step client
326 Step lastStep_;
327 bool lastStepRegistrationDone_;
328 //! ILastStepSignallerClient implementation
329 std::optional<SignallerCallback> registerLastStepCallback() override;
332 /*! \internal
333 * \ingroup module_modularsimulator
334 * \brief Element signalling energy related special steps
336 * This element informs its clients via callbacks
337 * of the following events:
338 * - energy calculation step
339 * - virial calculation step
340 * - free energy calculation step
342 class EnergySignaller final : public ISignaller, public ITrajectorySignallerClient, public ILoggingSignallerClient
344 public:
345 /*! \brief Run the signaller at a specific step / time
347 * Informs callbacks of energy / virial / free energy special steps
349 * \param step The current time step
350 * \param time The current time
352 void signal(Step step, Time time) override;
354 //! Check that necessary registration was done
355 void setup() override;
357 //! Allow builder to construct
358 friend class SignallerBuilder<EnergySignaller>;
359 //! Define client type
360 typedef IEnergySignallerClient Client;
362 private:
363 /*! \brief Constructor
365 * \param calculateEnergyCallbacks A vector of pointers to callbacks (energy steps)
366 * \param calculateVirialCallbacks A vector of pointers to callbacks (virial steps)
367 * \param calculateFreeEnergyCallbacks A vector of pointers to callbacks (free energy steps)
368 * \param nstcalcenergy The energy calculation frequency
369 * \param nstcalcfreeenergy The free energy calculation frequency
370 * \param nstcalcvirial The free energy calculation frequency
372 EnergySignaller(std::vector<SignallerCallback> calculateEnergyCallbacks,
373 std::vector<SignallerCallback> calculateVirialCallbacks,
374 std::vector<SignallerCallback> calculateFreeEnergyCallbacks,
375 int nstcalcenergy,
376 int nstcalcfreeenergy,
377 int nstcalcvirial);
379 //! Client callbacks
380 //! {
381 std::vector<SignallerCallback> calculateEnergyCallbacks_;
382 std::vector<SignallerCallback> calculateVirialCallbacks_;
383 std::vector<SignallerCallback> calculateFreeEnergyCallbacks_;
384 //! }
386 //! The energy calculation frequency
387 const int nstcalcenergy_;
388 //! The free energy calculation frequency
389 const int nstcalcfreeenergy_;
390 //! The virial calculation frequency
391 const int nstcalcvirial_;
393 //! ITrajectorySignallerClient implementation
394 std::optional<SignallerCallback> registerTrajectorySignallerCallback(TrajectoryEvent event) override;
395 //! The energy writing step (notified by signaller)
396 Step energyWritingStep_;
397 //! Whether we registered to the trajectory signaller
398 bool trajectoryRegistrationDone_;
400 //! ILoggingSignallerClient implementation
401 std::optional<SignallerCallback> registerLoggingCallback() override;
402 //! The next logging step (notified by signaller)
403 Step loggingStep_;
404 //! Whether we registered to the logging signaller
405 bool loggingRegistrationDone_;
408 //! Allows clients to register to the signaller
409 template<class Signaller>
410 void SignallerBuilder<Signaller>::registerSignallerClient(typename Signaller::Client* client)
412 if (client)
414 if (state_ == ModularSimulatorBuilderState::NotAcceptingClientRegistrations)
416 throw SimulationAlgorithmSetupError(
417 "Tried to register to signaller after it was built.");
419 signallerClients_.emplace_back(client);
423 /*! \brief Build the signaller
425 * General version - for NeighborSearchSignaller, LastStepSignaller, LoggingSignaller
427 template<class Signaller>
428 template<typename... Args>
429 std::unique_ptr<Signaller> SignallerBuilder<Signaller>::build(Args&&... args)
431 state_ = ModularSimulatorBuilderState::NotAcceptingClientRegistrations;
432 auto callbacks = buildCallbackVector();
433 // NOLINTNEXTLINE(modernize-make-unique): make_unique does not work with private constructor
434 return std::unique_ptr<Signaller>(new Signaller(std::move(callbacks), std::forward<Args>(args)...));
437 /*! \brief Build the signaller
439 * Specialized version - TrajectorySignaller has a different build process
441 template<>
442 template<typename... Args>
443 std::unique_ptr<TrajectorySignaller> SignallerBuilder<TrajectorySignaller>::build(Args&&... args)
445 state_ = ModularSimulatorBuilderState::NotAcceptingClientRegistrations;
446 auto signalEnergyCallbacks = buildCallbackVector(TrajectoryEvent::EnergyWritingStep);
447 auto signalStateCallbacks = buildCallbackVector(TrajectoryEvent::StateWritingStep);
448 // NOLINTNEXTLINE(modernize-make-unique): make_unique does not work with private constructor
449 return std::unique_ptr<TrajectorySignaller>(new TrajectorySignaller(
450 std::move(signalEnergyCallbacks), std::move(signalStateCallbacks), std::forward<Args>(args)...));
453 /*! \brief Build the signaller
455 * Specialized version - EnergySignaller has a significantly different build process
457 template<>
458 template<typename... Args>
459 std::unique_ptr<EnergySignaller> SignallerBuilder<EnergySignaller>::build(Args&&... args)
461 state_ = ModularSimulatorBuilderState::NotAcceptingClientRegistrations;
462 auto calculateEnergyCallbacks = buildCallbackVector(EnergySignallerEvent::EnergyCalculationStep);
463 auto calculateVirialCallbacks = buildCallbackVector(EnergySignallerEvent::VirialCalculationStep);
464 auto calculateFreeEnergyCallbacks =
465 buildCallbackVector(EnergySignallerEvent::FreeEnergyCalculationStep);
466 // NOLINTNEXTLINE(modernize-make-unique): make_unique does not work with private constructor
467 return std::unique_ptr<EnergySignaller>(new EnergySignaller(
468 std::move(calculateEnergyCallbacks), std::move(calculateVirialCallbacks),
469 std::move(calculateFreeEnergyCallbacks), std::forward<Args>(args)...));
472 //! Helper function to get the callbacks from the clients
473 template<typename Signaller>
474 template<typename... Args>
475 std::vector<SignallerCallback> SignallerBuilder<Signaller>::buildCallbackVector(Args&&... args)
477 std::vector<SignallerCallback> callbacks;
478 // Allow clients to register their callbacks
479 for (auto& client : signallerClients_)
481 if (auto callback = getSignallerCallback(client, std::forward<Args>(args)...)) // don't register nullptr
483 callbacks.emplace_back(std::move(*callback));
486 return callbacks;
489 //! Get a callback from a single client - NeighborSearchSignaller
490 template<>
491 template<typename... Args>
492 std::optional<SignallerCallback> SignallerBuilder<NeighborSearchSignaller>::getSignallerCallback(
493 typename NeighborSearchSignaller::Client* client,
494 Args&&... args)
496 return client->registerNSCallback(std::forward<Args>(args)...);
499 //! Get a callback from a single client - LastStepSignaller
500 template<>
501 template<typename... Args>
502 std::optional<SignallerCallback>
503 SignallerBuilder<LastStepSignaller>::getSignallerCallback(typename LastStepSignaller::Client* client,
504 Args&&... args)
506 return client->registerLastStepCallback(std::forward<Args>(args)...);
509 //! Get a callback from a single client - LoggingSignaller
510 template<>
511 template<typename... Args>
512 std::optional<SignallerCallback>
513 SignallerBuilder<LoggingSignaller>::getSignallerCallback(typename LoggingSignaller::Client* client,
514 Args&&... args)
516 return client->registerLoggingCallback(std::forward<Args>(args)...);
519 //! Get a callback from a single client - TrajectorySignaller
520 template<>
521 template<typename... Args>
522 std::optional<SignallerCallback>
523 SignallerBuilder<TrajectorySignaller>::getSignallerCallback(typename TrajectorySignaller::Client* client,
524 Args&&... args)
526 return client->registerTrajectorySignallerCallback(std::forward<Args>(args)...);
529 //! Get a callback from a single client - EnergySignaller
530 template<>
531 template<typename... Args>
532 std::optional<SignallerCallback>
533 SignallerBuilder<EnergySignaller>::getSignallerCallback(typename EnergySignaller::Client* client,
534 Args&&... args)
536 return client->registerEnergyCallback(std::forward<Args>(args)...);
539 } // namespace gmx
541 #endif // GMX_MODULARSIMULATOR_SIGNALLERS_H