Move responsibility for checking bondeds in reverse topology
[alexxy/gromacs.git] / src / gromacs / modularsimulator / simulatoralgorithm.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 Provides the modular simulator algorithm.
37  *
38  * Defines the ModularSimulatorAlgorithm class and its builder.
39  *
40  * \author Pascal Merz <pascal.merz@me.com>
41  * \ingroup module_modularsimulator
42  *
43  * This header is only used within the modular simulator module.
44  * Moving forward, the ModularSimulatorAlgorithmBuilder could be exposed to allow users to
45  * create custom algorithm - currently algorithms are only created an used by the ModularSimulator,
46  * meaning that this file is not exposed outside of the modular simulator module.
47  */
48 #ifndef GROMACS_MODULARSIMULATOR_SIMULATORALGORITHM_H
49 #define GROMACS_MODULARSIMULATOR_SIMULATORALGORITHM_H
50
51 #include <any>
52 #include <map>
53 #include <optional>
54 #include <string>
55 #include <typeinfo>
56
57 #include "gromacs/mdrun/isimulator.h"
58 #include "gromacs/mdtypes/state.h"
59 #include "gromacs/utility/exceptions.h"
60
61 #include "checkpointhelper.h"
62 #include "domdechelper.h"
63 #include "freeenergyperturbationdata.h"
64 #include "modularsimulatorinterfaces.h"
65 #include "pmeloadbalancehelper.h"
66 #include "signallers.h"
67 #include "topologyholder.h"
68 #include "trajectoryelement.h"
69
70 namespace gmx
71 {
72 enum class IntegrationStep;
73 class EnergyData;
74 class ModularSimulator;
75 class ResetHandler;
76 template<IntegrationStep algorithm>
77 class Propagator;
78 class TopologyHolder;
79
80 /*! \internal
81  * \ingroup module_modularsimulator
82  * \brief The modular simulator
83  *
84  * Based on the input given, this simulator builds independent elements and
85  * signallers and stores them in a respective vector. The run function
86  * runs the simulation by, in turn, building a task list from the elements
87  * for a predefined number of steps, then running the task list, and repeating
88  * until the stop criterion is fulfilled.
89  *
90  * The simulator algorithm is owning all elements involved in the simulation
91  * and is hence controlling their lifetime. This ensures that pointers and
92  * callbacks exchanged between elements remain valid throughout the duration
93  * of the simulation run.
94  */
95 class ModularSimulatorAlgorithm final
96 {
97 public:
98     //! Get next task in queue
99     [[nodiscard]] const SimulatorRunFunction* getNextTask();
100
101     // Only builder can construct
102     friend class ModularSimulatorAlgorithmBuilder;
103
104 private:
105     //! Constructor
106     ModularSimulatorAlgorithm(std::string              topologyName,
107                               FILE*                    fplog,
108                               t_commrec*               cr,
109                               const MDLogger&          mdlog,
110                               const MdrunOptions&      mdrunOptions,
111                               const t_inputrec*        inputrec,
112                               t_nrnb*                  nrnb,
113                               gmx_wallcycle*           wcycle,
114                               t_forcerec*              fr,
115                               gmx_walltime_accounting* walltime_accounting);
116
117     //! Calls setup for the simulator and all elements and signallers
118     void setup();
119     //! Calls teardown for the simulator and all elements
120     void teardown();
121     //! Re-populate task queue
122     void updateTaskQueue();
123
124     /*! \brief A function called once before the simulation run
125      *
126      * This function collects calls which need to be made once at the
127      * beginning of the simulation run. These are mainly printing information
128      * to the user and starting the wall time.
129      */
130     void simulatorSetup();
131
132     /*! \brief A function called once after the simulation run
133      *
134      * This function collects calls which need to be made once at the
135      * end of the simulation run.
136      */
137     void simulatorTeardown();
138
139     /*! \brief A function called before every step
140      *
141      * This function collects calls which need to be made before every
142      * simulation step. The need for this function could likely be
143      * eliminated by rewriting the stop and reset handler to fit the
144      * modular simulator approach.
145      */
146     void preStep(Step step, Time time, bool isNeighborSearchingStep);
147
148     /*! \brief A function called after every step
149      *
150      * This function collects calls which need to be made after every
151      * simulation step.
152      */
153     void postStep(Step step, Time time);
154
155     /*! \brief Add run functions to the task queue
156      *
157      * This function calls PME load balancing and domain decomposition first,
158      * and then queries the elements for all steps of the life time of the
159      * current neighbor list for their run functions. When the next step
160      * would be a neighbor-searching step, this function returns, such that
161      * the simulator can run the pre-computed task list before updating the
162      * neighbor list and re-filling the task list.
163      *
164      * TODO: In the current approach, unique pointers to tasks are created
165      *       repeatedly. While preliminary measurements do not seem to indicate
166      *       performance problems, a pool allocator would be an obvious
167      *       optimization possibility, as would reusing the task queue if
168      *       the task queue is periodic around the rebuild point (i.e. the
169      *       task queue is identical between rebuilds).
170      */
171     void populateTaskQueue();
172
173     //! The run queue
174     std::vector<SimulatorRunFunction> taskQueue_;
175     //! The task iterator
176     std::vector<SimulatorRunFunction>::const_iterator taskIterator_;
177
178     /* Note that the Simulator is owning the signallers and elements.
179      * The ownership list and the call list of the elements are kept
180      * separate, to allow to have elements more than once in the call
181      * lists - for example, using velocity verlet, the compute globals
182      * element needs to be scheduled more than once per step. For the
183      * signallers, no distinction between ownership and call list is
184      * made, all signallers are called exactly once per scheduling step.
185      *
186      * Objects being both a signaller and an element are not supported.
187      */
188     //! List of signalers (ownership and calling sequence)
189     std::vector<std::unique_ptr<ISignaller>> signallerList_;
190     //! List of schedulerElements (ownership)
191     std::vector<std::unique_ptr<ISimulatorElement>> elementsOwnershipList_;
192     //! List of schedulerElements (calling sequence)
193     std::vector<ISimulatorElement*> elementCallList_;
194
195     // Infrastructure elements
196     //! The domain decomposition element
197     std::unique_ptr<DomDecHelper> domDecHelper_;
198     //! The PME load balancing element
199     std::unique_ptr<PmeLoadBalanceHelper> pmeLoadBalanceHelper_;
200     //! The checkpoint helper
201     std::unique_ptr<CheckpointHelper> checkpointHelper_;
202     //! The stop handler
203     std::unique_ptr<StopHandler> stopHandler_;
204     //! The reset handler
205     std::unique_ptr<ResetHandler> resetHandler_;
206     //! Signal vector (used by stop / reset / checkpointing signaller)
207     std::unique_ptr<SimulationSignals> signals_;
208     //! The topology
209     std::unique_ptr<TopologyHolder> topologyHolder_;
210
211     // Data structures
212     //! The state propagator data
213     std::unique_ptr<StatePropagatorData> statePropagatorData_;
214     //! The energy data
215     std::unique_ptr<EnergyData> energyData_;
216     //! The free energy data
217     std::unique_ptr<FreeEnergyPerturbationData> freeEnergyPerturbationData_;
218
219     //! The current step
220     Step step_;
221     //! Whether the simulation run is finished
222     bool runFinished_;
223
224     /*! \internal
225      * \brief Signal helper
226      *
227      * The simulator needs to know about the last step and the
228      * neighbor searching step, which are determined in signallers.
229      * This object can be registered to the signals and accessed by
230      * the methods of the simulator.
231      */
232     class SignalHelper : public ILastStepSignallerClient, public INeighborSearchSignallerClient
233     {
234     public:
235         //! The last step
236         Step lastStep_ = std::numeric_limits<Step>::max();
237         //! The next NS step
238         Step nextNSStep_ = -1;
239         //! ILastStepSignallerClient implementation
240         std::optional<SignallerCallback> registerLastStepCallback() override;
241         //! INeighborSearchSignallerClient implementation
242         std::optional<SignallerCallback> registerNSCallback() override;
243     };
244     //! The signal helper object
245     std::unique_ptr<SignalHelper> signalHelper_;
246
247     // TODO: This is a hack for stop handler - needs to go once StopHandler
248     //       is adapted to the modular simulator
249     //! Whether this is a neighbor-searching step
250     bool stophandlerIsNSStep_ = false;
251     //! The current step
252     Step stophandlerCurrentStep_ = -1;
253
254     // ISimulator data - used for setup, teardown, pre step and post step
255     //! Name of the topology
256     std::string topologyName_;
257     //! Handles logging.
258     FILE* fplog;
259     //! Handles communication.
260     t_commrec* cr;
261     //! Handles logging.
262     const MDLogger& mdlog;
263     //! Contains command-line options to mdrun.
264     const MdrunOptions& mdrunOptions;
265     //! Contains user input mdp options.
266     const t_inputrec* inputrec;
267     //! Manages flop accounting.
268     t_nrnb* nrnb;
269     //! Manages wall cycle accounting.
270     gmx_wallcycle* wcycle;
271     //! Parameters for force calculations.
272     t_forcerec* fr;
273     //! Manages wall time accounting.
274     gmx_walltime_accounting* walltime_accounting;
275 };
276
277 /*! \internal
278  * \brief Helper container with data connected to global communication
279  *
280  * This includes data that needs to be shared between elements involved in
281  * global communication. This will become obsolete as soon as global
282  * communication is moved to a client system (#3421).
283  */
284 class GlobalCommunicationHelper
285 {
286 public:
287     //! Constructor
288     GlobalCommunicationHelper(int nstglobalcomm, SimulationSignals* simulationSignals);
289
290     //! Get the compute globals communication period
291     [[nodiscard]] int nstglobalcomm() const;
292     //! Get a pointer to the signals vector
293     [[nodiscard]] SimulationSignals* simulationSignals();
294
295 private:
296     //! Compute globals communication period
297     const int nstglobalcomm_;
298     //! Signal vector (used by stop / reset / checkpointing signaller)
299     SimulationSignals* simulationSignals_;
300 };
301
302 class ModularSimulatorAlgorithmBuilder;
303
304 /*! \internal
305  * \brief Helper for element addition
306  *
307  * Such an object will be given to each invocation of getElementPointer
308  *
309  * Note: It would be nicer to define this as a member type of
310  * ModularSimulatorAlgorithmBuilder, but this would break forward declaration.
311  * This object is therefore defined as friend class.
312  */
313 class ModularSimulatorAlgorithmBuilderHelper
314 {
315 public:
316     //! Constructor
317     ModularSimulatorAlgorithmBuilderHelper(ModularSimulatorAlgorithmBuilder* builder);
318     //! Store an element to the ModularSimulatorAlgorithmBuilder
319     ISimulatorElement* storeElement(std::unique_ptr<ISimulatorElement> element);
320     //! Check if an element is stored in the ModularSimulatorAlgorithmBuilder
321     bool elementIsStored(const ISimulatorElement* element) const;
322     //! Set arbitrary data in the ModularSimulatorAlgorithmBuilder. Helpful for stateful elements.
323     template<typename ValueType>
324     void storeValue(const std::string& key, const ValueType& value);
325     //! Get previously stored data. Returns std::nullopt if key is not found.
326     std::optional<std::any> getStoredValue(const std::string& key) const;
327     //! Register a thermostat that accepts propagator registrations
328     void registerThermostat(std::function<void(const PropagatorThermostatConnection&)> registrationFunction);
329     //! Register a barostat that accepts propagator registrations
330     void registerBarostat(std::function<void(const PropagatorBarostatConnection&)> registrationFunction);
331     //! Register a propagator to the thermostat used
332     void registerWithThermostat(PropagatorThermostatConnection connectionData);
333     //! Register a propagator to the barostat used
334     void registerWithBarostat(PropagatorBarostatConnection connectionData);
335
336 private:
337     //! Pointer to the associated ModularSimulatorAlgorithmBuilder
338     ModularSimulatorAlgorithmBuilder* builder_;
339     std::map<std::string, std::any>   values_;
340 };
341
342 /*!\internal
343  * \brief Builder for ModularSimulatorAlgorithm objects
344  *
345  * This builds a ModularSimulatorAlgorithm.
346  *
347  * Users can add elements and define their call order by calling the templated
348  * add<Element> function. Note that only elements that have a static
349  * getElementPointerImpl factory method can be built in that way.
350  *
351  * Note that each ModularSimulatorAlgorithmBuilder can only be used to build
352  * one ModularSimulatorAlgorithm object, i.e. build() can only be called once.
353  * During the call to build, all elements and other infrastructure objects will
354  * be moved to the built ModularSimulatorAlgorithm object, such that further use
355  * of the builder would not make sense.
356  * Any access to the build or add<> methods after the first call to
357  * build() will result in an exception being thrown.
358  */
359 class ModularSimulatorAlgorithmBuilder final
360 {
361 public:
362     //! Constructor
363     ModularSimulatorAlgorithmBuilder(compat::not_null<LegacySimulatorData*>    legacySimulatorData,
364                                      std::unique_ptr<ReadCheckpointDataHolder> checkpointDataHolder);
365     //! Build algorithm
366     ModularSimulatorAlgorithm build();
367
368     /*! \brief  Add element to the modular simulator algorithm builder
369      *
370      * This function has a general implementation, which will call the getElementPointer(...)
371      * factory function.
372      *
373      * \tparam Element  The element type
374      * \tparam Args     A variable number of argument types
375      * \param args      A variable number of arguments
376      */
377     template<typename Element, typename... Args>
378     void add(Args&&... args);
379
380     //! Allow access from helper
381     friend class ModularSimulatorAlgorithmBuilderHelper;
382
383 private:
384     //! The state of the builder
385     bool algorithmHasBeenBuilt_ = false;
386
387     // Data structures
388     //! The state propagator data
389     std::unique_ptr<StatePropagatorData> statePropagatorData_;
390     //! The energy data
391     std::unique_ptr<EnergyData> energyData_;
392     //! The free energy data
393     std::unique_ptr<FreeEnergyPerturbationData> freeEnergyPerturbationData_;
394
395     //! Pointer to the LegacySimulatorData object
396     compat::not_null<LegacySimulatorData*> legacySimulatorData_;
397
398     // Helper objects
399     //! Signal vector (used by stop / reset / checkpointing signaller)
400     std::unique_ptr<SimulationSignals> signals_;
401     //! Helper object passed to element factory functions
402     ModularSimulatorAlgorithmBuilderHelper elementAdditionHelper_;
403     //! Container for global computation data
404     GlobalCommunicationHelper globalCommunicationHelper_;
405
406     /*! \brief  Register an element to all applicable signallers and infrastructure elements
407      *
408      * \tparam Element  Type of the Element
409      * \param element   Pointer to the element
410      */
411     template<typename Element>
412     void registerWithInfrastructureAndSignallers(Element* element);
413
414     /*! \brief Take ownership of element
415      *
416      * This function returns a non-owning pointer to the new location of that
417      * element, allowing further usage (e.g. adding the element to the call list).
418      * Note that simply addin an element using this function will not call it
419      * during the simulation - it needs to be added to the call list separately.
420      * Note that generally, users will want to add elements to the call list, but
421      * it might not be practical to do this in the same order.
422      *
423      * \param element  A unique pointer to the element
424      * \return  A non-owning (raw) pointer to the element for further usage
425      */
426     ISimulatorElement* addElementToSimulatorAlgorithm(std::unique_ptr<ISimulatorElement> element);
427
428     /*! \brief Check if element is owned by *this
429      *
430      * \param element  Pointer to the element
431      * \return  Bool indicating whether element is owned by *this
432      */
433     [[nodiscard]] bool elementExists(const ISimulatorElement* element) const;
434
435     /*! \brief Add element to setupAndTeardownList_ if it's not already there
436      *
437      * \param element  Element pointer to be added
438      */
439     void addElementToSetupTeardownList(ISimulatorElement* element);
440
441     //! Vector to store elements, allowing the SimulatorAlgorithm to control their lifetime
442     std::vector<std::unique_ptr<ISimulatorElement>> elements_;
443     /*! \brief List defining in which order elements are called every step
444      *
445      * Elements may be referenced more than once if they should be called repeatedly
446      */
447     std::vector<ISimulatorElement*> callList_;
448     /*! \brief  List defining in which order elements are set up and torn down
449      *
450      * Elements should only appear once in this list
451      */
452     std::vector<ISimulatorElement*> setupAndTeardownList_;
453
454     //! Builder for the NeighborSearchSignaller
455     SignallerBuilder<NeighborSearchSignaller> neighborSearchSignallerBuilder_;
456     //! Builder for the LastStepSignaller
457     SignallerBuilder<LastStepSignaller> lastStepSignallerBuilder_;
458     //! Builder for the LoggingSignaller
459     SignallerBuilder<LoggingSignaller> loggingSignallerBuilder_;
460     //! Builder for the EnergySignaller
461     SignallerBuilder<EnergySignaller> energySignallerBuilder_;
462     //! Builder for the TrajectorySignaller
463     SignallerBuilder<TrajectorySignaller> trajectorySignallerBuilder_;
464     //! Builder for the TrajectoryElementBuilder
465     TrajectoryElementBuilder trajectoryElementBuilder_;
466     //! Builder for the TopologyHolder
467     TopologyHolder::Builder topologyHolderBuilder_;
468     //! Builder for the CheckpointHelper
469     CheckpointHelperBuilder checkpointHelperBuilder_;
470
471     /*! \brief List of clients for the CheckpointHelper
472      *
473      * \todo Replace this by proper builder (#3422)
474      */
475     std::vector<ICheckpointHelperClient*> checkpointClients_;
476
477     //! List of thermostat registration functions
478     std::vector<std::function<void(const PropagatorThermostatConnection&)>> thermostatRegistrationFunctions_;
479     //! List of barostat registration functions
480     std::vector<std::function<void(const PropagatorBarostatConnection&)>> barostatRegistrationFunctions_;
481     //! List of data to connect propagators to thermostats
482     std::vector<PropagatorThermostatConnection> propagatorThermostatConnections_;
483     //! List of data to connect propagators to barostats
484     std::vector<PropagatorBarostatConnection> propagatorBarostatConnections_;
485 };
486
487 /*! \internal
488  * \brief Factory function for elements that can be added via ModularSimulatorAlgorithmBuilder:
489  *        Get a pointer to an object of type \c Element to add to the call list
490  *
491  * This allows elements to be built via the templated ModularSimulatorAlgorithmBuilder::add<Element>
492  * method. Elements buildable throught this factor function are required to implement a static
493  * function with minimal signature
494  *
495  *     static ISimulatorElement* getElementPointerImpl(
496  *             LegacySimulatorData*                    legacySimulatorData,
497  *             ModularSimulatorAlgorithmBuilderHelper* builderHelper,
498  *             StatePropagatorData*                    statePropagatorData,
499  *             EnergyData*                             energyData,
500  *             FreeEnergyPerturbationData*             freeEnergyPerturbationData,
501  *             GlobalCommunicationHelper*              globalCommunicationHelper)
502  *
503  * This function may also accept additional parameters which are passed using the variadic
504  * template parameter pack forwarded in getElementPointer.
505  *
506  * This function returns a pointer to an object of the Element type. Note that the caller will
507  * check whether the returned object has previously been stored using the `storeElement`
508  * function, and throw an exception if the element is not found.
509  * The function can check whether a previously stored pointer is valid using
510  * the `checkElementExistence` function. Most implementing functions will simply want
511  * to create an object, store it using `storeElement`, and then use the return value of
512  * `storeElement` as a return value to the caller. However, this setup allows the function
513  * to store a created element (using a static pointer inside the function) and return it
514  * in case that the factory function is called repeatedly. This allows to create an element
515  * once, but have it called multiple times during the simulation run.
516  *
517  * \see ModularSimulatorAlgorithmBuilder::add
518  *      Function using this functionality
519  * \see ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>::getElementPointerImpl
520  *      Implementation using the single object / multiple call sites functionality
521  *
522  * \tparam Element The type of the element
523  * \tparam Args  Variable number of argument types allowing specific implementations to have
524  *               additional arguments
525  *
526  * \param legacySimulatorData  Pointer allowing access to simulator level data
527  * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
528  * \param statePropagatorData  Pointer to the \c StatePropagatorData object
529  * \param energyData  Pointer to the \c EnergyData object
530  * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
531  * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
532  * \param args  Variable number of additional parameters to be forwarded
533  *
534  * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
535  */
536 template<typename Element, typename... Args>
537 ISimulatorElement* getElementPointer(LegacySimulatorData*                    legacySimulatorData,
538                                      ModularSimulatorAlgorithmBuilderHelper* builderHelper,
539                                      StatePropagatorData*                    statePropagatorData,
540                                      EnergyData*                             energyData,
541                                      FreeEnergyPerturbationData* freeEnergyPerturbationData,
542                                      GlobalCommunicationHelper*  globalCommunicationHelper,
543                                      Args&&... args)
544 {
545     return Element::getElementPointerImpl(legacySimulatorData,
546                                           builderHelper,
547                                           statePropagatorData,
548                                           energyData,
549                                           freeEnergyPerturbationData,
550                                           globalCommunicationHelper,
551                                           std::forward<Args>(args)...);
552 }
553
554 template<typename Element, typename... Args>
555 void ModularSimulatorAlgorithmBuilder::add(Args&&... args)
556 {
557     if (algorithmHasBeenBuilt_)
558     {
559         throw SimulationAlgorithmSetupError(
560                 "Tried to add an element after ModularSimulationAlgorithm was built.");
561     }
562
563     // Get element from factory method
564     auto* element = static_cast<Element*>(getElementPointer<Element>(legacySimulatorData_,
565                                                                      &elementAdditionHelper_,
566                                                                      statePropagatorData_.get(),
567                                                                      energyData_.get(),
568                                                                      freeEnergyPerturbationData_.get(),
569                                                                      &globalCommunicationHelper_,
570                                                                      std::forward<Args>(args)...));
571
572     // Make sure returned element pointer is owned by *this
573     // Ensuring this makes sure we can control the life time
574     if (!elementExists(element))
575     {
576         throw ElementNotFoundError("Tried to append non-existing element to call list.");
577     }
578     // Add to call list
579     callList_.emplace_back(element);
580     // Add to setup / teardown list if element hasn't been added yet
581     addElementToSetupTeardownList(element);
582     // Register element to all applicable signallers
583     registerWithInfrastructureAndSignallers(element);
584 }
585
586 //! Returns a pointer casted to type Base if the Element is derived from Base
587 template<typename Base, typename Element>
588 static std::enable_if_t<std::is_base_of<Base, Element>::value, Base*> castOrNull(Element* element)
589 {
590     return static_cast<Base*>(element);
591 }
592
593 //! Returns a nullptr of type Base if Element is not derived from Base
594 template<typename Base, typename Element>
595 static std::enable_if_t<!std::is_base_of<Base, Element>::value, Base*> castOrNull(Element gmx_unused* element)
596 {
597     return nullptr;
598 }
599
600 template<typename Element>
601 void ModularSimulatorAlgorithmBuilder::registerWithInfrastructureAndSignallers(Element* element)
602 {
603     // Register element to all applicable signallers
604     neighborSearchSignallerBuilder_.registerSignallerClient(
605             castOrNull<INeighborSearchSignallerClient, Element>(element));
606     lastStepSignallerBuilder_.registerSignallerClient(castOrNull<ILastStepSignallerClient, Element>(element));
607     loggingSignallerBuilder_.registerSignallerClient(castOrNull<ILoggingSignallerClient, Element>(element));
608     energySignallerBuilder_.registerSignallerClient(castOrNull<IEnergySignallerClient, Element>(element));
609     trajectorySignallerBuilder_.registerSignallerClient(
610             castOrNull<ITrajectorySignallerClient, Element>(element));
611     // Register element to trajectory element (if applicable)
612     trajectoryElementBuilder_.registerWriterClient(castOrNull<ITrajectoryWriterClient, Element>(element));
613     // Register element to topology holder (if applicable)
614     topologyHolderBuilder_.registerClient(castOrNull<ITopologyHolderClient, Element>(element));
615     // Register element to checkpoint client (if applicable)
616     checkpointHelperBuilder_.registerClient(castOrNull<ICheckpointHelperClient, Element>(element));
617 }
618
619
620 template<typename ValueType>
621 void ModularSimulatorAlgorithmBuilderHelper::storeValue(const std::string& key, const ValueType& value)
622 {
623     values_[key] = std::any(value);
624 }
625
626
627 } // namespace gmx
628
629 #endif // GROMACS_MODULARSIMULATOR_SIMULATORALGORITHM_H