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.
37 * \defgroup module_awh Accelerated weight histogram (AWH) method
38 * \ingroup group_mdrun
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
54 * Declares the Awh class.
56 * \author Viveca Lindahl
57 * \author Berk Hess <hess@kth.se>
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/utility/basedefinitions.h"
74 struct gmx_multisim_t
;
90 struct BiasCoupledToSystem
;
91 class ForceWithVirial
;
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.
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.
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
,
136 /*! \brief Peform an AWH update, to be called every MD step.
138 * An update has two tasks: apply the bias force and improve
139 * the bias and the free energy estimate that AWH keeps internally.
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.
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.
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
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.
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.
165 real
applyBiasForcesAndUpdateBias(int ePBC
,
166 const t_mdatoms
&mdatoms
,
168 gmx::ForceWithVirial
*forceWithVirial
,
171 gmx_wallcycle
*wallcycle
,
175 * Update the AWH history in preparation for writing to checkpoint file.
177 * Should be called at least on the master rank at checkpoint steps.
179 * Should be called with a valid \p awhHistory (is checked).
181 * \param[in,out] awhHistory AWH history to set.
183 void updateHistory(AwhHistory
*awhHistory
) const;
186 * Allocate and initialize an AWH history with the given AWH state.
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.
193 * \returns a shared pointer to the AWH history on the rank that does I/O, nullptr otherwise.
195 std::shared_ptr
<AwhHistory
> initHistoryFromState() const;
197 /*! \brief Restore the AWH state from the given history.
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.
204 * \param[in] awhHistory AWH history to restore from.
206 void restoreStateFromHistory(const AwhHistory
*awhHistory
);
208 /*! \brief Fills the AWH data block of an energy frame with data at certain steps.
210 * \param[in] step The current MD step.
211 * \param[in,out] fr Energy data frame.
213 void writeToEnergyFrame(int64_t step
,
214 t_enxframe
*fr
) const;
216 /*! \brief Returns string "AWH" for registering AWH as an external potential provider with the pull module.
218 static const char *externalPotentialString();
220 /*! \brief Register the AWH biased coordinates with pull.
222 * This function is public because it needs to be called by grompp
223 * (and is otherwise only called by Awh()).
224 * Pull requires all external potentials to register themselves
225 * before the end of pre-processing and before the first MD step.
226 * If this has not happened, pull with throw an error.
228 * \param[in] awhParams The AWH parameters.
229 * \param[in,out] pull_work Pull struct which AWH will register the bias into.
231 static void registerAwhWithPull(const AwhParams
&awhParams
,
235 /*! \brief Returns whether we need to write output at the current step.
237 * \param[in] step The current MD step.
239 bool isOutputStep(int64_t step
) const;
242 std::vector
<BiasCoupledToSystem
> biasCoupledToSystem_
; /**< AWH biases and definitions of their coupling to the system. */
243 const int64_t seed_
; /**< Random seed for MC jumping with umbrella type bias potential. */
244 const int nstout_
; /**< Interval in steps for writing to energy file. */
245 const t_commrec
*commRecord_
; /**< Pointer to the communication record. */
246 const gmx_multisim_t
*multiSimRecord_
; /**< Handler for multi-simulations. */
247 pull_t
*pull_
; /**< Pointer to the pull working data. */
248 double potentialOffset_
; /**< The offset of the bias potential which changes due to bias updates. */
251 /*! \brief Makes an Awh and prepares to use it if the user input
254 * Restores state from history in checkpoint if needed.
256 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
257 * \param[in] inputRecord General input parameters (as set up by grompp).
258 * \param[in] stateGlobal A pointer to the global state structure.
259 * \param[in] commRecord Struct for communication, can be nullptr.
260 * \param[in] multiSimRecord Multi-sim handler
261 * \param[in] startingFromCheckpoint Whether the simulation is starting from a checkpoint
262 * \param[in] usingShellParticles Whether the user requested shell particles (which is unsupported)
263 * \param[in] biasInitFilename Name of file to read PMF and target from.
264 * \param[in,out] pull_work Pointer to a pull struct which AWH will couple to, has to be initialized,
265 * is assumed not to change during the lifetime of the Awh object.
266 * \returns An initialized Awh module, or nullptr if none was requested.
267 * \throws InvalidInputError If another active module is not supported.
270 prepareAwhModule(FILE *fplog
,
271 const t_inputrec
&inputRecord
,
272 t_state
*stateGlobal
,
273 const t_commrec
*commRecord
,
274 const gmx_multisim_t
*multiSimRecord
,
275 bool startingFromCheckpoint
,
276 bool usingShellParticles
,
277 const std::string
&biasInitFilename
,
282 #endif /* GMX_AWH_H */