Adaptive force scaling for densityfitting
[gromacs.git] / src / gromacs / utility / mdmodulenotification.h
blob0d8ecb39532cb4159afee4604208f32588695ba5
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 \file
36 * \brief
37 * Declares gmx::MdModuleNotification.
39 * \author Christian Blau <blau@kth.se>
40 * \inlibraryapi
41 * \ingroup module_utility
44 #ifndef GMX_MDRUNUTILITY_MDMODULENOTIFICATION_H
45 #define GMX_MDRUNUTILITY_MDMODULENOTIFICATION_H
47 #include <functional>
48 #include <string>
49 #include <vector>
51 struct t_commrec;
53 namespace gmx
56 /*! \libinternal \brief
57 * Subscribe and trigger notification functions.
59 * Extends MdModuleNotificationBase with new notification function and routine
60 * to subscribe new listeners.
62 * To create a class of this type that provides callbacks, e.g., for events
63 * EventA, and EventB use registerMdModuleNotification<EventA, EventB>::type.
65 * \tparam CallParameter of the function to be notified
66 * \tparam MdModuleNotificationBase class to be extended with a notification
67 * with CallParameter
69 \msc
70 wordwraparcs=true,
71 hscale="2";
73 runner [label="runner:\nMdrunner"],
74 CallParameter [label = "eventA:\nCallParameter"],
75 MOD [label = "mdModules_:\nMdModules"],
76 ModuleA [label="moduleA"],
77 ModuleB [label="moduleB"],
78 MdModuleNotification [label="notifier_:\nMdModuleNotification"];
80 MOD box MdModuleNotification [label = "mdModules_ owns notifier_ and moduleA/B"];
81 MOD =>> ModuleA [label="instantiates(notifier_)"];
82 ModuleA =>> MdModuleNotification [label="subscribe(otherfunc)"];
83 ModuleA =>> MOD;
84 MOD =>> ModuleB [label="instantiates(notifier_)"];
85 ModuleB =>> MdModuleNotification [label="subscribe(func)"];
86 ModuleB =>> MOD;
87 runner =>> CallParameter [label="instantiate"];
88 CallParameter =>> runner ;
89 runner =>> MOD [label="notify(eventA)"];
90 MOD =>> MdModuleNotification [label="notify(eventA)"];
91 MdModuleNotification =>> ModuleA [label="notify(eventA)"];
92 ModuleA -> ModuleA [label="func(eventA)"];
93 MdModuleNotification =>> ModuleB [label="notify(eventA)"];
94 ModuleB -> ModuleB [label="otherfunc(eventA)"];
96 \endmsc
98 * \note All added subscribers are required to out-live the MdModuleNotification
101 template <class CallParameter, class MdModuleNotificationBase>
102 class MdModuleNotification : public MdModuleNotificationBase
104 public:
105 //! Make base class notification trigger available to this class
106 using MdModuleNotificationBase::notify;
107 //! Make base class subscription available to this class
108 using MdModuleNotificationBase::subscribe;
110 /*! \brief Trigger the subscribed notifications.
111 * \param[in] callParameter of the function to be called back
113 void notify(CallParameter callParameter) const
115 for (auto &callBack : callBackFunctions_)
117 callBack(callParameter);
121 /*! \brief
122 * Add callback function to be called when notification is triggered.
124 * Notifications are distinguished by their call signature.
126 * \param[in] callBackFunction to be called from this class
128 void subscribe(std::function<void(CallParameter)> callBackFunction)
130 callBackFunctions_.emplace_back(callBackFunction);
133 private:
134 std::vector < std::function < void(CallParameter)>> callBackFunctions_;
137 /*! \internal
138 * \brief Aide to avoid nested MdModuleNotification definition.
140 * Instead of
141 * MdModuleNotification<CallParameterA, MdModuleNotification<CallParameterB, etc ... >>
142 * this allows to write
143 * registerMdModuleNotification<CallParameterA, CallParameterB, ...>::type
145 * \tparam CallParameter all the event types to be registered
147 template <class ... CallParameter> struct registerMdModuleNotification;
149 /*! \internal \brief Template specialization to end parameter unpacking recursion.
151 template <>
152 struct registerMdModuleNotification<>
154 /*! \internal
155 * \brief Do nothing but be base class of MdModuleNotification.
157 * Required so that using MdModuleNotificationBase::notify and
158 * MdModuleNotificationBase::subscribe are valid in derived class.
160 class NoCallParameter
162 public:
163 //! Do nothing but provide MdModuleNotification::notify to derived class
164 void notify() {}
165 //! Do nothing but provide MdModuleNotification::subscribe to derived class
166 void subscribe() {}
168 /*! \brief Defines a type if no notifications are managed.
170 * This ensures that code works with MdModuleCallParameterManagement that
171 * does not manage any notifications.
173 using type = NoCallParameter;
176 /*! \libinternal
177 * \brief Template specialization to assemble MdModuleNotification.
179 * Assembly of MdModuleNotification is performed by recursively taking off the
180 * front of the CallParameter parameter pack and constructing the nested type
181 * definition of MdModuleNotification base classes.
183 * \tparam CurrentCallParameter front of the template parameter pack
184 * \tparam CallParameter rest of the event types
186 template <class CurrentCallParameter, class ... CallParameter>
187 struct registerMdModuleNotification<CurrentCallParameter, CallParameter...>
189 // private:
190 //! The next type with rest of the arguments with the front parameter removed.
191 using next_type = typename registerMdModuleNotification<CallParameter...>::type;
192 //! The type of the MdModuleNotification
193 using type = MdModuleNotification<CurrentCallParameter, next_type>;
196 class KeyValueTreeObject;
197 class KeyValueTreeObjectBuilder;
198 class LocalAtomSetManager;
199 class IndexGroupsAndNames;
200 struct MdModulesCheckpointReadingDataOnMaster;
201 struct MdModulesCheckpointReadingBroadcast;
202 struct MdModulesWriteCheckpointData;
203 struct PeriodicBoundaryConditionType
205 int pbcType;
208 struct MdModulesEnergyOutputToDensityFittingRequestChecker
210 bool energyOutputToDensityFitting_ = false;
213 /*! \libinternal
214 * \brief Collect errors for the energy calculation frequency.
216 * Collect errors regarding energy calculation frequencies as strings that then
217 * may be used to issue errors.
219 * \note The mdp option "nstcalcenergy" is altered after reading the .mdp input
220 * and only used in certain integrators, thus this class is to be used
221 * only after all these operations are done.
223 class EnergyCalculationFrequencyErrors
225 public:
226 //! Construct by setting the energy calculation frequency
227 EnergyCalculationFrequencyErrors(int64_t energyCalculationIntervalInSteps) :
228 energyCalculationIntervalInSteps_(energyCalculationIntervalInSteps){}
229 //! Return the number of steps of an energy calculation interval
230 std::int64_t energyCalculationIntervalInSteps() const
232 return energyCalculationIntervalInSteps_;
234 //! Collect error messages
235 void addError(const std::string &errorMessage)
237 errorMessages_.push_back(errorMessage);
239 //! Return error messages
240 const std::vector<std::string> &errorMessages() const
242 return errorMessages_;
244 private:
245 //! The frequency of energy calculations
246 const std::int64_t energyCalculationIntervalInSteps_;
247 //! The error messages
248 std::vector<std::string> errorMessages_;
251 struct SimulationTimeStep
253 //! Time step (ps)
254 const double delta_t;
257 struct MdModulesNotifier
259 //! Register callback function types for MdModule
260 registerMdModuleNotification<
261 const t_commrec &,
262 EnergyCalculationFrequencyErrors *,
263 IndexGroupsAndNames,
264 KeyValueTreeObjectBuilder,
265 const KeyValueTreeObject &,
266 LocalAtomSetManager *,
267 MdModulesEnergyOutputToDensityFittingRequestChecker *,
268 MdModulesCheckpointReadingDataOnMaster,
269 MdModulesCheckpointReadingBroadcast,
270 MdModulesWriteCheckpointData,
271 PeriodicBoundaryConditionType,
272 const SimulationTimeStep &>::type notifier_;
275 } // namespace gmx
277 #endif