5c57dfabe137574859c1e32cdeaa1b7a23a9c97f
[alexxy/gromacs.git] / src / gromacs / awh / awh.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,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.
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
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.
41  *
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.
46  *
47  * \author Viveca Lindahl
48  * \author Berk Hess <hess@kth.se>
49  */
50
51 /*! \libinternal \file
52  *
53  * \brief
54  * Declares the Awh class.
55  *
56  * \author Viveca Lindahl
57  * \author Berk Hess <hess@kth.se>
58  * \inlibraryapi
59  * \ingroup module_awh
60  */
61
62 #ifndef GMX_AWH_H
63 #define GMX_AWH_H
64
65 #include <cstdio>
66
67 #include <memory>
68 #include <string>
69 #include <vector>
70
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/utility/basedefinitions.h"
73
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;
83
84 namespace gmx
85 {
86
87 struct AwhHistory;
88 struct AwhParams;
89 class Bias;
90 struct BiasCoupledToSystem;
91 class ForceWithVirial;
92
93 /*! \libinternal
94  * \brief Coupling of the accelerated weight histogram method (AWH) with the system.
95  *
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.
100  *
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.
105  *
106  * \todo Update parameter reading and checkpointing, when general C++ framework is ready.
107  */
108 class Awh
109 {
110 public:
111     /*! \brief Construct an AWH at the start of a simulation.
112      *
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.
117      *
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.
125      */
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);
133
134     ~Awh();
135
136     /*! \brief Peform an AWH update, to be called every MD step.
137      *
138      * An update has two tasks: apply the bias force and improve
139      * the bias and the free energy estimate that AWH keeps internally.
140      *
141      * For the first task, AWH retrieves the pull coordinate values from the pull struct.
142      * With these, the bias potential and forces are calculated.
143      * The bias force together with the atom forces and virial
144      * are passed on to pull which applies the bias force to the atoms.
145      * This is done at every step.
146      *
147      * Secondly, coordinate values are regularly sampled and kept by AWH.
148      * Convergence of the bias and free energy estimate is achieved by
149      * updating the AWH bias state after a certain number of samples has been collected.
150      *
151      * \note Requires that pull_potential from pull.h has been called first
152      * since AWH needs the current coordinate values (the pull code checks
153      * for this).
154      *
155      * \param[in]     mdatoms          Atom properties.
156      * \param[in]     ePBC             Type of periodic boundary conditions.
157      * \param[in]     box              Box vectors.
158      * \param[in,out] forceWithVirial  Force and virial buffers, should cover at least the local atoms.
159      * \param[in]     t                Time.
160      * \param[in]     step             The current MD step.
161      * \param[in,out] wallcycle        Wallcycle counter, can be nullptr.
162      * \param[in,out] fplog            General output file, normally md.log, can be nullptr.
163      * \returns the potential energy for the bias.
164      */
165     real applyBiasForcesAndUpdateBias(int                   ePBC,
166                                       const t_mdatoms&      mdatoms,
167                                       const matrix          box,
168                                       gmx::ForceWithVirial* forceWithVirial,
169                                       double                t,
170                                       int64_t               step,
171                                       gmx_wallcycle*        wallcycle,
172                                       FILE*                 fplog);
173
174     /*! \brief
175      * Update the AWH history in preparation for writing to checkpoint file.
176      *
177      * Should be called at least on the master rank at checkpoint steps.
178      *
179      * Should be called with a valid \p awhHistory (is checked).
180      *
181      * \param[in,out] awhHistory  AWH history to set.
182      */
183     void updateHistory(AwhHistory* awhHistory) const;
184
185     /*! \brief
186      * Allocate and initialize an AWH history with the given AWH state.
187      *
188      * This function should be called at the start of a new simulation
189      * at least on the master rank.
190      * Note that only constant data will be initialized here.
191      * History data is set by \ref Awh::updateHistory.
192      *
193      * \returns a shared pointer to the AWH history on the rank that does I/O, nullptr otherwise.
194      */
195     std::shared_ptr<AwhHistory> initHistoryFromState() const;
196
197     /*! \brief Restore the AWH state from the given history.
198      *
199      * Should be called on all ranks (for internal MPI broadcast).
200      * Should pass a point to an AwhHistory on the master rank that
201      * is compatible with the AWH setup in this simulation. Will throw
202      * an exception if it is not compatible.
203      *
204      * \param[in] awhHistory  AWH history to restore from.
205      */
206     void restoreStateFromHistory(const AwhHistory* awhHistory);
207
208     /*! \brief Fills the AWH data block of an energy frame with data at certain steps.
209      *
210      * \param[in]     step  The current MD step.
211      * \param[in,out] fr    Energy data frame.
212      */
213     void writeToEnergyFrame(int64_t step, t_enxframe* fr) const;
214
215     /*! \brief Returns string "AWH" for registering AWH as an external potential provider with the pull module.
216      */
217     static const char* externalPotentialString();
218
219     /*! \brief Register the AWH biased coordinates with pull.
220      *
221      * This function is public because it needs to be called by grompp
222      * (and is otherwise only called by Awh()).
223      * Pull requires all external potentials to register themselves
224      * before the end of pre-processing and before the first MD step.
225      * If this has not happened, pull with throw an error.
226      *
227      * \param[in]     awhParams  The AWH parameters.
228      * \param[in,out] pull_work  Pull struct which AWH will register the bias into.
229      */
230     static void registerAwhWithPull(const AwhParams& awhParams, pull_t* pull_work);
231
232 private:
233     /*! \brief Returns whether we need to write output at the current step.
234      *
235      * \param[in]     step             The current MD step.
236      */
237     bool isOutputStep(int64_t step) const;
238
239 private:
240     std::vector<BiasCoupledToSystem> biasCoupledToSystem_; /**< AWH biases and definitions of their coupling to the system. */
241     const int64_t    seed_;   /**< Random seed for MC jumping with umbrella type bias potential. */
242     const int        nstout_; /**< Interval in steps for writing to energy file. */
243     const t_commrec* commRecord_;          /**< Pointer to the communication record. */
244     const gmx_multisim_t* multiSimRecord_; /**< Handler for multi-simulations. */
245     pull_t*               pull_;           /**< Pointer to the pull working data. */
246     double                potentialOffset_; /**< The offset of the bias potential which changes due to bias updates. */
247 };
248
249 /*! \brief Makes an Awh and prepares to use it if the user input
250  * requests that
251  *
252  * Restores state from history in checkpoint if needed.
253  *
254  * \param[in,out] fplog                   General output file, normally md.log, can be nullptr.
255  * \param[in]     inputRecord             General input parameters (as set up by grompp).
256  * \param[in]     stateGlobal             A pointer to the global state structure.
257  * \param[in]     commRecord              Struct for communication, can be nullptr.
258  * \param[in]     multiSimRecord          Multi-sim handler
259  * \param[in]     startingFromCheckpoint  Whether the simulation is starting from a checkpoint
260  * \param[in]     usingShellParticles     Whether the user requested shell particles (which is unsupported)
261  * \param[in]     biasInitFilename        Name of file to read PMF and target from.
262  * \param[in,out] pull_work               Pointer to a pull struct which AWH will couple to, has to be initialized,
263  *                                        is assumed not to change during the lifetime of the Awh object.
264  * \returns       An initialized Awh module, or nullptr if none was requested.
265  * \throws        InvalidInputError       If another active module is not supported.
266  */
267 std::unique_ptr<Awh> prepareAwhModule(FILE*                 fplog,
268                                       const t_inputrec&     inputRecord,
269                                       t_state*              stateGlobal,
270                                       const t_commrec*      commRecord,
271                                       const gmx_multisim_t* multiSimRecord,
272                                       bool                  startingFromCheckpoint,
273                                       bool                  usingShellParticles,
274                                       const std::string&    biasInitFilename,
275                                       pull_t*               pull_work);
276
277 } // namespace gmx
278
279 #endif /* GMX_AWH_H */