Action: WHAM
Module | wham |
---|---|
Description | Usage |
Calculate the weights for configurations using the weighted histogram analysis method. | |
output value | type |
the vector of WHAM weights to use for reweighting the elements in a time series | vector |
Input
The arguments that serve as the input for this action are specified using one or more of the keywords in the following table.
Keyword | Type | Description |
---|---|---|
ARG | scalar/vector/matrix | the stored values for the bias |
Further details and examples
Calculate the weights for configurations using the weighted histogram analysis method.
The example input below shows how this command is used.
#SETTINGS NREPLICAS=4 phi: TORSIONCalculate one or multiple torsional angles. More details ATOMSthe four atoms involved in the torsional angle=5,7,9,15 psi: TORSIONCalculate one or multiple torsional angles. More details ATOMSthe four atoms involved in the torsional angle=7,9,15,17 rp: RESTRAINTAdds harmonic and/or linear restraints on one or more variables. More details ARGthe values the harmonic restraint acts upon=phi KAPPA specifies that the restraint is harmonic and what the values of the force constants on each of the variables are=50.0 ATthe position of the restraint=@replicas:This keyword specifies that different replicas have different values for this quantity. See here for more details.{-3.00,-1.45,0.10,1.65} # Get the bias on each of the replicas rep: GATHER_REPLICASCreate a vector that contains the copies of the input quantities from all replicas More details ARGthe argument from the various replicas that you would like to gather=rp.bias # Merge the biases on each of the replicas into a single vector all: CONCATENATEJoin vectors or matrices together More details ARGthe values that should be concatenated together to form the output vector=rep.* # Collect all the bias values col: COLLECTCollect data from the trajectory for later analysis More details TYPE required if you are collecting an object with rank>0=vector ARGthe label of the value whose time series is being stored for later analysis=all STRIDE the frequency with which the data should be collected and added to the quantity being averaged=1 wham: WHAMCalculate the weights for configurations using the weighted histogram analysis method. More details ARGthe stored values for the bias=col TEMPthe system temperature=300 DUMPVECTORPrint a vector to a file More details ARGthe labels of vectors/matrices that should be output in the file=wham FILE the file on which to write the vetors=wham_data
As illustrated in the example above this command is used when you have run N trajectories each of which was computed by integrating a different biased Hamiltonian. The input above calculates the probability of observing the set of configurations during the N trajectories that we ran using the product of multinomial distributions using the formula shown below:
P(T_)∝M∏j=1N∏k=1(ckwkjpj)tkj
In this expression the second product runs over the biases that were used when calculating the N trajectories. The first product then runs over the M bins in our histogram. The pj variable that is inside this product is the quantity we wish to extract; namely, the unbiased probability of having a set of CV values that lie within the range for the jth bin.
The quantity that we can easily extract from our simulations, tkj, however, measures the number of frames from trajectory k that are inside the jth bin. To interpret this quantity we must consider the bias that acts on each of the replicas so the wkj term is introduced. This quantity is calculated as:
wkj=e+βVk(xj)
and is essentially the factor that we have to multiply the unbiased probability of being in the bin by in order to get the probability that we would be inside this same bin in the kth of our biased simulations. Obviously, these wkj values depend on the value that the CVs take and also on the particular trajectory that we are investigating all of which, remember, have different simulation biases. Finally, ck is a free parameter that ensures that, for each k, the biased probability is normalized.
We can use the equation for the probability that was given above to find a set of values for pj that maximizes the likelihood of observing the trajectory. This constrained optimization must be performed using a set of Lagrange multipliers, λk, that ensure that each of the biased probability distributions are normalized, that is ∑jckwkjpj=1. Furthermore, as the problem is made easier if the quantity in the equation above is replaced by its logarithm we actually chose to minimize
L=M∑j=1N∑k=1tkjlnckwkjpj+∑kλk(N∑j=1ckwkjpj−1)
After some manipulations the following (WHAM) equations emerge:
pj∝∑Nk=1tkj∑kckwkjck=1∑Mj=1wkjpj
which can be solved by computing the pj values using the first of the two equations above with an initial guess for the ck values and by then refining these pj values using the ck values that are obtained by inserting the pj values obtained into the second of the two equations above.
Notice that only ∑ktkj, which is the total number of configurations from all the replicas that enter the jth bin, enters the WHAM equations above. There is thus no need to record which replica generated each of the frames. One can thus simply gather the trajectories from all the replicas together at the outset. This observation is important as it is the basis of the binless formulation of WHAM that is implemented within PLUMED.
Syntax
The following table describes the keywords and options that can be used with this action
Keyword | Type | Default | Description |
---|---|---|---|
ARG | input | none | the stored values for the bias |
MAXITER | compulsory | 1000 | maximum number of iterations for WHAM algorithm |
WHAMTOL | compulsory | 1e-10 | threshold for convergence of WHAM algorithm |
TEMP | optional | not used | the system temperature |