2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019,2020, 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.
39 * Declares the BiasState class.
41 * The data members of this class are the state variables of the bias.
42 * All interaction from the outside happens through the Bias class, which
43 * holds important helper classes such as DimParams and BiasGrid.
44 * This class holds many methods, but more are const methods that compute
45 * properties of the state.
47 * \author Viveca Lindahl
48 * \author Berk Hess <hess@kth.se>
52 #ifndef GMX_AWH_BIASSTATE_H
53 #define GMX_AWH_BIASSTATE_H
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/utility/alignedallocator.h"
61 #include "gromacs/utility/arrayref.h"
62 #include "gromacs/utility/basedefinitions.h"
63 #include "gromacs/utility/gmxassert.h"
65 #include "coordstate.h"
66 #include "dimparams.h"
67 #include "histogramsize.h"
69 struct gmx_multisim_t;
75 struct AwhBiasHistory;
83 * \brief The state of a bias.
85 * The bias state has the current coordinate state: its value and the grid point
86 * it maps to (the grid point of the umbrella potential if needed). It contains
87 * a vector with the state for each point on the grid. It also
88 * counts the number of updates issued and tracks which points have been sampled
89 * since last update. Finally, the convergence state is a global property set
90 * ultimately by the histogram size histogramSize in the sub-class HistogramSize,
91 * since the update sizes are ~ 1/histogramSize.
96 /*! \brief Constructor.
98 * Constructs the global state and the point states on a provided
99 * geometric grid passed in \p grid.
101 * \param[in] awhBiasParams The Bias parameters from inputrec.
102 * \param[in] histogramSizeInitial The estimated initial histogram size.
103 * This is floating-point, since histograms use weighted
104 * entries and grow by a floating-point scaling factor.
105 * \param[in] dimParams The dimension parameters.
106 * \param[in] grid The bias grid.
108 BiasState(const AwhBiasParams& awhBiasParams,
109 double histogramSizeInitial,
110 const std::vector<DimParams>& dimParams,
111 const BiasGrid& grid);
114 * Restore the bias state from history.
116 * \param[in] biasHistory Bias history struct.
117 * \param[in] grid The bias grid.
119 void restoreFromHistory(const AwhBiasHistory& biasHistory, const BiasGrid& grid);
122 * Broadcast the bias state over the MPI ranks in this simulation.
124 * \param[in] commRecord Struct for communication.
126 void broadcast(const t_commrec* commRecord);
129 * Allocate and initialize a bias history with the given bias state.
131 * This function will be called at the start of a new simulation.
132 * Note that this only sets the correct size and does produce
133 * a valid history object, but with all data set to zero.
134 * Actual history data is set by \ref updateHistory.
136 * \param[in,out] biasHistory AWH history to initialize.
138 void initHistoryFromState(AwhBiasHistory* biasHistory) const;
141 * Update the bias state history with the current state.
143 * \param[out] biasHistory Bias history struct.
144 * \param[in] grid The bias grid.
146 void updateHistory(AwhBiasHistory* biasHistory, const BiasGrid& grid) const;
149 /*! \brief Convolves the given PMF using the given AWH bias.
151 * \note: The PMF is in single precision, because it is a statistical
152 * quantity and therefore never reaches full float precision.
154 * \param[in] dimParams The bias dimensions parameters
155 * \param[in] grid The grid.
156 * \param[in,out] convolvedPmf Array returned will be of the same length as the AWH grid to store the convolved PMF in.
158 void calcConvolvedPmf(const std::vector<DimParams>& dimParams,
159 const BiasGrid& grid,
160 std::vector<float>* convolvedPmf) const;
163 * Convolves the PMF and sets the initial free energy to its convolution.
165 * \param[in] dimParams The bias dimensions parameters
166 * \param[in] grid The bias grid.
168 void setFreeEnergyToConvolvedPmf(const std::vector<DimParams>& dimParams, const BiasGrid& grid);
171 * Normalize the PMF histogram.
173 * \param[in] numSharingSims The number of simulations sharing the bias.
175 void normalizePmf(int numSharingSims);
179 * Initialize the state of grid coordinate points.
181 * \param[in] awhBiasParams Bias parameters from inputrec.
182 * \param[in] dimParams The dimension parameters.
183 * \param[in] grid The grid.
184 * \param[in] params The bias parameters.
185 * \param[in] filename Name of file to read PMF and target from.
186 * \param[in] numBias The number of biases.
188 void initGridPointState(const AwhBiasParams& awhBiasParams,
189 const std::vector<DimParams>& dimParams,
190 const BiasGrid& grid,
191 const BiasParams& params,
192 const std::string& filename,
196 * Performs statistical checks on the collected histograms and warns if issues are detected.
198 * \param[in] grid The grid.
199 * \param[in] biasIndex The index of the bias we are checking for.
201 * \param[in,out] fplog Output file for warnings.
202 * \param[in] maxNumWarnings Don't issue more than this number of warnings.
203 * \returns the number of warnings issued.
205 int warnForHistogramAnomalies(const BiasGrid& grid, int biasIndex, double t, FILE* fplog, int maxNumWarnings) const;
208 * Calculates and sets the force the coordinate experiences from an umbrella centered at the given point.
210 * The umbrella potential is an harmonic potential given by 0.5k(coord value - point value)^2. This
211 * value is also returned.
213 * \param[in] dimParams The bias dimensions parameters.
214 * \param[in] grid The grid.
215 * \param[in] point Point for umbrella center.
216 * \param[in,out] force Force vector to set.
217 * Returns the umbrella potential.
219 double calcUmbrellaForceAndPotential(const std::vector<DimParams>& dimParams,
220 const BiasGrid& grid,
222 gmx::ArrayRef<double> force) const;
225 * Calculates and sets the convolved force acting on the coordinate.
227 * The convolved force is the weighted sum of forces from umbrellas
228 * located at each point in the grid.
230 * \param[in] dimParams The bias dimensions parameters.
231 * \param[in] grid The grid.
232 * \param[in] probWeightNeighbor Probability weights of the neighbors.
233 * \param[in] forceWorkBuffer Force work buffer, values only used internally.
234 * \param[in,out] force Bias force vector to set.
236 void calcConvolvedForce(const std::vector<DimParams>& dimParams,
237 const BiasGrid& grid,
238 gmx::ArrayRef<const double> probWeightNeighbor,
239 gmx::ArrayRef<double> forceWorkBuffer,
240 gmx::ArrayRef<double> force) const;
243 * Move the center point of the umbrella potential.
245 * A new umbrella center is sampled from the biased distibution. Also, the bias
246 * force is updated and the new potential is return.
248 * This function should only be called when the bias force is not being convolved.
249 * It is assumed that the probability distribution has been updated.
251 * \param[in] dimParams Bias dimension parameters.
252 * \param[in] grid The grid.
253 * \param[in] probWeightNeighbor Probability weights of the neighbors.
254 * \param[in,out] biasForce The AWH bias force.
255 * \param[in] step Step number, needed for the random number generator.
256 * \param[in] seed Random seed.
257 * \param[in] indexSeed Second random seed, should be the bias Index.
258 * \returns the new potential value.
260 double moveUmbrella(const std::vector<DimParams>& dimParams,
261 const BiasGrid& grid,
262 gmx::ArrayRef<const double> probWeightNeighbor,
263 gmx::ArrayRef<double> biasForce,
270 * Gets the histogram rescaling factors needed for skipped updates.
272 * \param[in] params The bias parameters.
273 * \param[out] weighthistScaling Scaling factor for the reference weight histogram.
274 * \param[out] logPmfsumScaling Log of the scaling factor for the PMF histogram.
276 void getSkippedUpdateHistogramScaleFactors(const BiasParams& params,
277 double* weighthistScaling,
278 double* logPmfsumScaling) const;
282 * Do all previously skipped updates.
283 * Public for use by tests.
285 * \param[in] params The bias parameters.
287 void doSkippedUpdatesForAllPoints(const BiasParams& params);
290 * Do previously skipped updates in this neighborhood.
292 * \param[in] params The bias parameters.
293 * \param[in] grid The grid.
295 void doSkippedUpdatesInNeighborhood(const BiasParams& params, const BiasGrid& grid);
299 * Reset the range used to make the local update list.
301 * \param[in] grid The grid.
303 void resetLocalUpdateRange(const BiasGrid& grid);
306 * Returns the new size of the reference weight histogram in the initial stage.
308 * This function also takes care resetting the histogram used for covering checks
309 * and for exiting the initial stage.
311 * \param[in] params The bias parameters.
313 * \param[in] detectedCovering True if we detected that the sampling interval has been sufficiently covered.
314 * \param[in,out] fplog Log file.
315 * \returns the new histogram size.
317 double newHistogramSizeInitialStage(const BiasParams& params, double t, bool detectedCovering, FILE* fplog);
320 * Check if the sampling region has been covered "enough" or not.
322 * A one-dimensional interval is defined as covered if each point has
323 * accumulated the same weight as is in the peak of a discretized normal
324 * distribution. For multiple dimensions, the weights are simply projected
325 * onto each dimension and the multidimensional space is covered if each
328 * \note The covering criterion for multiple dimensions could improved, e.g.
329 * by using a path finding algorithm.
331 * \param[in] params The bias parameters.
332 * \param[in] dimParams Bias dimension parameters.
333 * \param[in] grid The grid.
334 * \param[in] commRecord Struct for intra-simulation communication.
335 * \param[in] multiSimComm Struct for multi-simulation communication.
336 * \returns true if covered.
338 bool isSamplingRegionCovered(const BiasParams& params,
339 const std::vector<DimParams>& dimParams,
340 const BiasGrid& grid,
341 const t_commrec* commRecord,
342 const gmx_multisim_t* multiSimComm) const;
345 * Return the new reference weight histogram size for the current update.
347 * This function also takes care of checking for covering in the initial stage.
349 * \param[in] params The bias parameters.
351 * \param[in] covered True if the sampling interval has been covered enough.
352 * \param[in,out] fplog Log file.
353 * \returns the new histogram size.
355 double newHistogramSize(const BiasParams& params, double t, bool covered, FILE* fplog);
359 * Update the reaction coordinate value.
361 * \param[in] grid The bias grid.
362 * \param[in] coordValue The current reaction coordinate value (there are no limits on allowed values).
364 void setCoordValue(const BiasGrid& grid, const awh_dvec coordValue)
366 coordState_.setCoordValue(grid, coordValue);
370 * Performs an update of the bias.
372 * The objective of the update is to use collected samples (probability weights)
373 * to improve the free energy estimate. For sake of efficiency, the update is
374 * local whenever possible, meaning that only points that have actually been sampled
375 * are accessed and updated here. For certain AWH settings or at certain steps
376 * however, global need to be performed. Besides the actual free energy update, this
377 * function takes care of ensuring future convergence of the free energy. Convergence
378 * is obtained by increasing the size of the reference weight histogram in a controlled
379 * (sometimes dynamic) manner. Also, there are AWH variables that are direct functions
380 * of the free energy or sampling history that need to be updated here, namely the target
381 * distribution and the bias function.
383 * \param[in] dimParams The dimension parameters.
384 * \param[in] grid The grid.
385 * \param[in] params The bias parameters.
386 * \param[in] commRecord Struct for intra-simulation communication.
387 * \param[in] ms Struct for multi-simulation communication.
389 * \param[in] step Time step.
390 * \param[in,out] fplog Log file.
391 * \param[in,out] updateList Work space to store a temporary list.
393 void updateFreeEnergyAndAddSamplesToHistogram(const std::vector<DimParams>& dimParams,
394 const BiasGrid& grid,
395 const BiasParams& params,
396 const t_commrec* commRecord,
397 const gmx_multisim_t* ms,
401 std::vector<int>* updateList);
404 * Update the probability weights and the convolved bias.
406 * Given a coordinate value, each grid point is assigned a probability
407 * weight, w(point|value), that depends on the current bias function. The sum
408 * of these weights is needed for normalizing the probability sum to 1 but
409 * also equals the effective, or convolved, biasing weight for this coordinate
410 * value. The convolved bias is needed e.g. for extracting the PMF, so we save
411 * it here since this saves us from doing extra exponential function evaluations
414 * \param[in] dimParams The bias dimensions parameters
415 * \param[in] grid The grid.
416 * \param[out] weight Probability weights of the neighbors, SIMD aligned.
417 * \returns the convolved bias.
420 double updateProbabilityWeightsAndConvolvedBias(const std::vector<DimParams>& dimParams,
421 const BiasGrid& grid,
422 std::vector<double, AlignedAllocator<double>>* weight) const;
425 * Take samples of the current probability weights for future updates and analysis.
427 * Points in the current neighborhood will now have data meaning they
428 * need to be included in the local update list of the next update.
429 * Therefore, the local update range is also update here.
431 * \param[in] grid The grid.
432 * \param[in] probWeightNeighbor Probability weights of the neighbors.
434 void sampleProbabilityWeights(const BiasGrid& grid, gmx::ArrayRef<const double> probWeightNeighbor);
437 * Sample the reaction coordinate and PMF for future updates or analysis.
439 * These samples do not affect the (future) sampling and are thus
440 * pure observables. Statisics of these are stored in the energy file.
442 * \param[in] grid The grid.
443 * \param[in] probWeightNeighbor Probability weights of the neighbors.
444 * \param[in] convolvedBias The convolved bias.
446 void sampleCoordAndPmf(const BiasGrid& grid,
447 gmx::ArrayRef<const double> probWeightNeighbor,
448 double convolvedBias);
450 * Calculates the convolved bias for a given coordinate value.
452 * The convolved bias is the effective bias acting on the coordinate.
453 * Since the bias here has arbitrary normalization, this only makes
454 * sense as a relative, to other coordinate values, measure of the bias.
456 * \note If it turns out to be costly to calculate this pointwise
457 * the convolved bias for the whole grid could be returned instead.
459 * \param[in] dimParams The bias dimensions parameters
460 * \param[in] grid The grid.
461 * \param[in] coordValue Coordinate value.
462 * \returns the convolved bias >= -GMX_FLOAT_MAX.
464 double calcConvolvedBias(const std::vector<DimParams>& dimParams,
465 const BiasGrid& grid,
466 const awh_dvec& coordValue) const;
469 * Fills the given array with PMF values.
471 * Points outside of the biasing target region will get PMF = GMX_FLOAT_MAX.
472 * \note: The PMF is in single precision, because it is a statistical
473 * quantity and therefore never reaches full float precision.
475 * \param[out] pmf Array(ref) to be filled with the PMF values, should have the same size as the bias grid.
477 void getPmf(gmx::ArrayRef<float> /*pmf*/) const;
479 /*! \brief Returns the current coordinate state.
481 const CoordState& coordState() const { return coordState_; }
483 /*! \brief Returns a const reference to the point state.
485 const std::vector<PointState>& points() const { return points_; }
487 /*! \brief Returns true if we are in the initial stage.
489 bool inInitialStage() const { return histogramSize_.inInitialStage(); }
491 /*! \brief Returns the current histogram size.
493 inline HistogramSize histogramSize() const { return histogramSize_; }
497 CoordState coordState_; /**< The Current coordinate state */
499 /* The grid point state */
500 std::vector<PointState> points_; /**< Vector of state of the grid points */
502 /* Covering values for each point on the grid */
503 std::vector<double> weightSumCovering_; /**< Accumulated weights for covering checks */
505 HistogramSize histogramSize_; /**< Global histogram size related values. */
507 /* Track the part of the grid sampled since the last update. */
508 awh_ivec originUpdatelist_; /**< The origin of the rectangular region that has been sampled since last update. */
509 awh_ivec endUpdatelist_; /**< The end of the rectangular region that has been sampled since last update. */
512 //! Linewidth used for warning output
513 static const int c_linewidth = 80 - 2;
515 //! Indent used for warning output
516 static const int c_indent = 0;
520 #endif /* GMX_AWH_BIASSTATE_H */