This is part of the ves module | |
It is only available if you configure PLUMED with ./configure –enable-modules=ves . Furthermore, this feature is still being developed so take care when using it and report any problems on the mailing list. |
Linear basis set expansion bias.
This VES bias action takes the bias potential to be a linear expansion in some basis set that is written as a product of one-dimensional basis functions. For example, for one CV the bias would be written as
\[ V(s_{1};\boldsymbol{\alpha}) = \sum_{i_{1}} \alpha_{i_{1}} \, f_{i_{1}}(s_{1}), \]
while for two CVs it is written as
\[ V(s_{1},s_{2};\boldsymbol{\alpha}) = \sum_{i_{1},i_{2}} \alpha_{i_{1},i_{2}} \, f_{i_{1}}(s_{1}) \, f_{i_{2}}(s_{2}) \]
where \(\boldsymbol{\alpha}\) is the set of expansion coefficients that are optimized within VES. With an appropriate choice of the basis functions it is possible to represent any generic free energy surface. The relationship between the bias and the free energy surface is given by
\[ V(\mathbf{s}) = - F(\mathbf{s}) - \frac{1}{\beta} \log p(\mathbf{s}). \]
where \(p(\mathbf{s})\) is the target distribution that is employed in the VES simulation.
Various one-dimensional basis functions are available in the VES code, see the complete list here. At the current moment we recommend to use Legendre polynomials (BF_LEGENDRE) for non-periodic CVs and Fourier basis functions (BF_FOURIER) for periodic CV (e.g. dihedral angles).
To use basis functions within VES_LINEAR_EXPANSION you first need to define them in the input file before the VES_LINEAR_EXPANSION action and then give their labels using the BASIS_FUNCTIONS keyword.
Various target distributions \(p(\mathbf{s})\) are available in the VES code, see the complete list here.
To use a target distribution within VES_LINEAR_EXPANSION you first need to define it in the input file before the VES_LINEAR_EXPANSION action and then give its label using the TARGET_DISTRIBUTION keyword. The default behavior if no TARGET_DISTRIBUTION is given is to employ a uniform target distribution.
Some target distribution, like the well-tempered one (TD_WELLTEMPERED), are dynamic and need to be iteratively updated during the optimization.
In order to optimize the coefficients you will need to use VES_LINEAR_EXPANSION in combination with an optimizer, see the list of optimizers available in the VES code here. At the current moment we recommend to use the averaged stochastic gradient decent optimizer (OPT_AVERAGED_SGD).
The optimizer should be defined after the VES_LINEAR_EXPANSION action.
Internally the code uses grids to calculate the basis set averages over the target distribution that is needed for the gradient. The same grid is also used for the output files (see next section). The size of the grid is determined by the GRID_BINS keyword. By default it has 100 grid points in each dimension, and generally this value should be sufficient.
It is possible to output on-the-fly during the simulation the free energy surface estimated from the bias potential. How often this is done is specified within the optimizer by using the FES_OUTPUT keyword. The filename is specified by the FES_FILE keyword, but by default is it fes.LABEL.data, with an added suffix indicating the iteration number (iter-#).
For multi-dimensional case is it possible to also output projections of the free energy surfaces. The arguments for which to do these projections is specified using the numbered PROJ_ARG keywords. For these files a suffix indicating the projection (proj-#) will be added to the filenames. You will also need to specify the frequency of the output by using the FES_PROJ_OUTPUT keyword within the optimizer.
It is also possible to output the bias potential itself, for this the relevant keyword is BIAS_OUTPUT within the optimizer. The filename is specified by the BIAS_FILE keyword, but by default is it bias.LABEL.data, with an added suffix indicating the iteration number (iter-#).
Furthermore is it possible to output the target distribution, and its projections (i.e. marginal distributions). The filenames of these files are specified with the TARGETDIST_FILE, but by default is it targetdist.LABEL.data. The logarithm of the target distribution will also be outputted to file that has the added suffix log. For static target distribution these files will be outputted in the beginning of the simulation while for dynamic ones you will need to specify the frequency of the output by using the TARGETDIST_OUTPUT and TARGETDIST_PROJ_OUTPUT keywords within the optimizer.
It is also possible to output free energy surfaces and bias in post processing by using the VES_OUTPUT_FES action. However, be aware that this action does does not support dynamic target distribution (e.g. well-tempered).
It is also possible to use VES_LINEAR_EXPANSION as a static bias that uses previously obtained coefficients. In this case the coefficients should be read in from the coefficient file given in the COEFFS keyword.
It is possible to impose a cutoff on the bias potential using the procedure introduced in [70] such that the free energy surface is only flooded up to a certain value. The bias that results from this procedure can then be used as a static bias for obtaining kinetic rates. The value of the cutoff is given by the BIAS_CUTOFF keyword. To impose the cutoff the code uses a Fermi switching function \(1/(1+e^{\lambda x})\) where the parameter \(\lambda\) controls how sharply the switchingfunction goes to zero. The default value is \(\lambda=10\) but this can be changed by using the BIAS_CUTOFF_FERMI_LAMBDA keyword.
By default this Action calculates the following quantities. These quantities can be referenced elsewhere in the input by using this Action's label followed by a dot and the name of the quantity required from the list below.
Quantity | Description |
bias | the instantaneous value of the bias potential |
force2 | the instantaneous value of the squared force due to this bias potential. |
BASIS_FUNCTIONS | the label of the one dimensional basis functions that should be used. |
NUMERICAL_DERIVATIVES | ( default=off ) calculate the derivatives for these quantities numerically |
TEMP | the system temperature - this is needed if the MD code does not pass the temperature to PLUMED. |
BIAS_FILE | filename of the file on which the bias should be written out. By default it is bias.LABEL.data. Note that suffixes indicating the iteration number (iter-#) are added to the filename when optimizing coefficients. |
FES_FILE | filename of the file on which the FES should be written out. By default it is fes.LABEL.data. Note that suffixes indicating the iteration number (iter-#) are added to the filename when optimizing coefficients. |
TARGETDIST_FILE | filename of the file on which the target distribution should be written out. By default it is targetdist.LABEL.data. Note that suffixes indicating the iteration number (iter-#) are added to the filename when optimizing coefficients and the target distribution is dynamic. |
COEFFS | read in the coefficients from files. |
TARGET_DISTRIBUTION | the label of the target distribution to be used. |
BIAS_CUTOFF | cutoff the bias such that it only fills the free energy surface up to certain level F_cutoff, here you should give the value of the F_cutoff. |
BIAS_CUTOFF_FERMI_LAMBDA | the lambda value used in the Fermi switching function for the bias cutoff (BIAS_CUTOFF), the default value is 10.0. |
GRID_BINS | the number of bins used for the grid. The default value is 100 bins per dimension. |
PROJ_ARG | arguments for doing projections of the FES or the target distribution. You can use multiple instances of this keyword i.e. PROJ_ARG1, PROJ_ARG2, PROJ_ARG3... |
ARG | the input for this action is the scalar output from one or more other actions. The particular scalars that you will use are referenced using the label of the action. If the label appears on its own then it is assumed that the Action calculates a single scalar value. The value of this scalar is thus used as the input to this new action. If * or *.* appears the scalars calculated by all the proceeding actions in the input file are taken. Some actions have multi-component outputs and each component of the output has a specific label. For example a DISTANCE action labelled dist may have three components x, y and z. To take just the x component you should use dist.x, if you wish to take all three components then use dist.*.More information on the referencing of Actions can be found in the section of the manual on the PLUMED Getting Started. Scalar values can also be referenced using POSIX regular expressions as detailed in the section on Regular Expressions. To use this feature you you must compile PLUMED with the appropriate flag. You can use multiple instances of this keyword i.e. ARG1, ARG2, ARG3... |
In the following example we run a VES_LINEAR_EXPANSION for one CV using a Legendre basis functions (BF_LEGENDRE) and a uniform target distribution as no target distribution is specified. The coefficients are optimized using averaged stochastic gradient descent optimizer (OPT_AVERAGED_SGD). Within the optimizer we specify that the FES should be outputted to file every 500 coefficients iterations (the FES_OUTPUT keyword). Parameters that are very specific to the problem at hand, like the order of the basis functions, the interval on which the basis functions are defined, and the step size used in the optimizer, are left unfilled.
bf1: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ VES_LINEAR_EXPANSION ... ARG=d1 BASIS_FUNCTIONS=bf1 TEMP=__ GRID_BINS=200 LABEL=b1 ... VES_LINEAR_EXPANSION OPT_AVERAGED_SGD ... BIAS=b1 STRIDE=1000 LABEL=o1 STEPSIZE=__ FES_OUTPUT=500 COEFFS_OUTPUT=10 ... OPT_AVERAGED_SGD
In the following example we employ VES_LINEAR_EXPANSION for two CVs, The first CV is periodic and therefore we employ a Fourier basis functions (BF_LEGENDRE) while the second CV is non-periodic so we employ a Legendre polynomials as in the previous example. For the target distribution we employ a well-tempered target distribution (TD_WELLTEMPERED), which is dynamic and needs to be iteratively updated with a stride that is given using the TARGETDIST_STRIDE within the optimizer.
bf1: BF_FOURIER ORDER=__ MINIMUM=__ MAXIMUM=__ bf2: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ td_wt: TD_WELLTEMPERED BIASFACTOR=10.0 VES_LINEAR_EXPANSION ... ARG=cv1,cv2 BASIS_FUNCTIONS=bf1,bf2 TEMP=__ GRID_BINS=100 LABEL=b1 TARGET_DISTRIBUTION=td_wt ... VES_LINEAR_EXPANSION OPT_AVERAGED_SGD ... BIAS=b1 STRIDE=1000 LABEL=o1 STEPSIZE=__ FES_OUTPUT=500 COEFFS_OUTPUT=10 TARGETDIST_STRIDE=500 ... OPT_AVERAGED_SGD
In the following example we employ a bias cutoff such that the bias only fills the free energy landscape up a certain level. In this case the target distribution is also dynamic and needs to iteratively updated.
bf1: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ bf2: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ VES_LINEAR_EXPANSION ... ARG=cv1,cv2 BASIS_FUNCTIONS=bf1,bf2 TEMP=__ GRID_BINS=100 LABEL=b1 BIAS_CUTOFF=20.0 ... VES_LINEAR_EXPANSION OPT_AVERAGED_SGD ... BIAS=b1 STRIDE=1000 LABEL=o1 STEPSIZE=__ FES_OUTPUT=500 COEFFS_OUTPUT=10 TARGETDIST_STRIDE=500 ... OPT_AVERAGED_SGD
The optimized bias potential can then be used as a static bias for obtaining kinetics. For this you need read in the final coefficients from file (e.g. coeffs_final.data in this case) by using the COEFFS keyword (also, no optimizer should be defined in the input)
bf1: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ bf2: BF_LEGENDRE ORDER=__ MINIMUM=__ MAXIMUM=__ VES_LINEAR_EXPANSION ... ARG=cv1,cv2 BASIS_FUNCTIONS=bf1,bf2 TEMP=__ GRID_BINS=100 LABEL=b1 BIAS_CUTOFF=20.0 COEFFS=coeffs_final.data ... VES_LINEAR_EXPANSION