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 /*! \defgroup module_modularsimulator The modular simulator
36 * \ingroup group_mdrun
38 * The modular simulator improves extensibility, adds Monte Carlo capabilities,
39 * promotes data locality and communication via interfaces, supports
40 * multi-stepping integrators, and paves the way for some task parallelism.
42 * For more information, see page_modularsimulator
43 * \todo Can we link to `docs/doxygen/lib/modularsimulator.md`?
45 * \author Pascal Merz <pascal.merz@me.com>
49 * Declares the main interfaces used by the modular simulator
51 * \author Pascal Merz <pascal.merz@me.com>
52 * \ingroup module_modularsimulator
54 * This header is only used within the modular simulator module
56 #ifndef GMX_MODULARSIMULATOR_MODULARSIMULATORINTERFACES_H
57 #define GMX_MODULARSIMULATOR_MODULARSIMULATORINTERFACES_H
63 #include "gromacs/math/vectypes.h"
64 #include "gromacs/mdtypes/checkpointdata.h"
65 #include "gromacs/utility/basedefinitions.h"
66 #include "gromacs/utility/exceptions.h"
68 struct gmx_localtop_t
;
77 template<class Signaller
>
78 class SignallerBuilder
;
79 class NeighborSearchSignaller
;
80 class LastStepSignaller
;
81 class LoggingSignaller
;
82 class TrajectorySignaller
;
83 class EnergySignaller
;
85 //! \addtogroup module_modularsimulator
88 // This will make signatures more legible
94 //! The function type that can be scheduled to be run during the simulator run
95 typedef std::function
<void()> SimulatorRunFunction
;
97 //! The function type that allows to register run functions
98 typedef std::function
<void(SimulatorRunFunction
)> RegisterRunFunction
;
101 * \brief The general interface for elements of the modular simulator
103 * Setup and teardown are run once at the beginning of the simulation
104 * (after all elements were set up, before the first step) and at the
105 * end of the simulation (after the last step, before elements are
106 * destructed), respectively. `registerRun` is periodically called
107 * during the run, at which point elements can decide to register one
108 * or more run functions to be run at a specific time / step. Registering
109 * more than one function is especially valuable for collective elements
110 * consisting of more than one single element.
112 class ISimulatorElement
115 /*! \brief Query whether element wants to run at step / time
117 * Element can register one or more functions to be run at that step through
118 * the registration pointer.
120 virtual void scheduleTask(Step
, Time
, const RegisterRunFunction
&) = 0;
121 //! Method guaranteed to be called after construction, before simulator run
122 virtual void elementSetup() = 0;
123 //! Method guaranteed to be called after simulator run, before deconstruction
124 virtual void elementTeardown() = 0;
125 //! Standard virtual destructor
126 virtual ~ISimulatorElement() = default;
130 * \brief The general Signaller interface
132 * Signallers are run at the beginning of Simulator steps, informing
133 * their clients about the upcoming step. This allows clients to
134 * decide if they need to be activated at this time step, and what
135 * functions they will have to run. Examples for signallers
136 * include the neighbor-search signaller (informs its clients when a
137 * neighbor-searching step is about to happen) or the logging
138 * signaller (informs its clients when a logging step is about to
141 * We expect all signallers to accept registration from clients, but
142 * different signallers might handle that differently, so we don't
143 * include this in the interface.
148 //! Function run before every step of scheduling
149 virtual void signal(Step
, Time
) = 0;
150 //! Method guaranteed to be called after construction, before simulator run
151 virtual void setup() = 0;
152 //! Standard virtual destructor
153 virtual ~ISignaller() = default;
156 //! The function type that can be registered to signallers for callback
157 typedef std::function
<void(Step
, Time
)> SignallerCallback
;
160 * \brief Interface for clients of the NeighborSearchSignaller
162 * Defining registerNSCallback allows clients to register an arbitrary callback
163 * for notification by the signaller.
165 class INeighborSearchSignallerClient
169 // (doxygen doesn't like these...)
170 //! Allow builder of NeighborSearchSignaller to ask for callback registration
171 friend class SignallerBuilder
<NeighborSearchSignaller
>;
173 //! Standard virtual destructor
174 virtual ~INeighborSearchSignallerClient() = default;
177 //! Return callback to NeighborSearchSignaller
178 virtual std::optional
<SignallerCallback
> registerNSCallback() = 0;
182 * \brief Interface for clients of the LastStepSignaller
184 * Defining registerLastStepCallback allows clients to register an arbitrary callback
185 * for notification by the signaller.
187 class ILastStepSignallerClient
191 // (doxygen doesn't like these...)
192 //! Allow builder of LastStepSignaller to ask for callback registration
193 friend class SignallerBuilder
<LastStepSignaller
>;
195 //! Standard virtual destructor
196 virtual ~ILastStepSignallerClient() = default;
199 //! Return callback to LastStepSignaller
200 virtual std::optional
<SignallerCallback
> registerLastStepCallback() = 0;
204 * \brief Interface for clients of the LoggingSignaller
206 * Defining registerLoggingCallback allows clients to register an arbitrary callback
207 * for notification by the signaller.
209 class ILoggingSignallerClient
213 // (doxygen doesn't like these...)
214 //! Allow builder of LoggingSignaller to ask for callback registration
215 friend class SignallerBuilder
<LoggingSignaller
>;
217 //! Standard virtual destructor
218 virtual ~ILoggingSignallerClient() = default;
221 //! Return callback to LoggingSignaller
222 virtual std::optional
<SignallerCallback
> registerLoggingCallback() = 0;
225 //! The energy events signalled by the EnergySignaller
226 enum class EnergySignallerEvent
228 EnergyCalculationStep
,
229 VirialCalculationStep
,
230 FreeEnergyCalculationStep
234 * \brief Interface for clients of the EnergySignaller
236 * Defining registerEnergyCallback allows clients to register an arbitrary callback
237 * for notification by the signaller for every EnergySignallerEvent separately.
239 class IEnergySignallerClient
243 // (doxygen doesn't like these...)
244 //! Allow builder of EnergySignaller to ask for callback registration
245 friend class SignallerBuilder
<EnergySignaller
>;
247 //! Standard virtual destructor
248 virtual ~IEnergySignallerClient() = default;
251 //! Return callback to EnergySignaller
252 virtual std::optional
<SignallerCallback
> registerEnergyCallback(EnergySignallerEvent
) = 0;
255 //! The trajectory writing events
256 enum class TrajectoryEvent
263 * \brief Interface for signaller clients of the TrajectoryElement
265 * Defining registerTrajectorySignallerCallback allows clients to register an arbitrary
266 * callback for notification by the signaller for every TrajectoryEvent separately.
268 class ITrajectorySignallerClient
272 // (doxygen doesn't like these...)
273 //! Allow builder of TrajectorySignaller to ask for callback registration
274 friend class SignallerBuilder
<TrajectorySignaller
>;
276 //! Standard virtual destructor
277 virtual ~ITrajectorySignallerClient() = default;
280 //! Return callback to TrajectoryElement
281 virtual std::optional
<SignallerCallback
> registerTrajectorySignallerCallback(TrajectoryEvent
) = 0;
284 /* Trajectory writing clients are handed a pointer to the output file handler,
285 * allowing them to write their own trajectory contribution.
287 * As trajectory writing and log writing cannot currently be separated for the
288 * energy, clients also get informed whether this is a trajectory-writing step
289 * and / or a log-writing step.
291 //! Function type for trajectory writing clients
292 typedef std::function
<void(gmx_mdoutf
*, Step
, Time
, bool, bool)> ITrajectoryWriterCallback
;
295 * \brief Interface for writer clients of the TrajectoryElement
297 * Defining registerTrajectoryWriterCallback allows clients to register an arbitrary
298 * callback called by the TrajectoryElement when trajectory writing happens.
300 * Setup and teardown methods allow clients to perform tasks with a valid output pointer.
302 class ITrajectoryWriterClient
306 // (doxygen doesn't like these...)
307 //! Allow TrajectoryElement to ask for callback registration
308 friend class TrajectoryElement
;
310 //! Standard virtual destructor
311 virtual ~ITrajectoryWriterClient() = default;
314 //! Setup method with valid output pointer.
315 virtual void trajectoryWriterSetup(gmx_mdoutf
* outf
) = 0;
316 //! Teardown method with valid output pointer.
317 virtual void trajectoryWriterTeardown(gmx_mdoutf
* outf
) = 0;
319 //! Return callback to TrajectoryElement
320 virtual std::optional
<ITrajectoryWriterCallback
> registerTrajectoryWriterCallback(TrajectoryEvent
) = 0;
324 * \brief Client requiring read access to the local topology
327 class ITopologyHolderClient
331 // (doxygen doesn't like these...)
332 //! Allow TopologyHolder to set new topology
333 friend class TopologyHolder
;
335 //! Standard virtual destructor
336 virtual ~ITopologyHolderClient() = default;
339 //! Pass pointer to new local topology
340 virtual void setTopology(const gmx_localtop_t
*) = 0;
344 * \brief Client that needs to store data during checkpointing
346 * Clients receive a CheckpointData object for reading and writing.
347 * Note that `ReadCheckpointData` is a typedef for
348 * `CheckpointData<CheckpointDataOperation::Read>`, and
349 * `WriteCheckpointData` is a typedef for
350 * `CheckpointData<CheckpointDataOperation::Write>`. This allows clients
351 * to write a single templated function, e.g.
352 * template<CheckpointDataOperation operation>
353 * void doCheckpointData(CheckpointData<operation>* checkpointData,
354 * const t_commrec* cr)
356 * checkpointData->scalar("important value", &value_);
358 * for both checkpoint reading and writing. This function can then be
359 * dispatched from the interface functions,
360 * void writeCheckpoint(WriteCheckpointData checkpointData, const t_commrec* cr)
362 * doCheckpointData<CheckpointDataOperation::Write>(&checkpointData, cr);
364 * void readCheckpoint(ReadCheckpointData checkpointData, const t_commrec* cr)
366 * doCheckpointData<CheckpointDataOperation::Read>(&checkpointData, cr);
368 * This reduces code duplication and ensures that reading and writing
369 * operations will not get out of sync.
371 class ICheckpointHelperClient
374 //! Standard virtual destructor
375 virtual ~ICheckpointHelperClient() = default;
377 //! Write checkpoint (CheckpointData object only passed on master rank)
378 virtual void saveCheckpointState(std::optional
<WriteCheckpointData
> checkpointData
,
379 const t_commrec
* cr
) = 0;
380 //! Read checkpoint (CheckpointData object only passed on master rank)
381 virtual void restoreCheckpointState(std::optional
<ReadCheckpointData
> checkpointData
,
382 const t_commrec
* cr
) = 0;
383 //! Get unique client id
384 [[nodiscard
]] virtual const std::string
& clientID() = 0;
388 * Exception class signalling that a requested element was not found.
392 class ElementNotFoundError final
: public ModularSimulatorError
395 //! \copydoc FileIOError::FileIOError()
396 explicit ElementNotFoundError(const ExceptionInitializer
& details
) :
397 ModularSimulatorError(details
)
403 * Exception class signalling that elements were not connected properly.
407 class MissingElementConnectionError final
: public ModularSimulatorError
410 //! \copydoc FileIOError::FileIOError()
411 explicit MissingElementConnectionError(const ExceptionInitializer
& details
) :
412 ModularSimulatorError(details
)
418 * Exception class signalling that the ModularSimulatorAlgorithm was set up
419 * in an incompatible way.
423 class SimulationAlgorithmSetupError final
: public ModularSimulatorError
426 //! \copydoc FileIOError::FileIOError()
427 explicit SimulationAlgorithmSetupError(const ExceptionInitializer
& details
) :
428 ModularSimulatorError(details
)
434 * Exception class signalling an error in reading or writing modular checkpoints.
438 class CheckpointError final
: public ModularSimulatorError
441 //! \copydoc FileIOError::FileIOError()
442 explicit CheckpointError(const ExceptionInitializer
& details
) : ModularSimulatorError(details
)
447 //! Enum allowing builders to store whether they can accept client registrations
448 enum class ModularSimulatorBuilderState
450 AcceptingClientRegistrations
,
451 NotAcceptingClientRegistrations
454 //! Generic callback to the propagator
455 typedef std::function
<void(Step
)> PropagatorCallback
;
458 * \brief Information needed to connect a propagator to a thermostat
460 struct PropagatorThermostatConnection
462 //! Function variable for setting velocity scaling variables.
463 std::function
<void(int)> setNumVelocityScalingVariables
;
464 //! Function variable for receiving view on velocity scaling.
465 std::function
<ArrayRef
<real
>()> getViewOnVelocityScaling
;
466 //! Function variable for callback.
467 std::function
<PropagatorCallback()> getVelocityScalingCallback
;
471 * \brief Information needed to connect a propagator to a barostat
473 struct PropagatorBarostatConnection
475 //! Function variable for receiving view on pressure scaling matrix.
476 std::function
<ArrayRef
<rvec
>()> getViewOnPRScalingMatrix
;
477 //! Function variable for callback.
478 std::function
<PropagatorCallback()> getPRScalingCallback
;
484 #endif // GMX_MODULARSIMULATOR_MODULARSIMULATORINTERFACES_H