Update Awh initialization and lifetime management
[gromacs.git] / src / gromacs / awh / awh.h
blob858da6b7566d989db61bff50c39f931d3339f802
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018, 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.
36 /*! \libinternal
37 * \defgroup module_awh Accelerated weight histogram (AWH) method
38 * \ingroup group_mdrun
39 * \brief
40 * Implements the "accelerated weight histogram" sampling method.
42 * This class provides the interface between the AWH module and
43 * other modules using it. Currently AWH can only act on COM pull
44 * reaction coordinates, but this can easily be extended to other
45 * types of reaction coordinates.
47 * \author Viveca Lindahl
48 * \author Berk Hess <hess@kth.se>
51 /*! \libinternal \file
53 * \brief
54 * Declares the Awh class.
56 * \author Viveca Lindahl
57 * \author Berk Hess <hess@kth.se>
58 * \inlibraryapi
59 * \ingroup module_awh
62 #ifndef GMX_AWH_H
63 #define GMX_AWH_H
65 #include <cstdio>
67 #include <memory>
68 #include <string>
69 #include <vector>
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/utility/basedefinitions.h"
74 struct gmx_multisim_t;
75 struct gmx_wallcycle;
76 struct pull_work_t;
77 struct pull_t;
78 class t_state;
79 struct t_commrec;
80 struct t_enxframe;
81 struct t_inputrec;
82 struct t_mdatoms;
84 namespace gmx
87 struct AwhHistory;
88 struct AwhParams;
89 class Bias;
90 struct BiasCoupledToSystem;
91 class ForceWithVirial;
93 /*! \libinternal
94 * \brief Coupling of the accelerated weight histogram method (AWH) with the system.
96 * AWH calculates the free energy along order parameters of the system.
97 * Free energy barriers are overcome by adaptively tuning a bias potential along
98 * the order parameter such that the biased distribution along the parameter
99 * converges toward a chosen target distribution.
101 * The Awh class takes care of the coupling between the system and the AWH
102 * bias(es). The Awh class contains one or more BiasCoupledToSystem objects.
103 * The BiasCoupledToSystem class takes care of the reaction coordinate input
104 * and force output for the single Bias object it containts.
106 * \todo Update parameter reading and checkpointing, when general C++ framework is ready.
108 class Awh
110 public:
111 /*! \brief Construct an AWH at the start of a simulation.
113 * AWH will here also register itself with the pull struct as the
114 * potential provider for the pull coordinates given as AWH coordinates
115 * in the user input. This allows AWH to later apply the bias force to
116 * these coordinate in \ref Awh::applyBiasForcesAndUpdateBias.
118 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
119 * \param[in] inputRecord General input parameters (as set up by grompp).
120 * \param[in] commRecord Struct for communication, can be nullptr.
121 * \param[in] multiSimRecord Multi-sim handler
122 * \param[in] awhParams AWH input parameters, consistent with the relevant parts of \p inputRecord (as set up by grompp).
123 * \param[in] biasInitFilename Name of file to read PMF and target from.
124 * \param[in,out] pull_work Pointer to a pull struct which AWH will couple to, has to be initialized, is assumed not to change during the lifetime of the Awh object.
126 Awh(FILE *fplog,
127 const t_inputrec &inputRecord,
128 const t_commrec *commRecord,
129 const gmx_multisim_t *multiSimRecord,
130 const AwhParams &awhParams,
131 const std::string &biasInitFilename,
132 pull_t *pull_work);
134 /*! \brief Destructor. */
135 ~Awh();
137 /*! \brief Peform an AWH update, to be called every MD step.
139 * An update has two tasks: apply the bias force and improve
140 * the bias and the free energy estimate that AWH keeps internally.
142 * For the first task, AWH retrieves the pull coordinate values from the pull struct.
143 * With these, the bias potential and forces are calculated.
144 * The bias force together with the atom forces and virial
145 * are passed on to pull which applies the bias force to the atoms.
146 * This is done at every step.
148 * Secondly, coordinate values are regularly sampled and kept by AWH.
149 * Convergence of the bias and free energy estimate is achieved by
150 * updating the AWH bias state after a certain number of samples has been collected.
152 * \note Requires that pull_potential from pull.h has been called first
153 * since AWH needs the current coordinate values (the pull code checks
154 * for this).
156 * \param[in] mdatoms Atom properties.
157 * \param[in] ePBC Type of periodic boundary conditions.
158 * \param[in] box Box vectors.
159 * \param[in,out] forceWithVirial Force and virial buffers, should cover at least the local atoms.
160 * \param[in] t Time.
161 * \param[in] step The current MD step.
162 * \param[in,out] wallcycle Wallcycle counter, can be nullptr.
163 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
164 * \returns the potential energy for the bias.
166 real applyBiasForcesAndUpdateBias(int ePBC,
167 const t_mdatoms &mdatoms,
168 const matrix box,
169 gmx::ForceWithVirial *forceWithVirial,
170 double t,
171 gmx_int64_t step,
172 gmx_wallcycle *wallcycle,
173 FILE *fplog);
175 /*! \brief
176 * Update the AWH history in preparation for writing to checkpoint file.
178 * Should be called at least on the master rank at checkpoint steps.
180 * Should be called with a valid \p awhHistory (is checked).
182 * \param[in,out] awhHistory AWH history to set.
184 void updateHistory(AwhHistory *awhHistory) const;
186 /*! \brief
187 * Allocate and initialize an AWH history with the given AWH state.
189 * This function should be called at the start of a new simulation
190 * at least on the master rank.
191 * Note that only constant data will be initialized here.
192 * History data is set by \ref Awh::updateHistory.
194 * \returns a shared pointer to the AWH history on the rank that does I/O, nullptr otherwise.
196 std::shared_ptr<AwhHistory> initHistoryFromState() const;
198 /*! \brief Restore the AWH state from the given history.
200 * Should be called on all ranks (for internal MPI broadcast).
201 * Should pass a point to an AwhHistory on the master rank that
202 * is compatible with the AWH setup in this simulation. Will throw
203 * an exception if it is not compatible.
205 * \param[in] awhHistory AWH history to restore from.
207 void restoreStateFromHistory(const AwhHistory *awhHistory);
209 /*! \brief Fills the AWH data block of an energy frame with data at certain steps.
211 * \param[in] step The current MD step.
212 * \param[in,out] fr Energy data frame.
214 void writeToEnergyFrame(gmx_int64_t step,
215 t_enxframe *fr) const;
217 /*! \brief Returns string "AWH" for registering AWH as an external potential provider with the pull module.
219 static const char *externalPotentialString();
221 /*! \brief Register the AWH biased coordinates with pull.
223 * This function is public because it needs to be called by grompp
224 * (and is otherwise only called by Awh()).
225 * Pull requires all external potentials to register themselves
226 * before the end of pre-processing and before the first MD step.
227 * If this has not happened, pull with throw an error.
229 * \param[in] awhParams The AWH parameters.
230 * \param[in,out] pull_work Pull struct which AWH will register the bias into.
232 static void registerAwhWithPull(const AwhParams &awhParams,
233 pull_t *pull_work);
235 private:
236 /*! \brief Returns whether we need to write output at the current step.
238 * \param[in] step The current MD step.
240 bool isOutputStep(gmx_int64_t step) const;
242 private:
243 std::vector<BiasCoupledToSystem> biasCoupledToSystem_; /**< AWH biases and definitions of their coupling to the system. */
244 const gmx_int64_t seed_; /**< Random seed for MC jumping with umbrella type bias potential. */
245 const int nstout_; /**< Interval in steps for writing to energy file. */
246 const t_commrec *commRecord_; /**< Pointer to the communication record. */
247 const gmx_multisim_t *multiSimRecord_; /**< Handler for multi-simulations. */
248 pull_t *pull_; /**< Pointer to the pull working data. */
249 double potentialOffset_; /**< The offset of the bias potential which changes due to bias updates. */
252 /*! \brief Makes an Awh and prepares to use it if the user input
253 * requests that
255 * Restores state from history in checkpoint if needed.
257 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
258 * \param[in] inputRecord General input parameters (as set up by grompp).
259 * \param[in] stateGlobal A pointer to the global state structure.
260 * \param[in] commRecord Struct for communication, can be nullptr.
261 * \param[in] multiSimRecord Multi-sim handler
262 * \param[in] startingFromCheckpoint Whether the simulation is starting from a checkpoint
263 * \param[in] usingShellParticles Whether the user requested shell particles (which is unsupported)
264 * \param[in] biasInitFilename Name of file to read PMF and target from.
265 * \param[in,out] pull_work Pointer to a pull struct which AWH will couple to, has to be initialized,
266 * is assumed not to change during the lifetime of the Awh object.
267 * \returns An initialized Awh module, or nullptr if none was requested.
268 * \throws InvalidInputError If another active module is not supported.
270 std::unique_ptr<Awh>
271 prepareAwhModule(FILE *fplog,
272 const t_inputrec &inputRecord,
273 t_state *stateGlobal,
274 const t_commrec *commRecord,
275 const gmx_multisim_t *multiSimRecord,
276 const bool startingFromCheckpoint,
277 const bool usingShellParticles,
278 const std::string &biasInitFilename,
279 pull_t *pull_work);
281 } // namespace gmx
283 #endif /* GMX_AWH_H */