Gromacs
2018.8
|
#include <gromacs/awh/biasstate.h>
The state of a bias.
The bias state has the current coordinate state: its value and the grid point it maps to (the grid point of the umbrella potential if needed). It contains a vector with the state for each point on the grid. It also counts the number of updates issued and tracks which points have been sampled since last update. Finally, the convergence state is a global property set ultimately by the histogram size histogramSize in the sub-class HistogramSize, since the update sizes are ~ 1/histogramSize.
Public Member Functions | |
BiasState (const AwhBiasParams &awhBiasParams, double histogramSizeInitial, const std::vector< DimParams > &dimParams, const Grid &grid) | |
Constructor. More... | |
void | restoreFromHistory (const AwhBiasHistory &biasHistory, const Grid &grid) |
Restore the bias state from history. More... | |
void | broadcast (const t_commrec *commRecord) |
Broadcast the bias state over the MPI ranks in this simulation. More... | |
void | initHistoryFromState (AwhBiasHistory *biasHistory) const |
Allocate and initialize a bias history with the given bias state. More... | |
void | updateHistory (AwhBiasHistory *biasHistory, const Grid &grid) const |
Update the bias state history with the current state. More... | |
void | initGridPointState (const AwhBiasParams &awhBiasParams, const std::vector< DimParams > &dimParams, const Grid &grid, const BiasParams ¶ms, const std::string &filename, int numBias) |
Initialize the state of grid coordinate points. More... | |
int | warnForHistogramAnomalies (const Grid &grid, int biasIndex, double t, FILE *fplog, int maxNumWarnings) const |
Performs statistical checks on the collected histograms and warns if issues are detected. More... | |
double | calcUmbrellaForceAndPotential (const std::vector< DimParams > &dimParams, const Grid &grid, int point, gmx::ArrayRef< double > force) const |
Calculates and sets the force the coordinate experiences from an umbrella centered at the given point. More... | |
void | calcConvolvedForce (const std::vector< DimParams > &dimParams, const Grid &grid, gmx::ArrayRef< const double > probWeightNeighbor, gmx::ArrayRef< double > forceWorkBuffer, gmx::ArrayRef< double > force) const |
Calculates and sets the convolved force acting on the coordinate. More... | |
double | moveUmbrella (const std::vector< DimParams > &dimParams, const Grid &grid, gmx::ArrayRef< const double > probWeightNeighbor, gmx::ArrayRef< double > biasForce, gmx_int64_t step, gmx_int64_t seed, int indexSeed) |
Move the center point of the umbrella potential. More... | |
void | doSkippedUpdatesForAllPoints (const BiasParams ¶ms) |
Do all previously skipped updates. Public for use by tests. More... | |
void | doSkippedUpdatesInNeighborhood (const BiasParams ¶ms, const Grid &grid) |
Do previously skipped updates in this neighborhood. More... | |
void | setCoordValue (const Grid &grid, const awh_dvec coordValue) |
Update the reaction coordinate value. More... | |
void | updateFreeEnergyAndAddSamplesToHistogram (const std::vector< DimParams > &dimParams, const Grid &grid, const BiasParams ¶ms, const t_commrec *commRecord, const gmx_multisim_t *ms, double t, gmx_int64_t step, FILE *fplog, std::vector< int > *updateList) |
Performs an update of the bias. More... | |
double | updateProbabilityWeightsAndConvolvedBias (const std::vector< DimParams > &dimParams, const Grid &grid, std::vector< double, AlignedAllocator< double >> *weight) const |
Update the probability weights and the convolved bias. More... | |
void | sampleProbabilityWeights (const Grid &grid, gmx::ArrayRef< const double > probWeightNeighbor) |
Take samples of the current probability weights for future updates and analysis. More... | |
void | sampleCoordAndPmf (const Grid &grid, gmx::ArrayRef< const double > probWeightNeighbor, double convolvedBias) |
Sample the reaction coordinate and PMF for future updates or analysis. More... | |
double | calcConvolvedBias (const std::vector< DimParams > &dimParams, const Grid &grid, const awh_dvec &coordValue) const |
Calculates the convolved bias for a given coordinate value. More... | |
void | getPmf (gmx::ArrayRef< float >) const |
Fills the given array with PMF values. More... | |
const CoordState & | coordState () const |
Returns the current coordinate state. | |
const std::vector< PointState > & | points () const |
Returns a const reference to the point state. | |
bool | inInitialStage () const |
Returns true if we are in the initial stage. | |
HistogramSize | histogramSize () const |
Returns the current histogram size. | |
gmx::BiasState::BiasState | ( | const AwhBiasParams & | awhBiasParams, |
double | histogramSizeInitial, | ||
const std::vector< DimParams > & | dimParams, | ||
const Grid & | grid | ||
) |
Constructor.
Constructs the global state and the point states on a provided geometric grid passed in grid
.
[in] | awhBiasParams | The Bias parameters from inputrec. |
[in] | histogramSizeInitial | The estimated initial histogram size. This is floating-point, since histograms use weighted entries and grow by a floating-point scaling factor. |
[in] | dimParams | The dimension parameters. |
[in] | grid | The bias grid. |
void gmx::BiasState::broadcast | ( | const t_commrec * | commRecord | ) |
Broadcast the bias state over the MPI ranks in this simulation.
[in] | commRecord | Struct for communication. |
double gmx::BiasState::calcConvolvedBias | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
const awh_dvec & | coordValue | ||
) | const |
Calculates the convolved bias for a given coordinate value.
The convolved bias is the effective bias acting on the coordinate. Since the bias here has arbitrary normalization, this only makes sense as a relative, to other coordinate values, measure of the bias.
[in] | dimParams | The bias dimensions parameters |
[in] | grid | The grid. |
[in] | coordValue | Coordinate value. |
void gmx::BiasState::calcConvolvedForce | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
gmx::ArrayRef< const double > | probWeightNeighbor, | ||
gmx::ArrayRef< double > | forceWorkBuffer, | ||
gmx::ArrayRef< double > | force | ||
) | const |
Calculates and sets the convolved force acting on the coordinate.
The convolved force is the weighted sum of forces from umbrellas located at each point in the grid.
[in] | dimParams | The bias dimensions parameters. |
[in] | grid | The grid. |
[in] | probWeightNeighbor | Probability weights of the neighbors. |
[in] | forceWorkBuffer | Force work buffer, values only used internally. |
[in,out] | force | Bias force vector to set. |
double gmx::BiasState::calcUmbrellaForceAndPotential | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
int | point, | ||
gmx::ArrayRef< double > | force | ||
) | const |
Calculates and sets the force the coordinate experiences from an umbrella centered at the given point.
The umbrella potential is an harmonic potential given by 0.5k(coord value - point value)^2. This value is also returned.
[in] | dimParams | The bias dimensions parameters. |
[in] | grid | The grid. |
[in] | point | Point for umbrella center. |
[in,out] | force | Force vector to set. Returns the umbrella potential. |
void gmx::BiasState::doSkippedUpdatesForAllPoints | ( | const BiasParams & | params | ) |
Do all previously skipped updates. Public for use by tests.
[in] | params | The bias parameters. |
void gmx::BiasState::doSkippedUpdatesInNeighborhood | ( | const BiasParams & | params, |
const Grid & | grid | ||
) |
Do previously skipped updates in this neighborhood.
[in] | params | The bias parameters. |
[in] | grid | The grid. |
void gmx::BiasState::getPmf | ( | gmx::ArrayRef< float > | pmf | ) | const |
Fills the given array with PMF values.
Points outside of the biasing target region will get PMF = GMX_FLOAT_MAX.
[out] | pmf | Array(ref) to be filled with the PMF values, should have the same size as the bias grid. |
void gmx::BiasState::initGridPointState | ( | const AwhBiasParams & | awhBiasParams, |
const std::vector< DimParams > & | dimParams, | ||
const Grid & | grid, | ||
const BiasParams & | params, | ||
const std::string & | filename, | ||
int | numBias | ||
) |
Initialize the state of grid coordinate points.
[in] | awhBiasParams | Bias parameters from inputrec. |
[in] | dimParams | The dimension parameters. |
[in] | grid | The grid. |
[in] | params | The bias parameters. |
[in] | filename | Name of file to read PMF and target from. |
[in] | numBias | The number of biases. |
void gmx::BiasState::initHistoryFromState | ( | AwhBiasHistory * | biasHistory | ) | const |
Allocate and initialize a bias history with the given bias state.
This function will be called at the start of a new simulation. Note that this only sets the correct size and does produce a valid history object, but with all data set to zero. Actual history data is set by updateHistory.
[in,out] | biasHistory | AWH history to initialize. |
double gmx::BiasState::moveUmbrella | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
gmx::ArrayRef< const double > | probWeightNeighbor, | ||
gmx::ArrayRef< double > | biasForce, | ||
gmx_int64_t | step, | ||
gmx_int64_t | seed, | ||
int | indexSeed | ||
) |
Move the center point of the umbrella potential.
A new umbrella center is sampled from the biased distibution. Also, the bias force is updated and the new potential is return.
This function should only be called when the bias force is not being convolved. It is assumed that the probability distribution has been updated.
[in] | dimParams | Bias dimension parameters. |
[in] | grid | The grid. |
[in] | probWeightNeighbor | Probability weights of the neighbors. |
[in,out] | biasForce | The AWH bias force. |
[in] | step | Step number, needed for the random number generator. |
[in] | seed | Random seed. |
[in] | indexSeed | Second random seed, should be the bias Index. |
void gmx::BiasState::restoreFromHistory | ( | const AwhBiasHistory & | biasHistory, |
const Grid & | grid | ||
) |
Restore the bias state from history.
[in] | biasHistory | Bias history struct. |
[in] | grid | The bias grid. |
void gmx::BiasState::sampleCoordAndPmf | ( | const Grid & | grid, |
gmx::ArrayRef< const double > | probWeightNeighbor, | ||
double | convolvedBias | ||
) |
Sample the reaction coordinate and PMF for future updates or analysis.
These samples do not affect the (future) sampling and are thus pure observables. Statisics of these are stored in the energy file.
[in] | grid | The grid. |
[in] | probWeightNeighbor | Probability weights of the neighbors. |
[in] | convolvedBias | The convolved bias. |
void gmx::BiasState::sampleProbabilityWeights | ( | const Grid & | grid, |
gmx::ArrayRef< const double > | probWeightNeighbor | ||
) |
Take samples of the current probability weights for future updates and analysis.
Points in the current neighborhood will now have data meaning they need to be included in the local update list of the next update. Therefore, the local update range is also update here.
[in] | grid | The grid. |
[in] | probWeightNeighbor | Probability weights of the neighbors. |
Update the reaction coordinate value.
[in] | grid | The bias grid. |
[in] | coordValue | The current reaction coordinate value (there are no limits on allowed values). |
void gmx::BiasState::updateFreeEnergyAndAddSamplesToHistogram | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
const BiasParams & | params, | ||
const t_commrec * | commRecord, | ||
const gmx_multisim_t * | ms, | ||
double | t, | ||
gmx_int64_t | step, | ||
FILE * | fplog, | ||
std::vector< int > * | updateList | ||
) |
Performs an update of the bias.
The objective of the update is to use collected samples (probability weights) to improve the free energy estimate. For sake of efficiency, the update is local whenever possible, meaning that only points that have actually been sampled are accessed and updated here. For certain AWH settings or at certain steps however, global need to be performed. Besides the actual free energy update, this function takes care of ensuring future convergence of the free energy. Convergence is obtained by increasing the size of the reference weight histogram in a controlled (sometimes dynamic) manner. Also, there are AWH variables that are direct functions of the free energy or sampling history that need to be updated here, namely the target distribution and the bias function.
[in] | dimParams | The dimension parameters. |
[in] | grid | The grid. |
[in] | params | The bias parameters. |
[in] | commRecord | Struct for intra-simulation communication. |
[in] | ms | Struct for multi-simulation communication. |
[in] | t | Time. |
[in] | step | Time step. |
[in,out] | fplog | Log file. |
[in,out] | updateList | Work space to store a temporary list. |
void gmx::BiasState::updateHistory | ( | AwhBiasHistory * | biasHistory, |
const Grid & | grid | ||
) | const |
Update the bias state history with the current state.
[out] | biasHistory | Bias history struct. |
[in] | grid | The bias grid. |
double gmx::BiasState::updateProbabilityWeightsAndConvolvedBias | ( | const std::vector< DimParams > & | dimParams, |
const Grid & | grid, | ||
std::vector< double, AlignedAllocator< double >> * | weight | ||
) | const |
Update the probability weights and the convolved bias.
Given a coordinate value, each grid point is assigned a probability weight, w(point|value), that depends on the current bias function. The sum of these weights is needed for normalizing the probability sum to 1 but also equals the effective, or convolved, biasing weight for this coordinate value. The convolved bias is needed e.g. for extracting the PMF, so we save it here since this saves us from doing extra exponential function evaluations later on.
[in] | dimParams | The bias dimensions parameters |
[in] | grid | The grid. |
[out] | weight | Probability weights of the neighbors, SIMD aligned. |
int gmx::BiasState::warnForHistogramAnomalies | ( | const Grid & | grid, |
int | biasIndex, | ||
double | t, | ||
FILE * | fplog, | ||
int | maxNumWarnings | ||
) | const |
Performs statistical checks on the collected histograms and warns if issues are detected.
[in] | grid | The grid. |
[in] | biasIndex | The index of the bias we are checking for. |
[in] | t | Time. |
[in,out] | fplog | Output file for warnings. |
[in] | maxNumWarnings | Don't issue more than this number of warnings. |