2 * This file is part of the GROMACS molecular simulation package.
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.
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 \file
37 * Declares gmx::MDModulesNotifiers.
39 * \author Christian Blau <blau@kth.se>
41 * \ingroup module_utility
44 #ifndef GMX_UTILITY_MDMODULESNOTIFIERS_H
45 #define GMX_UTILITY_MDMODULESNOTIFIERS_H
50 #include "gromacs/math/arrayrefwithpadding.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/utility/mdmodulesnotifier.h"
56 enum class PbcType : int;
61 class KeyValueTreeObject;
62 class KeyValueTreeObjectBuilder;
63 class LocalAtomSetManager;
65 class IndexGroupsAndNames;
66 class SeparatePmeRanksPermitted;
67 struct MDModulesCheckpointReadingDataOnMaster;
68 struct MDModulesCheckpointReadingBroadcast;
69 struct MDModulesWriteCheckpointData;
71 /*! \libinternal \brief Check if module outputs energy to a specific field.
73 * Ensures that energy is output for this module.
75 struct MDModulesEnergyOutputToDensityFittingRequestChecker
77 //! Trigger output to density fitting energy field
78 bool energyOutputToDensityFitting_ = false;
82 * \brief Collect errors for the energy calculation frequency.
84 * Collect errors regarding energy calculation frequencies as strings that then
85 * may be used to issue errors.
87 * \note The mdp option "nstcalcenergy" is altered after reading the .mdp input
88 * and only used in certain integrators, thus this class is to be used
89 * only after all these operations are done.
91 class EnergyCalculationFrequencyErrors
94 //! Construct by setting the energy calculation frequency
95 EnergyCalculationFrequencyErrors(int64_t energyCalculationIntervalInSteps) :
96 energyCalculationIntervalInSteps_(energyCalculationIntervalInSteps)
99 //! Return the number of steps of an energy calculation interval
100 std::int64_t energyCalculationIntervalInSteps() const
102 return energyCalculationIntervalInSteps_;
104 //! Collect error messages
105 void addError(const std::string& errorMessage) { errorMessages_.push_back(errorMessage); }
106 //! Return error messages
107 const std::vector<std::string>& errorMessages() const { return errorMessages_; }
110 //! The frequency of energy calculations
111 const std::int64_t energyCalculationIntervalInSteps_;
112 //! The error messages
113 std::vector<std::string> errorMessages_;
116 /*! \libinternal \brief Provides the simulation time step in ps.
118 struct SimulationTimeStep
121 const double delta_t;
124 /*! \libinternal \brief Provides coordinates and simulation box.
126 struct CoordinatesAndBoxPreprocessed
128 ArrayRefWithPadding<RVec> coordinates_;
133 /*! \libinternal \brief Mdrun input filename.
135 struct MdRunInputFilename
137 //! The name of the run input file (.tpr) as output by grompp
138 std::string mdRunFilename_;
142 * \brief Group of notifers to organize that MDModules
143 * can receive callbacks they subscribe to.
145 * MDModules use members of this struct to subscribe to notifications
146 * of particular events. When the event occurs, the callback provided
147 * by a particular MDModule will be passed a parameter of the
148 * particular type they are interested in.
150 * Typically, during the setup phase, modules subscribe to notifiers
151 * that interest them by passing callbacks that expect a single parameter
152 * that describes the event. These are stored for later use. See the
153 * sequence diagram that follows:
158 modules [label = "mdModules:\nMDModules"],
159 notifiers [label="notifiers\nMDModulesNotifiers"],
160 notifier [label="exampleNotifier:\nBuildMDModulesNotifier\n<EventX, EventY>::type"],
161 moduleA [label="moduleA"],
162 moduleB [label="moduleB"],
163 moduleC [label="moduleC"];
165 modules box moduleC [label = "mdModules creates and owns moduleA, moduleB, and moduleC"];
166 modules =>> notifiers [label="creates"];
167 notifiers =>> notifier [label="creates"];
168 notifier =>> notifiers [label="returns"];
169 notifiers =>> modules [label="returns"];
171 modules =>> moduleA [label="provides notifiers"];
172 moduleA =>> moduleA [label="unpacks\nnotifiers.exampleNotifier"];
173 moduleA =>> notifier [label="subscribes with\ncallback(EventX&)"];
174 notifier =>> notifier [label="records subscription\nto EventX"];
175 moduleA =>> notifier [label="subscribes with\ncallback(EventY&)"];
176 notifier =>> notifier [label="records subscription\nto EventY"];
177 moduleA =>> modules [label="returns"];
179 modules =>> moduleB [label="provides notifiers"];
180 moduleB =>> moduleB [label="unpacks\nnotifiers.exampleNotifier"];
181 moduleA =>> notifier [label="subscribes with\ncallback(EventY&)"];
182 notifier =>> notifier [label="records subscription\nto EventY"];
183 moduleB =>> modules [label="returns"];
185 modules =>> moduleC [label="provides notifiers"];
186 moduleC =>> moduleC [label="unpacks and keeps\nnotifiers.exampleNotifier"];
187 moduleC =>> modules [label="returns"];
191 * When the event occurs later on, the stored callbacks are used to
192 * allow the modules to react. See the following sequence diagram,
193 * which assumes that exampleNotifier was configured as in the
194 * previous sequence diagram.
200 moduleC [label="moduleC"],
201 notifier [label="exampleNotifier:\nBuildMDModulesNotifier\n<EventX, EventY>::type"],
202 moduleA [label="moduleA"],
203 moduleB [label="moduleB"];
205 moduleC box moduleB [label = "Later, when ModuleC is doing work"];
206 moduleC =>> moduleC [label="generates EventX"];
207 moduleC =>> moduleC [label="generates EventY"];
208 moduleC =>> notifier [label="calls notify(eventX)"];
209 notifier =>> moduleA [label="calls callback(eventX)"];
210 moduleA =>> moduleA [label="reacts to eventX"];
211 moduleA =>> notifier [label="returns"];
213 notifier =>> moduleC [label="returns"];
214 moduleC =>> notifier [label="calls notify(eventY)"];
215 notifier =>> moduleA [label="calls callback(eventY)"];
216 moduleA =>> moduleA [label="reacts to eventY"];
217 moduleA =>> notifier [label="returns"];
218 notifier =>> moduleB [label="calls callback(eventY)"];
219 moduleB =>> moduleB [label="reacts to eventY"];
220 moduleB =>> notifier [label="returns"];
221 notifier =>> moduleC [label="returns"];
224 * The template arguments to the members of this struct are the
225 * parameters passed to the callback functions, one type per
226 * callback. Arguments passed as pointers are always meant to be
227 * modified, but never meant to be stored (in line with the policy
231 struct MDModulesNotifiers
233 /*! \brief Pre-processing callback functions.
234 * const CoordinatesAndBoxPreprocessed Allows modules to access coordinates,
235 * box and pbc during grompp
236 * const MDLogger& Allows MdModule to use standard logging class for messages output
237 * EnergyCalculationFrequencyErrors* allows modules to check if they match
238 * their required calculation frequency
239 * and add their error message if needed
240 * to the collected error messages
241 * gmx_mtop_t * Allows modules to modify the topology during pre-processing
242 * IndexGroupsAndNames provides modules with atom indices and their names
243 * KeyValueTreeObjectBuilder enables writing of module internal data to
246 BuildMDModulesNotifier<const CoordinatesAndBoxPreprocessed,
248 EnergyCalculationFrequencyErrors*,
250 const IndexGroupsAndNames&,
251 KeyValueTreeObjectBuilder>::type preProcessingNotifier_;
253 /*! \brief Handles subscribing and calling checkpointing callback functions.
255 * MDModulesCheckpointReadingDataOnMaster provides modules with their
256 * checkpointed data on the master
257 * node and checkpoint file version
258 * MDModulesCheckpointReadingBroadcast provides modules with a communicator
259 * and the checkpoint file version to
260 * distribute their data
261 * MDModulesWriteCheckpointData provides the modules with a key-value-tree
262 * builder to store their checkpoint data and
263 * the checkpoint file version
265 BuildMDModulesNotifier<MDModulesCheckpointReadingDataOnMaster, MDModulesCheckpointReadingBroadcast, MDModulesWriteCheckpointData>::type
266 checkpointingNotifier_;
268 /*! \brief Handles subscribing and calling callbacks during simulation setup.
270 * const KeyValueTreeObject& provides modules with the internal data they
271 * wrote to .tpr files
272 * LocalAtomSetManager* enables modules to add atom indices to local atom sets
274 * const MDLogger& Allows MdModule to use standard logging class for messages output
275 * const gmx_mtop_t& provides the topology of the system to the modules
276 * MDModulesEnergyOutputToDensityFittingRequestChecker* enables modules to
277 * report if they want to write their energy output
278 * to the density fitting field in the energy files
279 * SeparatePmeRanksPermitted* enables modules to report if they want
280 * to disable dedicated PME ranks
281 * const PbcType& provides modules with the periodic boundary condition type
282 * that is used during the simulation
283 * const SimulationTimeStep& provides modules with the simulation time-step
284 * that allows them to interconvert between step
286 * const t_commrec& provides a communicator to the modules during simulation
288 * MdRunInputFilename Allows modules to know .tpr filename during mdrun
290 BuildMDModulesNotifier<const KeyValueTreeObject&,
291 LocalAtomSetManager*,
294 MDModulesEnergyOutputToDensityFittingRequestChecker*,
295 SeparatePmeRanksPermitted*,
297 const SimulationTimeStep&,
299 MdRunInputFilename>::type simulationSetupNotifier_;