PBMETAD
This is part of the bias module

Used to performed Parallel Bias metadynamics.

This action activate Parallel Bias Metadynamics (PBMetaD) [79], a version of metadynamics [60] in which multiple low-dimensional bias potentials are applied in parallel. In the current implementation, these have the form of mono-dimensional metadynamics bias potentials:

\[ {V(s_1,t), ..., V(s_N,t)} \]

where:

\[ V(s_i,t) = \sum_{ k \tau < t} W_i(k \tau) \exp\left( - \frac{(s_i-s_i^{(0)}(k \tau))^2}{2\sigma_i^2} \right). \]

To ensure the convergence of each mono-dimensional bias potential to the corresponding free energy, at each deposition step the Gaussian heights are multiplied by the so-called conditional term:

\[ W_i(k \tau)=W_0 \frac{\exp\left( - \frac{V(s_i,k \tau)}{k_B T} \right)}{\sum_{i=1}^N \exp\left( - \frac{V(s_i,k \tau)}{k_B T} \right)} \]

where \(W_0\) is the initial Gaussian height.

The PBMetaD bias potential is defined by:

\[ V_{PB}(\vec{s},t) = -k_B T \log{\sum_{i=1}^N \exp\left( - \frac{V(s_i,t)}{k_B T} \right)}. \]

Information on the Gaussian functions that build each bias potential are printed to multiple HILLS files, which are used both to restart the calculation and to reconstruct the mono-dimensional free energies as a function of the corresponding CVs. These can be reconstructed using the sum_hills utility because the final bias is given by:

\[ V(s_i) = -F(s_i) \]

Currently, only a subset of the METAD options are available in PBMetaD.

The bias potentials can be stored on a grid to increase performances of long PBMetaD simulations. You should provide either the number of bins for every collective variable (GRID_BIN) or the desired grid spacing (GRID_SPACING). In case you provide both PLUMED will use the most conservative choice (highest number of bins) for each dimension. In case you do not provide any information about bin size (neither GRID_BIN nor GRID_SPACING) and if Gaussian width is fixed PLUMED will use 1/5 of the Gaussian width as grid spacing. This default choice should be reasonable for most applications.

Another option that is available is well-tempered metadynamics [7]. In this variant of PBMetaD the heights of the Gaussian hills are scaled at each step by the additional well-tempered metadynamics term. This ensures that each bias converges more smoothly. It should be noted that, in the case of well-tempered metadynamics, in the output printed the Gaussian height is re-scaled using the bias factor. Also notice that with well-tempered metadynamics the HILLS files do not contain the bias, but the negative of the free-energy estimate. This choice has the advantage that one can restart a simulation using a different value for the \(\Delta T\). The applied bias will be scaled accordingly.

Note that you can use here also the flexible gaussian approach [22] in which you can adapt the gaussian to the extent of Cartesian space covered by a variable or to the space in collective variable covered in a given time. In this case the width of the deposited gaussian potential is denoted by one value only that is a Cartesian space (ADAPTIVE=GEOM) or a time (ADAPTIVE=DIFF). Note that a specific integration technique for the deposited Gaussian kernels should be used in this case. Check the documentation for utility sum_hills.

With the keyword INTERVAL one changes the metadynamics algorithm setting the bias force equal to zero outside boundary [6]. If, for example, metadynamics is performed on a CV s and one is interested only to the free energy for s > boundary, the history dependent potential is still updated according to the above equations but the metadynamics force is set to zero for s < boundary. Notice that Gaussians are added also if s < boundary, as the tails of these Gaussians influence VG in the relevant region s > boundary. In this way, the force on the system in the region s > boundary comes from both metadynamics and the force field, in the region s < boundary only from the latter. This approach allows obtaining a history-dependent bias potential VG that fluctuates around a stable estimator, equal to the negative of the free energy far enough from the boundaries. Note that:

  • It works only for one-dimensional biases;
  • It works both with and without GRID;
  • The interval limit boundary in a region where the free energy derivative is not large;
  • If in the region outside the limit boundary the system has a free energy minimum, the INTERVAL keyword should be used together with a UPPER_WALLS or LOWER_WALLS at boundary.

Multiple walkers [85] can also be used. See below the examples.

Description of components

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
Compulsory keywords
SIGMA the widths of the Gaussian hills
PACE the frequency for hill addition, one for all biases
Options
NUMERICAL_DERIVATIVES ( default=off ) calculate the derivatives for these quantities numerically
GRID_SPARSE ( default=off ) use a sparse grid to store hills
GRID_NOSPLINE ( default=off ) don't use spline interpolation with grids
WALKERS_MPI

( default=off ) Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR

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...
FILE files in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found
HEIGHT the height of the Gaussian hills, one for all biases. Compulsory unless TAU, TEMP and BIASFACTOR are given
FMT specify format for HILLS files (useful for decrease the number of digits in regtests)
BIASFACTOR use well tempered metadynamics with this bias factor, one for all biases. Please note you must also specify temp
TEMP the system temperature - this is only needed if you are doing well-tempered metadynamics
TAU in well tempered metadynamics, sets height to ( \(k_B \Delta T\)*pace*timestep)/tau
GRID_RFILES read grid for the bias
GRID_WSTRIDE frequency for dumping the grid
GRID_WFILES dump grid for the bias, default names are used if GRID_WSTRIDE is used without GRID_WFILES.
GRID_MIN the lower bounds for the grid
GRID_MAX the upper bounds for the grid
GRID_BIN the number of bins for the grid
GRID_SPACING the approximate grid spacing (to be used as an alternative or together with GRID_BIN)
SELECTOR add forces and do update based on the value of SELECTOR
SELECTOR_ID value of SELECTOR
WALKERS_ID walker id
WALKERS_N number of walkers
WALKERS_DIR shared directory with the hills files from all the walkers
WALKERS_RSTRIDE stride for reading hills files
INTERVAL_MIN one dimensional lower limits, outside the limits the system will not feel the biasing force.
INTERVAL_MAX one dimensional upper limits, outside the limits the system will not feel the biasing force.
ADAPTIVE use a geometric (=GEOM) or diffusion (=DIFF) based hills width scheme. Sigma is one number that has distance units or timestep dimensions
SIGMA_MAX the upper bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds
SIGMA_MIN the lower bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds
RESTART allows per-action setting of restart (YES/NO/AUTO)
UPDATE_FROM Only update this action from this time
UPDATE_UNTIL Only update this action until this time
Examples

The following input is for PBMetaD calculation using as collective variables the distance between atoms 3 and 5 and the distance between atoms 2 and 4. The value of the CVs and the PBMetaD bias potential are written to the COLVAR file every 100 steps.

DISTANCE ATOMS=3,5 LABEL=d1
DISTANCE ATOMS=2,4 LABEL=d2
PBMETAD ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3 PACE=500 LABEL=pb FILE=HILLS_d1,HILLS_d2
PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR

(See also DISTANCE and PRINT).

If you use well-tempered metadynamics, you should specify a single bias factor and initial Gaussian height.
DISTANCE ATOMS=3,5 LABEL=d1
DISTANCE ATOMS=2,4 LABEL=d2
PBMETAD ...
ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
PACE=500 BIASFACTOR=8 LABEL=pb
FILE=HILLS_d1,HILLS_d2
... PBMETAD
PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
The following input enables the MPI version of multiple-walkers.
DISTANCE ATOMS=3,5 LABEL=d1
DISTANCE ATOMS=2,4 LABEL=d2
PBMETAD ...
ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
PACE=500 BIASFACTOR=8 LABEL=pb
FILE=HILLS_d1,HILLS_d2
WALKERS_MPI
... PBMETAD
PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
The disk version of multiple-walkers can be enabled by setting the number of walker used, the id of the current walker which interprets the input file, the directory where the hills containing files resides, and the frequency to read the other walkers. Here is an example
DISTANCE ATOMS=3,5 LABEL=d1
DISTANCE ATOMS=2,4 LABEL=d2
PBMETAD ...
ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
PACE=500 BIASFACTOR=8 LABEL=pb
FILE=HILLS_d1,HILLS_d2
WALKERS_N=10
WALKERS_ID=3
WALKERS_DIR=../
WALKERS_RSTRIDE=100
... PBMETAD
PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
where WALKERS_N is the total number of walkers, WALKERS_ID is the id of the present walker (starting from 0 ) and the WALKERS_DIR is the directory where all the walkers are located. WALKERS_RSTRIDE is the number of step between one update and the other.