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.
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
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/utility/alignedallocator.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"
76 struct AwhBiasHistory;
84 * \brief The state of a bias.
86 * The bias state has the current coordinate state: its value and the grid point
87 * it maps to (the grid point of the umbrella potential if needed). It contains
88 * a vector with the state for each point on the grid. It also
89 * counts the number of updates issued and tracks which points have been sampled
90 * since last update. Finally, the convergence state is a global property set
91 * ultimately by the histogram size histogramSize in the sub-class HistogramSize,
92 * since the update sizes are ~ 1/histogramSize.
97 /*! \brief Constructor.
99 * Constructs the global state and the point states on a provided
100 * geometric grid passed in \p grid.
102 * \param[in] awhBiasParams The Bias parameters from inputrec.
103 * \param[in] histogramSizeInitial The estimated initial histogram size.
104 * This is floating-point, since histograms use weighted
105 * entries and grow by a floating-point scaling factor.
106 * \param[in] dimParams The dimension parameters.
107 * \param[in] grid The bias grid.
108 * \param[in] biasSharing Multisim bias sharing object, can be nullptrx
110 BiasState(const AwhBiasParams& awhBiasParams,
111 double histogramSizeInitial,
112 ArrayRef<const DimParams> dimParams,
113 const BiasGrid& grid,
114 const BiasSharing* biasSharing);
117 * Restore the bias state from history.
119 * \param[in] biasHistory Bias history struct.
120 * \param[in] grid The bias grid.
122 void restoreFromHistory(const AwhBiasHistory& biasHistory, const BiasGrid& grid);
125 * Broadcast the bias state over the MPI ranks in this simulation.
127 * \param[in] commRecord Struct for communication.
129 void broadcast(const t_commrec* commRecord);
132 * Allocate and initialize a bias history with the given bias state.
134 * This function will be called at the start of a new simulation.
135 * Note that this only sets the correct size and does produce
136 * a valid history object, but with all data set to zero.
137 * Actual history data is set by \ref updateHistory.
139 * \param[in,out] biasHistory AWH history to initialize.
141 void initHistoryFromState(AwhBiasHistory* biasHistory) const;
144 * Update the bias state history with the current state.
146 * \param[out] biasHistory Bias history struct.
147 * \param[in] grid The bias grid.
149 void updateHistory(AwhBiasHistory* biasHistory, const BiasGrid& grid) const;
152 /*! \brief Convolves the given PMF using the given AWH bias.
154 * \note: The PMF is in single precision, because it is a statistical
155 * quantity and therefore never reaches full float precision.
157 * \param[in] dimParams The bias dimensions parameters
158 * \param[in] grid The grid.
159 * \param[in,out] convolvedPmf Array returned will be of the same length as the AWH grid to store the convolved PMF in.
161 void calcConvolvedPmf(ArrayRef<const DimParams> dimParams,
162 const BiasGrid& grid,
163 std::vector<float>* convolvedPmf) const;
166 * Convolves the PMF and sets the initial free energy to its convolution.
168 * \param[in] dimParams The bias dimensions parameters
169 * \param[in] grid The bias grid.
171 void setFreeEnergyToConvolvedPmf(ArrayRef<const DimParams> dimParams, const BiasGrid& grid);
174 * Normalize the PMF histogram.
176 * \param[in] numSharingSims The number of simulations sharing the bias.
178 void normalizePmf(int numSharingSims);
182 * Initialize the state of grid coordinate points.
184 * \param[in] awhBiasParams Bias parameters from inputrec.
185 * \param[in] dimParams The dimension parameters.
186 * \param[in] grid The grid.
187 * \param[in] params The bias parameters.
188 * \param[in] filename Name of file to read PMF and target from.
189 * \param[in] numBias The number of biases.
191 void initGridPointState(const AwhBiasParams& awhBiasParams,
192 ArrayRef<const DimParams> dimParams,
193 const BiasGrid& grid,
194 const BiasParams& params,
195 const std::string& filename,
199 * Performs statistical checks on the collected histograms and warns if issues are detected.
201 * \param[in] grid The grid.
202 * \param[in] biasIndex The index of the bias we are checking for.
204 * \param[in,out] fplog Output file for warnings.
205 * \param[in] maxNumWarnings Don't issue more than this number of warnings.
206 * \returns the number of warnings issued.
208 int warnForHistogramAnomalies(const BiasGrid& grid, int biasIndex, double t, FILE* fplog, int maxNumWarnings) const;
211 * Calculates and sets the force the coordinate experiences from an umbrella centered at the given point.
213 * The umbrella potential is an harmonic potential given by 0.5k(coord value - point value)^2. This
214 * value is also returned.
216 * \param[in] dimParams The bias dimensions parameters.
217 * \param[in] grid The grid.
218 * \param[in] point Point for umbrella center.
219 * \param[in] neighborLambdaDhdl An array containing the dHdL at the neighboring lambda
220 * points. The array is of length numLambdas+1, where numLambdas is the number of free
221 * energy lambda states. Element 0 in the array is the dHdL
222 * of the current state and elements 1..numLambdas contain the dHdL of the system in the
223 * neighboring lambda states (also including the current state). When there are no free
224 * energy lambda state dimensions this can be empty.
225 * \param[in,out] force Force vector to set.
226 * Returns the umbrella potential.
228 double calcUmbrellaForceAndPotential(ArrayRef<const DimParams> dimParams,
229 const BiasGrid& grid,
231 ArrayRef<const double> neighborLambdaDhdl,
232 ArrayRef<double> force) const;
235 * Calculates and sets the convolved force acting on the coordinate.
237 * The convolved force is the weighted sum of forces from umbrellas
238 * located at each point in the grid.
240 * \param[in] dimParams The bias dimensions parameters.
241 * \param[in] grid The grid.
242 * \param[in] probWeightNeighbor Probability weights of the neighbors.
243 * \param[in] neighborLambdaDhdl An array containing the dHdL at the neighboring lambda
244 * points. The array is of length numLambdas+1, where numLambdas is the number of free
245 * energy lambda states. Element 0 in the array is the dHdL
246 * of the current state and elements 1..numLambdas contain the dHdL of the system in the
247 * neighboring lambda states (also including the current state). When there are no free
248 * energy lambda state dimensions this can be empty.
249 * \param[in] forceWorkBuffer Force work buffer, values only used internally.
250 * \param[in,out] force Bias force vector to set.
252 void calcConvolvedForce(ArrayRef<const DimParams> dimParams,
253 const BiasGrid& grid,
254 ArrayRef<const double> probWeightNeighbor,
255 ArrayRef<const double> neighborLambdaDhdl,
256 ArrayRef<double> forceWorkBuffer,
257 ArrayRef<double> force) const;
260 * Move the center point of the umbrella potential.
262 * A new umbrella center is sampled from the biased distibution. Also, the bias
263 * force is updated and the new potential is return.
265 * This function should only be called when the bias force is not being convolved.
266 * It is assumed that the probability distribution has been updated.
268 * \param[in] dimParams Bias dimension parameters.
269 * \param[in] grid The grid.
270 * \param[in] probWeightNeighbor Probability weights of the neighbors.
271 * \param[in] neighborLambdaDhdl An array containing the dHdL at the neighboring lambda
272 * points. The array is of length numLambdas+1, where numLambdas is the number of free
273 * energy lambda states. Element 0 in the array is the dHdL
274 * of the current state and elements 1..numLambdas contain the dHdL of the system in the
275 * neighboring lambda states (also including the current state). When there are no free
276 * energy lambda state dimensions this can be empty.
277 * \param[in,out] biasForce The AWH bias force.
278 * \param[in] step Step number, needed for the random number generator.
279 * \param[in] seed Random seed.
280 * \param[in] indexSeed Second random seed, should be the bias Index.
281 * \param[in] onlySampleUmbrellaGridpoint Only sample the umbrella gridpoint without calculating
282 * force and potential.
283 * \returns the new potential value.
285 double moveUmbrella(ArrayRef<const DimParams> dimParams,
286 const BiasGrid& grid,
287 ArrayRef<const double> probWeightNeighbor,
288 ArrayRef<const double> neighborLambdaDhdl,
289 ArrayRef<double> biasForce,
293 bool onlySampleUmbrellaGridpoint);
297 * Gets the histogram rescaling factors needed for skipped updates.
299 * \param[in] params The bias parameters.
300 * \param[out] weighthistScaling Scaling factor for the reference weight histogram.
301 * \param[out] logPmfsumScaling Log of the scaling factor for the PMF histogram.
303 void getSkippedUpdateHistogramScaleFactors(const BiasParams& params,
304 double* weighthistScaling,
305 double* logPmfsumScaling) const;
309 * Do all previously skipped updates.
310 * Public for use by tests.
312 * \param[in] params The bias parameters.
314 void doSkippedUpdatesForAllPoints(const BiasParams& params);
317 * Do previously skipped updates in this neighborhood.
319 * \param[in] params The bias parameters.
320 * \param[in] grid The grid.
322 void doSkippedUpdatesInNeighborhood(const BiasParams& params, const BiasGrid& grid);
326 * Reset the range used to make the local update list.
328 * \param[in] grid The grid.
330 void resetLocalUpdateRange(const BiasGrid& grid);
333 * Returns the new size of the reference weight histogram in the initial stage.
335 * This function also takes care resetting the histogram used for covering checks
336 * and for exiting the initial stage.
338 * \param[in] params The bias parameters.
340 * \param[in] detectedCovering True if we detected that the sampling interval has been sufficiently covered.
341 * \param[in,out] fplog Log file.
342 * \returns the new histogram size.
344 double newHistogramSizeInitialStage(const BiasParams& params, double t, bool detectedCovering, FILE* fplog);
347 * Check if the sampling region has been covered "enough" or not.
349 * A one-dimensional interval is defined as covered if each point has
350 * accumulated the same weight as is in the peak of a discretized normal
351 * distribution. For multiple dimensions, the weights are simply projected
352 * onto each dimension and the multidimensional space is covered if each
355 * \note The covering criterion for multiple dimensions could improved, e.g.
356 * by using a path finding algorithm.
358 * \param[in] params The bias parameters.
359 * \param[in] dimParams Bias dimension parameters.
360 * \param[in] grid The grid.
361 * \returns true if covered.
363 bool isSamplingRegionCovered(const BiasParams& params,
364 ArrayRef<const DimParams> dimParams,
365 const BiasGrid& grid) const;
368 * Return the new reference weight histogram size for the current update.
370 * This function also takes care of checking for covering in the initial stage.
372 * \param[in] params The bias parameters.
374 * \param[in] covered True if the sampling interval has been covered enough.
375 * \param[in,out] fplog Log file.
376 * \returns the new histogram size.
378 double newHistogramSize(const BiasParams& params, double t, bool covered, FILE* fplog);
382 * Update the reaction coordinate value.
384 * \param[in] grid The bias grid.
385 * \param[in] coordValue The current reaction coordinate value (there are no limits on allowed values).
387 void setCoordValue(const BiasGrid& grid, const awh_dvec coordValue)
389 coordState_.setCoordValue(grid, coordValue);
393 * Performs an update of the bias.
395 * The objective of the update is to use collected samples (probability weights)
396 * to improve the free energy estimate. For sake of efficiency, the update is
397 * local whenever possible, meaning that only points that have actually been sampled
398 * are accessed and updated here. For certain AWH settings or at certain steps
399 * however, global need to be performed. Besides the actual free energy update, this
400 * function takes care of ensuring future convergence of the free energy. Convergence
401 * is obtained by increasing the size of the reference weight histogram in a controlled
402 * (sometimes dynamic) manner. Also, there are AWH variables that are direct functions
403 * of the free energy or sampling history that need to be updated here, namely the target
404 * distribution and the bias function.
406 * \param[in] dimParams The dimension parameters.
407 * \param[in] grid The grid.
408 * \param[in] params The bias parameters.
410 * \param[in] step Time step.
411 * \param[in,out] fplog Log file.
412 * \param[in,out] updateList Work space to store a temporary list.
414 void updateFreeEnergyAndAddSamplesToHistogram(ArrayRef<const DimParams> dimParams,
415 const BiasGrid& grid,
416 const BiasParams& params,
420 std::vector<int>* updateList);
423 * Update the probability weights and the convolved bias.
425 * Given a coordinate value, each grid point is assigned a probability
426 * weight, w(point|value), that depends on the current bias function. The sum
427 * of these weights is needed for normalizing the probability sum to 1 but
428 * also equals the effective, or convolved, biasing weight for this coordinate
429 * value. The convolved bias is needed e.g. for extracting the PMF, so we save
430 * it here since this saves us from doing extra exponential function evaluations
433 * \param[in] dimParams The bias dimensions parameters
434 * \param[in] grid The grid.
435 * \param[in] neighborLambdaEnergies An array containing the energy of the system
436 * in neighboring lambdas. The array is of length numLambdas+1, where numLambdas is
437 * the number of free energy lambda states. Element 0 in the array is the energy
438 * of the current state and elements 1..numLambdas contain the energy of the system in the
439 * neighboring lambda states (also including the current state). When there are no free
440 * energy lambda state dimensions this can be empty.
441 * \param[out] weight Probability weights of the neighbors, SIMD aligned.
442 * \returns the convolved bias.
445 double updateProbabilityWeightsAndConvolvedBias(ArrayRef<const DimParams> dimParams,
446 const BiasGrid& grid,
447 ArrayRef<const double> neighborLambdaEnergies,
448 std::vector<double, AlignedAllocator<double>>* weight) const;
451 * Take samples of the current probability weights for future updates and analysis.
453 * Points in the current neighborhood will now have data meaning they
454 * need to be included in the local update list of the next update.
455 * Therefore, the local update range is also update here.
457 * \param[in] grid The grid.
458 * \param[in] probWeightNeighbor Probability weights of the neighbors.
460 void sampleProbabilityWeights(const BiasGrid& grid, ArrayRef<const double> probWeightNeighbor);
463 * Sample the reaction coordinate and PMF for future updates or analysis.
465 * These samples do not affect the (future) sampling and are thus
466 * pure observables. Statisics of these are stored in the energy file.
468 * \param[in] dimParams The bias dimensions parameters
469 * \param[in] grid The grid.
470 * \param[in] probWeightNeighbor Probability weights of the neighbors.
471 * \param[in] convolvedBias The convolved bias.
473 void sampleCoordAndPmf(const std::vector<DimParams>& dimParams,
474 const BiasGrid& grid,
475 ArrayRef<const double> probWeightNeighbor,
476 double convolvedBias);
478 * Calculates the convolved bias for a given coordinate value.
480 * The convolved bias is the effective bias acting on the coordinate.
481 * Since the bias here has arbitrary normalization, this only makes
482 * sense as a relative, to other coordinate values, measure of the bias.
484 * \note If it turns out to be costly to calculate this pointwise
485 * the convolved bias for the whole grid could be returned instead.
487 * \param[in] dimParams The bias dimensions parameters
488 * \param[in] grid The grid.
489 * \param[in] coordValue Coordinate value.
490 * \returns the convolved bias >= -GMX_FLOAT_MAX.
492 double calcConvolvedBias(ArrayRef<const DimParams> dimParams,
493 const BiasGrid& grid,
494 const awh_dvec& coordValue) const;
497 * Fills the given array with PMF values.
499 * Points outside of the biasing target region will get PMF = GMX_FLOAT_MAX.
500 * \note: The PMF is in single precision, because it is a statistical
501 * quantity and therefore never reaches full float precision.
503 * \param[out] pmf Array(ref) to be filled with the PMF values, should have the same size as the bias grid.
505 void getPmf(ArrayRef<float> /*pmf*/) const;
507 /*! \brief Returns the current coordinate state.
509 const CoordState& coordState() const { return coordState_; }
511 /*! \brief Returns a const reference to the point state.
513 const std::vector<PointState>& points() const { return points_; }
515 /*! \brief Returns true if we are in the initial stage.
517 bool inInitialStage() const { return histogramSize_.inInitialStage(); }
519 /*! \brief Returns the current histogram size.
521 inline HistogramSize histogramSize() const { return histogramSize_; }
523 /*! \brief Sets the umbrella grid point to the current grid point
525 void setUmbrellaGridpointToGridpoint() { coordState_.setUmbrellaGridpointToGridpoint(); }
529 CoordState coordState_; /**< The Current coordinate state */
531 /* The grid point state */
532 std::vector<PointState> points_; /**< Vector of state of the grid points */
534 /* Covering values for each point on the grid */
535 std::vector<double> weightSumCovering_; /**< Accumulated weights for covering checks */
537 HistogramSize histogramSize_; /**< Global histogram size related values. */
539 /* Track the part of the grid sampled since the last update. */
540 awh_ivec originUpdatelist_; /**< The origin of the rectangular region that has been sampled since last update. */
541 awh_ivec endUpdatelist_; /**< The end of the rectangular region that has been sampled since last update. */
543 //! Object for sharing biases over multiple simulations, can be nullptr
544 const BiasSharing* biasSharing_;
547 //! Linewidth used for warning output
548 static const int c_linewidth = 80 - 2;
550 //! Indent used for warning output
551 static const int c_indent = 0;
555 #endif /* GMX_AWH_BIASSTATE_H */