2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,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.
37 * \defgroup module_awh Accelerated weight histogram (AWH) method
38 * \ingroup module_applied_forces
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>
49 * \author Magnus Lundborg
52 /*! \libinternal \file
55 * Declares the Awh class.
57 * \author Viveca Lindahl
58 * \author Berk Hess <hess@kth.se>
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/utility/basedefinitions.h"
75 struct gmx_multisim_t;
83 enum class PbcType : int;
93 struct BiasCoupledToSystem;
95 class ForceWithVirial;
98 * \brief Coupling of the accelerated weight histogram method (AWH) with the system.
100 * AWH calculates the free energy along order parameters of the system.
101 * Free energy barriers are overcome by adaptively tuning a bias potential along
102 * the order parameter such that the biased distribution along the parameter
103 * converges toward a chosen target distribution.
105 * The Awh class takes care of the coupling between the system and the AWH
106 * bias(es). The Awh class contains one or more BiasCoupledToSystem objects.
107 * The BiasCoupledToSystem class takes care of the reaction coordinate input
108 * and force output for the single Bias object it containts.
110 * \todo Update parameter reading and checkpointing, when general C++ framework is ready.
115 /*! \brief Construct an AWH at the start of a simulation.
117 * AWH will here also register itself with the pull struct as the
118 * potential provider for the pull coordinates given as AWH coordinates
119 * in the user input. This allows AWH to later apply the bias force to
120 * these coordinate in \ref Awh::applyBiasForcesAndUpdateBias.
122 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
123 * \param[in] inputRecord General input parameters (as set up by grompp).
124 * \param[in] commRecord Struct for communication, can be nullptr.
125 * \param[in] multiSimRecord Multi-sim handler
126 * \param[in] awhParams AWH input parameters, consistent with the relevant
127 * parts of \p inputRecord (as set up by grompp).
128 * \param[in] biasInitFilename Name of file to read PMF and target from.
129 * \param[in,out] pull_work Pointer to a pull struct which AWH will
130 * couple to, has to be initialized, is assumed not to change during the
131 * lifetime of the Awh object.
132 * \param[in] numLambdaStates The number of free energy lambda states.
133 * \param[in] lambdaState The initial free energy lambda state of the system.
134 * \throws InvalidInputError If there is user input (or combinations thereof) that is not allowed.
137 const t_inputrec& inputRecord,
138 const t_commrec* commRecord,
139 const gmx_multisim_t* multiSimRecord,
140 const AwhParams& awhParams,
141 const std::string& biasInitFilename,
148 /*! \brief Peform an AWH update, to be called every MD step.
150 * An update has two tasks: apply the bias force and improve
151 * the bias and the free energy estimate that AWH keeps internally.
153 * For the first task, AWH retrieves the pull coordinate values from the pull struct.
154 * With these, the bias potential and forces are calculated.
155 * The bias force together with the atom forces and virial
156 * are passed on to pull which applies the bias force to the atoms.
157 * This is done at every step.
159 * Secondly, coordinate values are regularly sampled and kept by AWH.
160 * Convergence of the bias and free energy estimate is achieved by
161 * updating the AWH bias state after a certain number of samples has been collected.
163 * \note Requires that pull_potential from pull.h has been called first
164 * since AWH needs the current coordinate values (the pull code checks
167 * \param[in] pbcType Type of periodic boundary conditions.
168 * \param[in] masses Atoms masses.
169 * \param[in] neighborLambdaEnergies An array containing the energy of the system
170 * in neighboring lambdas. The array is of length numLambdas+1, where numLambdas is
171 * the number of free energy lambda states. Element 0 in the array is the energy
172 * of the current state and elements 1..numLambdas contain the energy of the system in the
173 * neighboring lambda states (also including the current state). When there are no free
174 * energy lambda state dimensions this can be empty.
175 * \param[in] neighborLambdaDhdl An array containing the dHdL at the neighboring lambda
176 * points. The array is of length numLambdas+1, where numLambdas is the number of free
177 * energy lambda states. Element 0 in the array is the dHdL
178 * of the current state and elements 1..numLambdas contain the dHdL of the system in the
179 * neighboring lambda states (also including the current state). When there are no free
180 * energy lambda state dimensions this can be empty.
181 * \param[in] box Box vectors.
182 * \param[in,out] forceWithVirial Force and virial buffers, should cover at least the local atoms.
184 * \param[in] step The current MD step.
185 * \param[in,out] wallcycle Wallcycle counter, can be nullptr.
186 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
187 * \returns the potential energy for the bias.
189 real applyBiasForcesAndUpdateBias(PbcType pbcType,
190 ArrayRef<const real> masses,
191 ArrayRef<const double> neighborLambdaEnergies,
192 ArrayRef<const double> neighborLambdaDhdl,
194 gmx::ForceWithVirial* forceWithVirial,
197 gmx_wallcycle* wallcycle,
201 * Update the AWH history in preparation for writing to checkpoint file.
203 * Should be called at least on the master rank at checkpoint steps.
205 * Should be called with a valid \p awhHistory (is checked).
207 * \param[in,out] awhHistory AWH history to set.
209 void updateHistory(AwhHistory* awhHistory) const;
212 * Allocate and initialize an AWH history with the given AWH state.
214 * This function should be called at the start of a new simulation
215 * at least on the master rank.
216 * Note that only constant data will be initialized here.
217 * History data is set by \ref Awh::updateHistory.
219 * \returns a shared pointer to the AWH history on the rank that does I/O, nullptr otherwise.
221 std::shared_ptr<AwhHistory> initHistoryFromState() const;
223 /*! \brief Restore the AWH state from the given history.
225 * Should be called on all ranks (for internal MPI broadcast).
226 * Should pass a point to an AwhHistory on the master rank that
227 * is compatible with the AWH setup in this simulation. Will throw
228 * an exception if it is not compatible.
230 * \param[in] awhHistory AWH history to restore from.
232 void restoreStateFromHistory(const AwhHistory* awhHistory);
234 /*! \brief Fills the AWH data block of an energy frame with data at certain steps.
236 * \param[in] step The current MD step.
237 * \param[in,out] fr Energy data frame.
239 void writeToEnergyFrame(int64_t step, t_enxframe* fr) const;
241 /*! \brief Returns string "AWH" for registering AWH as an external potential provider with the pull module.
243 static const char* externalPotentialString();
245 /*! \brief Register the AWH biased coordinates with pull.
247 * This function is public because it needs to be called by grompp
248 * (and is otherwise only called by Awh()).
249 * Pull requires all external potentials to register themselves
250 * before the end of pre-processing and before the first MD step.
251 * If this has not happened, pull with throw an error.
253 * \param[in] awhParams The AWH parameters.
254 * \param[in,out] pull_work Pull struct which AWH will register the bias into.
256 static void registerAwhWithPull(const AwhParams& awhParams, pull_t* pull_work);
258 /*! \brief Get the current free energy lambda state.
259 * \returns The value of lambdaState_.
261 int fepLambdaState() const { return fepLambdaState_; }
263 /*! \brief Returns if foreign energy differences are required during this step.
264 * \param[in] step The current MD step.
266 bool needForeignEnergyDifferences(int64_t step) const;
268 /*! \brief Returns true if AWH has a bias with a free energy lambda state dimension at all.
270 bool hasFepLambdaDimension() const;
273 /*! \brief Returns whether we need to write output at the current step.
275 * \param[in] step The current MD step.
277 bool isOutputStep(int64_t step) const;
279 std::vector<BiasCoupledToSystem> biasCoupledToSystem_; /**< AWH biases and definitions of their coupling to the system. */
280 const int64_t seed_; /**< Random seed for MC jumping with umbrella type bias potential. */
281 const int nstout_; /**< Interval in steps for writing to energy file. */
282 const t_commrec* commRecord_; /**< Pointer to the communication record. */
283 //! Object for sharing bias between simulations, only set when needed
284 std::unique_ptr<BiasSharing> biasSharing_;
285 pull_t* pull_; /**< Pointer to the pull working data. */
286 double potentialOffset_; /**< The offset of the bias potential which changes due to bias updates. */
287 const int numFepLambdaStates_; /**< The number of free energy lambda states of the system. */
288 int fepLambdaState_; /**< The current free energy lambda state. */
291 /*! \brief Makes an Awh and prepares to use it if the user input
294 * Restores state from history in checkpoint if needed.
296 * \param[in,out] fplog General output file, normally md.log, can be nullptr.
297 * \param[in] inputRecord General input parameters (as set up by grompp).
298 * \param[in] stateGlobal A pointer to the global state structure.
299 * \param[in] commRecord Struct for communication, can be nullptr.
300 * \param[in] multiSimRecord Multi-sim handler
301 * \param[in] startingFromCheckpoint Whether the simulation is starting from a checkpoint
302 * \param[in] usingShellParticles Whether the user requested shell particles (which is unsupported)
303 * \param[in] biasInitFilename Name of file to read PMF and target from.
304 * \param[in,out] pull_work Pointer to a pull struct which AWH will couple to, has to be initialized,
305 * is assumed not to change during the lifetime of the Awh object.
306 * \returns An initialized Awh module, or nullptr if none was requested.
307 * \throws InvalidInputError If another active module is not supported.
309 std::unique_ptr<Awh> prepareAwhModule(FILE* fplog,
310 const t_inputrec& inputRecord,
311 t_state* stateGlobal,
312 const t_commrec* commRecord,
313 const gmx_multisim_t* multiSimRecord,
314 bool startingFromCheckpoint,
315 bool usingShellParticles,
316 const std::string& biasInitFilename,
321 #endif /* GMX_AWH_H */