Add reading and writing to AWH module
[gromacs.git] / src / gromacs / awh / bias.h
blob841f08098525ae75522881db6d2161b1de4569df
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017, 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 /*! \internal \file
38 * \brief
39 * Declares the Bias class.
41 * This class is essentially a wrapper around the BiasState class.
42 * In addition to BiasState, it holds all data that BiasState needs
43 * to update the bias. Interaction of the outside world, such as updating
44 * BiasState or extracting bias data all happen through Bias.
46 * \author Viveca Lindahl
47 * \author Berk Hess <hess@kth.se>
48 * \ingroup module_awh
51 #ifndef GMX_AWH_BIAS_H
52 #define GMX_AWH_BIAS_H
54 #include <memory>
55 #include <vector>
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/utility/alignedallocator.h"
59 #include "gromacs/utility/basedefinitions.h"
60 #include "gromacs/utility/gmxassert.h"
62 #include "biasparams.h"
63 #include "biasstate.h"
64 #include "biaswriter.h"
65 #include "dimparams.h"
66 #include "grid.h"
68 struct gmx_multisim_t;
69 struct t_commrec;
70 struct t_enxsubblock;
72 namespace gmx
75 struct AwhBiasHistory;
76 struct AwhBiasParams;
77 struct AwhHistory;
78 struct AwhParams;
79 struct AwhPointStateHistory;
80 class Grid;
81 class GridAxis;
82 class PointState;
84 /*! \internal
85 * \brief A bias acting on a multidimensional coordinate.
87 * At each step AWH should provide its biases with updated
88 * values of their coordinates. Each bias provides AWH with an updated
89 * bias forces and the corresponding potential.
91 * See the user manual for details on the algorithm and equations.
93 * The bias is responsible for keeping and updating a free energy estimate
94 * along the coordinate. The bias potential is basically a function of the
95 * free energy estimate and so also changes by the update.
96 * The free energy update is based on information from coordinate samples
97 * collected at a constant bias potential, between updates.
99 * The bias keeps a grid with coordinate points that organizes spatial
100 * information about the coordinate. The grid has the the same geometry
101 * as the coordinate, i.e. they have the same dimensionality and periodicity
102 * (if any). The number of points in the grid sets the resolution of
103 * the collected data and its extent defines the sampling region of interest.
105 * Each coordinate point has further statistical properties and function values
106 * which a grid point does not know about. E.g., for the bias each coordinate point
107 * is associated with values of the bias, free energy and target distribution,
108 * accumulated sampling weight, etc. For this the bias attaches to each grid
109 * point a state. The grid + vector of point states are the bias coordinate points.
111 * The bias has a fairly complex global state keeping track of where
112 * the system (coordinate) currently is (CoordState), where it has
113 * sampled since the last update (BiasState) and controlling the free energy
114 * convergence rate (HistogramSize).
116 * Partly, the complexity comes from the bias having two convergence stages:
117 * an initial stage which in an heuristic, non-deterministic way restricts
118 * the early convergence rate for sake of robustness; and a final stage
119 * where the convergence rate is constant. The length of the initial stage
120 * depends on the sampling and is unknown beforehand.
122 * Another complexity comes from the fact that coordinate points,
123 * for sake of efficiency in the case of many grid points, are typically
124 * only accessed in recently sampled regions even though the free energy
125 * update is inherently global and affects all points.
126 * The bias allows points thay are non-local at the time the update
127 * was issued to postpone ("skip", as it is called in the code) the update.
128 * A non-local point is defined as a point which has not been sampled since
129 * the last update. Local points are points that have been sampled since
130 * the last update. The (current) set of local points are kept track of by
131 * the bias state and reset after every update. An update is called local
132 * if it only updates local points. Non-local points will temporarily "skip"
133 * the update until next time they are local (or when a global update
134 * is issued). For this to work, the bias keeps a global "clock"
135 * (in HistogramSize) of the number of issued updates. Each PointState
136 * also has its own local "clock" with the counting the number of updates
137 * it has pulled through. When a point updates its state it asserts
138 * that its local clock is synchronized with the global clock.
140 class Bias
142 public:
143 //! Enum for requesting Bias set up with(out) I/O on this rank.
144 enum class ThisRankWillDoIO
146 No, //!< This rank will not do I/O.
147 Yes //!< This rank will do I/O.
150 /*! \brief
151 * Constructor.
153 * \param[in] biasIndexInCollection Index of the bias in collection.
154 * \param[in] awhParams AWH parameters.
155 * \param[in] awhBiasParams Bias parameters.
156 * \param[in] dimParams Bias dimension parameters.
157 * \param[in] beta 1/(k_B T).
158 * \param[in] mdTimeStep The MD time step.
159 * \param[in] numSharingSimulations The number of simulations to share the bias across.
160 * \param[in] biasInitFilename Name of file to read PMF and target from.
161 * \param[in] thisRankWillDoIO Tells whether this MPI rank will do I/O (checkpointing, AWH output), normally (only) the master rank does I/O.
162 * \param[in] disableUpdateSkips If to disable update skips, useful for testing.
164 Bias(int biasIndexInCollection,
165 const AwhParams &awhParams,
166 const AwhBiasParams &awhBiasParams,
167 const std::vector<DimParams> &dimParams,
168 double beta,
169 double mdTimeStep,
170 int numSharingSimulations,
171 const std::string &biasInitFilename,
172 ThisRankWillDoIO thisRankWillDoIO,
173 BiasParams::DisableUpdateSkips disableUpdateSkips = BiasParams::DisableUpdateSkips::no);
175 /*! \brief
176 * Evolves the bias at every step.
178 * At each step the bias step needs to:
179 * - set the bias force and potential;
180 * - update the free energy and bias if needed;
181 * - reweight samples to extract the PMF.
183 * \param[in] coordValue The current coordinate value(s).
184 * \param[out] biasForce The bias force.
185 * \param[out] awhPotential Bias potential.
186 * \param[out] potentialJump Change in bias potential for this bias.
187 * \param[in] ms Struct for multi-simulation communication.
188 * \param[in] t Time.
189 * \param[in] step Time step.
190 * \param[in] seed Random seed.
191 * \param[in,out] fplog Log file.
193 void calcForceAndUpdateBias(const awh_dvec coordValue,
194 awh_dvec biasForce,
195 double *awhPotential,
196 double *potentialJump,
197 const gmx_multisim_t *ms,
198 double t,
199 gmx_int64_t step,
200 gmx_int64_t seed,
201 FILE *fplog);
203 /*! \brief
204 * Calculates the convolved bias for a given coordinate value.
206 * The convolved bias is the effective bias acting on the coordinate.
207 * Since the bias here has arbitrary normalization, this only makes
208 * sense as a relative, to other coordinate values, measure of the bias.
210 * \param[in] coordValue The coordinate value.
211 * \returns the convolved bias >= -GMX_FLOAT_MAX.
213 double calcConvolvedBias(const awh_dvec &coordValue) const
215 return state_.calcConvolvedBias(dimParams_, grid_, coordValue);
218 /*! \brief
219 * Restore the bias state from history on the master rank and broadcast it.
221 * \param[in] biasHistory Bias history struct, only allowed to be nullptr on non-master ranks.
222 * \param[in] cr The communication record.
224 void restoreStateFromHistory(const AwhBiasHistory *biasHistory,
225 const t_commrec *cr);
227 /*! \brief
228 * Allocate and initialize a bias history with the given bias state.
230 * This function will be called at the start of a new simulation.
231 * Note that only constant data will be initialized here.
232 * History data is set by \ref updateHistory.
234 * \param[in,out] biasHistory AWH history to initialize.
236 void initHistoryFromState(AwhBiasHistory *biasHistory) const;
238 /*! \brief
239 * Update the bias history with the current state.
241 * \param[out] biasHistory Bias history struct.
243 void updateHistory(AwhBiasHistory *biasHistory) const;
245 /*! \brief
246 * Do all previously skipped updates.
247 * Public for use by tests.
249 void doSkippedUpdatesForAllPoints();
251 //! Returns the dimensionality of the bias.
252 inline int ndim() const
254 return dimParams_.size();
257 /*! \brief Returns the dimension parameters.
259 inline const std::vector<DimParams> &dimParams() const
261 return dimParams_;
264 //! Returns the bias parameters
265 inline const BiasParams &params() const
267 return params_;
270 //! Returns the global state of the bias.
271 inline const BiasState &state() const
273 return state_;
276 //! Returns the index of the bias.
277 inline int biasIndex() const
279 return params_.biasIndex;
282 /*! \brief Return the coordinate value for a grid point.
284 * \param[in] gridPointIndex The index of the grid point.
286 inline const awh_dvec &getGridCoordValue(size_t gridPointIndex) const
288 GMX_ASSERT(gridPointIndex < grid_.numPoints(), "gridPointIndex should be in the range of the grid");
290 return grid_.point(gridPointIndex).coordValue;
293 private:
294 /*! \brief
295 * Performs statistical checks on the collected histograms and warns if issues are detected.
297 * \param[in] t Time.
298 * \param[in] step Time step.
299 * \param[in,out] fplog Output file for warnings.
301 void warnForHistogramAnomalies(double t,
302 gmx_int64_t step,
303 FILE *fplog);
305 public:
306 /*! \brief Return the number of data blocks that have been prepared for writing.
308 int numEnergySubblocksToWrite() const;
310 /*! \brief Write bias data blocks to energy subblocks.
312 * \param[in,out] subblock Energy subblocks to write to.
313 * \returns the number of subblocks written.
315 int writeToEnergySubblocks(t_enxsubblock *subblock) const;
317 /* Data members. */
318 private:
319 const std::vector<DimParams> dimParams_; /**< Parameters for each dimension. */
320 const Grid grid_; /**< The multidimensional grid organizing the coordinate point locations. */
322 const BiasParams params_; /**< Constant parameters for the method. */
324 BiasState state_; /**< The state, both global and of the grid points */
325 std::vector<int> updateList_; /**< List of points for update for temporary use (could be made another tempWorkSpace) */
327 const bool thisRankDoesIO_; /**< Tells whether this MPI rank will do I/O (checkpointing, AWH output) */
329 /* I/O */
330 std::unique_ptr<BiasWriter> writer_; /**< Takes care of AWH data output. */
332 /* Temporary working vector used during the update.
333 * This only here to avoid allocation at every MD step.
335 std::vector < double, AlignedAllocator < double>> alignedTempWorkSpace_; /**< Working vector of doubles. */
337 /* Run-local counter to avoid flooding log with warnings. */
338 int numWarningsIssued_; /**< The number of warning issued in the current run. */
341 } // namespace gmx
343 #endif /* GMX_AWH_BIAS_H */