Gromacs
2018.8
|
#include <gromacs/awh/bias.h>
A bias acting on a multidimensional coordinate.
At each step AWH should provide its biases with updated values of their coordinates. Each bias provides AWH with an updated bias forces and the corresponding potential.
See the user manual for details on the algorithm and equations.
The bias is responsible for keeping and updating a free energy estimate along the coordinate. The bias potential is basically a function of the free energy estimate and so also changes by the update. The free energy update is based on information from coordinate samples collected at a constant bias potential, between updates.
The bias keeps a grid with coordinate points that organizes spatial information about the coordinate. The grid has the the same geometry as the coordinate, i.e. they have the same dimensionality and periodicity (if any). The number of points in the grid sets the resolution of the collected data and its extent defines the sampling region of interest.
Each coordinate point has further statistical properties and function values which a grid point does not know about. E.g., for the bias each coordinate point is associated with values of the bias, free energy and target distribution, accumulated sampling weight, etc. For this the bias attaches to each grid point a state. The grid + vector of point states are the bias coordinate points.
The bias has a fairly complex global state keeping track of where the system (coordinate) currently is (CoordState), where it has sampled since the last update (BiasState) and controlling the free energy convergence rate (HistogramSize).
Partly, the complexity comes from the bias having two convergence stages: an initial stage which in an heuristic, non-deterministic way restricts the early convergence rate for sake of robustness; and a final stage where the convergence rate is constant. The length of the initial stage depends on the sampling and is unknown beforehand.
Another complexity comes from the fact that coordinate points, for sake of efficiency in the case of many grid points, are typically only accessed in recently sampled regions even though the free energy update is inherently global and affects all points. The bias allows points thay are non-local at the time the update was issued to postpone ("skip", as it is called in the code) the update. A non-local point is defined as a point which has not been sampled since the last update. Local points are points that have been sampled since the last update. The (current) set of local points are kept track of by the bias state and reset after every update. An update is called local if it only updates local points. Non-local points will temporarily "skip" the update until next time they are local (or when a global update is issued). For this to work, the bias keeps a global "clock" (in HistogramSize) of the number of issued updates. Each PointState also has its own local "clock" with the counting the number of updates it has pulled through. When a point updates its state it asserts that its local clock is synchronized with the global clock.
Public Types | |
enum | ThisRankWillDoIO { ThisRankWillDoIO::No, ThisRankWillDoIO::Yes } |
Enum for requesting Bias set up with(out) I/O on this rank. More... | |
Public Member Functions | |
Bias (int biasIndexInCollection, const AwhParams &awhParams, const AwhBiasParams &awhBiasParams, const std::vector< DimParams > &dimParams, double beta, double mdTimeStep, int numSharingSimulations, const std::string &biasInitFilename, ThisRankWillDoIO thisRankWillDoIO, BiasParams::DisableUpdateSkips disableUpdateSkips=BiasParams::DisableUpdateSkips::no) | |
Constructor. More... | |
void | printInitializationToLog (FILE *fplog) const |
Print information about initialization to log file. More... | |
gmx::ArrayRef< const double > | calcForceAndUpdateBias (const awh_dvec coordValue, double *awhPotential, double *potentialJump, const t_commrec *commRecord, const gmx_multisim_t *ms, double t, gmx_int64_t step, gmx_int64_t seed, FILE *fplog) |
Evolves the bias at every step. More... | |
double | calcConvolvedBias (const awh_dvec &coordValue) const |
Calculates the convolved bias for a given coordinate value. More... | |
void | restoreStateFromHistory (const AwhBiasHistory *biasHistory, const t_commrec *cr) |
Restore the bias state from history on the master rank and broadcast it. More... | |
void | initHistoryFromState (AwhBiasHistory *biasHistory) const |
Allocate and initialize a bias history with the given bias state. More... | |
void | updateHistory (AwhBiasHistory *biasHistory) const |
Update the bias history with the current state. More... | |
void | doSkippedUpdatesForAllPoints () |
Do all previously skipped updates. Public for use by tests. | |
int | ndim () const |
Returns the dimensionality of the bias. | |
const std::vector< DimParams > & | dimParams () const |
Returns the dimension parameters. | |
const BiasParams & | params () const |
Returns the bias parameters. | |
const BiasState & | state () const |
Returns the global state of the bias. | |
int | biasIndex () const |
Returns the index of the bias. | |
const awh_dvec & | getGridCoordValue (size_t gridPointIndex) const |
Return the coordinate value for a grid point. More... | |
const CorrelationGrid & | forceCorrelationGrid () const |
Return a const reference to the force correlation grid. | |
int | numEnergySubblocksToWrite () const |
Return the number of data blocks that have been prepared for writing. | |
int | writeToEnergySubblocks (t_enxsubblock *subblock) const |
Write bias data blocks to energy subblocks. More... | |
|
strong |
Enum for requesting Bias set up with(out) I/O on this rank.
Enumerator | |
---|---|
No |
This rank will not do I/O. |
Yes |
This rank will do I/O. |
gmx::Bias::Bias | ( | int | biasIndexInCollection, |
const AwhParams & | awhParams, | ||
const AwhBiasParams & | awhBiasParams, | ||
const std::vector< DimParams > & | dimParams, | ||
double | beta, | ||
double | mdTimeStep, | ||
int | numSharingSimulations, | ||
const std::string & | biasInitFilename, | ||
ThisRankWillDoIO | thisRankWillDoIO, | ||
BiasParams::DisableUpdateSkips | disableUpdateSkips = BiasParams::DisableUpdateSkips::no |
||
) |
Constructor.
[in] | biasIndexInCollection | Index of the bias in collection. |
[in] | awhParams | AWH parameters. |
[in] | awhBiasParams | Bias parameters. |
[in] | dimParams | Bias dimension parameters. |
[in] | beta | 1/(k_B T). |
[in] | mdTimeStep | The MD time step. |
[in] | numSharingSimulations | The number of simulations to share the bias across. |
[in] | biasInitFilename | Name of file to read PMF and target from. |
[in] | thisRankWillDoIO | Tells whether this MPI rank will do I/O (checkpointing, AWH output), normally (only) the master rank does I/O. |
[in] | disableUpdateSkips | If to disable update skips, useful for testing. |
|
inline |
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] | coordValue | The coordinate value. |
gmx::ArrayRef< const double > gmx::Bias::calcForceAndUpdateBias | ( | const awh_dvec | coordValue, |
double * | awhPotential, | ||
double * | potentialJump, | ||
const t_commrec * | commRecord, | ||
const gmx_multisim_t * | ms, | ||
double | t, | ||
gmx_int64_t | step, | ||
gmx_int64_t | seed, | ||
FILE * | fplog | ||
) |
Evolves the bias at every step.
At each step the bias step needs to:
[in] | coordValue | The current coordinate value(s). |
[out] | awhPotential | Bias potential. |
[out] | potentialJump | Change in bias potential for this bias. |
[in] | commRecord | Struct for intra-simulation communication. |
[in] | ms | Struct for multi-simulation communication. |
[in] | t | Time. |
[in] | step | Time step. |
[in] | seed | Random seed. |
[in,out] | fplog | Log file. |
|
inline |
Return the coordinate value for a grid point.
[in] | gridPointIndex | The index of the grid point. |
void gmx::Bias::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 only constant data will be initialized here. History data is set by updateHistory.
[in,out] | biasHistory | AWH history to initialize. |
void gmx::Bias::printInitializationToLog | ( | FILE * | fplog | ) | const |
Print information about initialization to log file.
Prints information about AWH variables that are set internally but might be of interest to the user.
[in,out] | fplog | Log file, can be nullptr. |
void gmx::Bias::restoreStateFromHistory | ( | const AwhBiasHistory * | biasHistory, |
const t_commrec * | cr | ||
) |
Restore the bias state from history on the master rank and broadcast it.
[in] | biasHistory | Bias history struct, only allowed to be nullptr on non-master ranks. |
[in] | cr | The communication record. |
void gmx::Bias::updateHistory | ( | AwhBiasHistory * | biasHistory | ) | const |
Update the bias history with the current state.
[out] | biasHistory | Bias history struct. |
int gmx::Bias::writeToEnergySubblocks | ( | t_enxsubblock * | subblock | ) | const |
Write bias data blocks to energy subblocks.
[in,out] | subblock | Energy subblocks to write to. |