Shortcut: WHAM_HISTOGRAM

Module wham
Description Usage
This can be used to output the a histogram using the weighted histogram technique used in 0 tutorialsused in 0 eggs
output value type
the histogram that was generated using the WHAM weights grid

Further details and examples

This can be used to output the a histogram using the weighted histogram technique

This shortcut action allows you to calculate a histogram using the weighted histogram analysis technique. The following input illustrates how this is used in practise to analyze the output from a series of umbrella sampling calculations. The trajectory from each of the simulations run with the different biases should be concatenated into a single trajectory before running the following analysis script on the concatenated trajectory using PLUMED driver. The umbrella sampling simulations that will be analyzed using the script below applied a harmonic restraint that restrained the torsional angle involving atoms 5, 7, 9 and 15 to particular values. The script below calculates the reweighting weights for each of the trajectories and then applies the binless WHAM algorithm to determine a weight for each configuration in the concatenated trajectory. A histogram is then constructed from the configurations visited and their weights. This histogram is then converted into a free energy surface and output to a file called fes.dat

Click on the labels of the actions for more information on what each action computes
tested on2.11
#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. This action has hidden defaults. 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}
hh: WHAM_HISTOGRAMThis can be used to output the a histogram using the weighted histogram technique This action is a shortcut and it has hidden defaults. More details ARGthe arguments that you would like to make the histogram for=phi BIAS the value of the biases to use when performing WHAM=rp.bias TEMPthe temperature at which the simulation was run=300 GRID_MINthe minimum to use for the grid=-pi GRID_MAXthe maximum to use for the grid=pi GRID_BINthe number of bins to use for the grid=50
fes: CONVERT_TO_FESConvert a histogram to a free energy surface. This action is a shortcut. More details ARGthe histogram that you would like to convert into a free energy surface=hh TEMPthe temperature at which you are operating=300
DUMPGRIDOutput the function on the grid to a file with the PLUMED grid format. More details ARGthe label for the grid that you would like to output=fes FILE the file on which to write the grid=fes.dat

The script above must be run with multiple replicas using the following command:

mpirun -np 4 plumed driver --mf_xtc alltraj.xtc --multi 4

For more details on how the weights for configurations are determined using the wham method see the documentation for the WHAM action.

Syntax

The following table describes the keywords and options that can be used with this action

Keyword Type Default Description
ARG compulsory none the arguments that you would like to make the histogram for
BIAS compulsory *.bias the value of the biases to use when performing WHAM
TEMP compulsory none the temperature at which the simulation was run
STRIDE compulsory 1 the frequency with which the data should be stored to perform WHAM
GRID_MIN compulsory none the minimum to use for the grid
GRID_MAX compulsory none the maximum to use for the grid
GRID_BIN compulsory none the number of bins to use for the grid
BANDWIDTH optional not used the bandwidth for kernel density estimation