Action: PBMETAD

Module bias
Description Usage
Used to performed Parallel Bias metadynamics. used in 2 tutorialsused in 34 eggs

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 the labels of the scalars on which the bias will act
PF scalar specify which CVs belong in a partitioned family

Further details and examples

Used to performed Parallel Bias metadynamics.

This action provides an implementation of the Parallel Bias Metadynamics (PBMetaD) method that was introduced in the first paper cited below. This variant of metadynamics (see METAD) uses multiple low-dimensional bias potentials that are applied in parallel. In the current implementation, these have the form of mono-dimensional metadynamics bias potentials:

where:

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:

where is the initial Gaussian height.

The PBMetaD bias potential is defined by:

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:

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 the well-tempered metadynamics that is discussed in the second paper cited below. 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 . The applied bias will be scaled accordingly.

Note that you can also use here the flexible Gaussian approach that is discussed in the third paper cited below. This methods allows you to 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 and sets the bias force equal to zero outside boundary as discussed in the fourth paper cited below. 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 Gaussian kernels are added also if s < boundary, as the tails of these Gaussian kernels 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.

For systems with multiple CVs that share identical properties, PBMetaD with partitioned families can be used to group them under one bias potential that each contributes to as discussed in the penultimate paper in the references below. This is done with a list of PF keywords, where each PF* argument contains the list of CVs from ARG to be placed in that family. Once invoked, each CV in ARG must be placed in exactly one PF, even if it results in families containing only one CV. Additionally, in cases where each of SIGMA or GRID entry would correspond to each ARG entry, they now correspond to each PF and must be adjusted accordingly.

The multiple walkers that is discussed in the final paper cited below can also be used. See below the examples.

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.

Click on the labels of the actions for more information on what each action computes
tested on2.11
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=3,5 LABELa label for the action so that its output can be referenced in the input to other actions=d1
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=2,4 LABELa label for the action so that its output can be referenced in the input to other actions=d2
PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d1,d2 SIGMAthe widths of the Gaussian hills=0.2,0.2 HEIGHTthe height of the Gaussian hills, one for all biases=0.3 PACEthe frequency for hill addition, one for all biases=500 LABELa label for the action so that its output can be referenced in the input to other actions=pb FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS_d1,HILLS_d2
PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=d1,d2,pb.bias STRIDE the frequency with which the quantities of interest should be output=100 FILEthe name of the file on which to output these quantities=COLVAR

(See also DISTANCE and PRINT).

If you use well-tempered metadynamics, you should specify a single bias factor and initial Gaussian height.

Click on the labels of the actions for more information on what each action computes
tested on2.11
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=3,5 LABELa label for the action so that its output can be referenced in the input to other actions=d1
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=2,4 LABELa label for the action so that its output can be referenced in the input to other actions=d2
PBMETADUsed to performed Parallel Bias metadynamics. More details ...
ARGthe labels of the scalars on which the bias will act=d1,d2 SIGMAthe widths of the Gaussian hills=0.2,0.2 HEIGHTthe height of the Gaussian hills, one for all biases=0.3
PACEthe frequency for hill addition, one for all biases=500 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 LABELa label for the action so that its output can be referenced in the input to other actions=pb
FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS_d1,HILLS_d2
... PBMETAD
PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=d1,d2,pb.bias STRIDE the frequency with which the quantities of interest should be output=100 FILEthe name of the file on which to output these quantities=COLVAR

Using partitioned families, each CV in ARG must be in exactly one family. Here, the distance between atoms 1,2 is degenerate with 2,4, but not with the distance between 3,5. Note that two SIGMA are provided to match the two PF.

Click on the labels of the actions for more information on what each action computes
tested on2.11
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=3,5 LABELa label for the action so that its output can be referenced in the input to other actions=d1
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=2,4 LABELa label for the action so that its output can be referenced in the input to other actions=d2
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=1,2 LABELa label for the action so that its output can be referenced in the input to other actions=d3
PBMETADUsed to performed Parallel Bias metadynamics. More details ...
ARGthe labels of the scalars on which the bias will act=d1,d2,d3 SIGMAthe widths of the Gaussian hills=0.2,0.2 HEIGHTthe height of the Gaussian hills, one for all biases=0.3
PF0specify which CVs belong in a partitioned family=d1 PF1specify which CVs belong in a partitioned family=d2,d3
PACEthe frequency for hill addition, one for all biases=500 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 LABELa label for the action so that its output can be referenced in the input to other actions=pb
FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS_d1,HILLS_d2
... PBMETAD
PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=d1,d2,d3,pb.bias STRIDE the frequency with which the quantities of interest should be output=100 FILEthe name of the file on which to output these quantities=COLVAR

The following input enables the MPI version of multiple-walkers.

Click on the labels of the actions for more information on what each action computes
tested on2.11
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=3,5 LABELa label for the action so that its output can be referenced in the input to other actions=d1
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=2,4 LABELa label for the action so that its output can be referenced in the input to other actions=d2
PBMETADUsed to performed Parallel Bias metadynamics. More details ...
ARGthe labels of the scalars on which the bias will act=d1,d2 SIGMAthe widths of the Gaussian hills=0.2,0.2 HEIGHTthe height of the Gaussian hills, one for all biases=0.3
PACEthe frequency for hill addition, one for all biases=500 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 LABELa label for the action so that its output can be referenced in the input to other actions=pb
FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS_d1,HILLS_d2
WALKERS_MPI Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR
... PBMETAD
PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=d1,d2,pb.bias STRIDE the frequency with which the quantities of interest should be output=100 FILEthe name of the file on which to output these quantities=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

Click on the labels of the actions for more information on what each action computes
tested on2.11
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=3,5 LABELa label for the action so that its output can be referenced in the input to other actions=d1
DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=2,4 LABELa label for the action so that its output can be referenced in the input to other actions=d2
PBMETADUsed to performed Parallel Bias metadynamics. More details ...
ARGthe labels of the scalars on which the bias will act=d1,d2 SIGMAthe widths of the Gaussian hills=0.2,0.2 HEIGHTthe height of the Gaussian hills, one for all biases=0.3
PACEthe frequency for hill addition, one for all biases=500 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 LABELa label for the action so that its output can be referenced in the input to other actions=pb
FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS_d1,HILLS_d2
WALKERS_Nnumber of walkers=10
WALKERS_IDwalker id=3
WALKERS_DIRshared directory with the hills files from all the walkers=../
WALKERS_RSTRIDEstride for reading hills files=100
... PBMETAD
PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=d1,d2,pb.bias STRIDE the frequency with which the quantities of interest should be output=100 FILEthe name of the file on which to output these quantities=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.

References

More information about how this action can be used is available in the following articles:

Syntax

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

Keyword Type Default Description
ARG input none the labels of the scalars on which the bias will act
PF input none specify which CVs belong in a partitioned family
SIGMA compulsory none the widths of the Gaussian hills
PACE compulsory none the frequency for hill addition, one for all biases
NUMERICAL_DERIVATIVES optional false calculate the derivatives for these quantities numerically
FILE optional not used files in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found
HEIGHT optional not used the height of the Gaussian hills, one for all biases
FMT optional not used specify format for HILLS files (useful for decrease the number of digits in regtests)
BIASFACTOR optional not used use well tempered metadynamics with this bias factor, one for all biases
TEMP optional not used the system temperature - this is only needed if you are doing well-tempered metadynamics
TAU optional not used in well tempered metadynamics, sets height to (k_B Delta Tpacetimestep)/tau
GRID_MIN optional not used the lower bounds for the grid
GRID_MAX optional not used the upper bounds for the grid
GRID_BIN optional not used the number of bins for the grid
GRID_SPACING optional not used the approximate grid spacing (to be used as an alternative or together with GRID_BIN)
GRID_SPARSE optional false use a sparse grid to store hills
GRID_NOSPLINE optional false don't use spline interpolation with grids
GRID_WSTRIDE optional not used frequency for dumping the grid
GRID_WFILES optional not used dump grid for the bias, default names are used if GRID_WSTRIDE is used without GRID_WFILES
GRID_RFILES optional not used read grid for the bias
ADAPTIVE optional not used use a geometric (=GEOM) or diffusion (=DIFF) based hills width scheme
SIGMA_MAX optional not used the upper bounds for the sigmas (in CV units) when using adaptive hills
SIGMA_MIN optional not used the lower bounds for the sigmas (in CV units) when using adaptive hills
SELECTOR optional not used add forces and do update based on the value of SELECTOR
SELECTOR_ID optional not used value of SELECTOR
WALKERS_ID optional not used walker id
WALKERS_N optional not used number of walkers
WALKERS_DIR optional not used shared directory with the hills files from all the walkers
WALKERS_RSTRIDE optional not used stride for reading hills files
WALKERS_MPI optional false Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR
INTERVAL_MIN optional not used one dimensional lower limits, outside the limits the system will not feel the biasing force
INTERVAL_MAX optional not used one dimensional upper limits, outside the limits the system will not feel the biasing force
RESTART optional not used allows per-action setting of restart (YES/NO/AUTO)
UPDATE_FROM optional not used Only update this action from this time
UPDATE_UNTIL optional not used Only update this action until this time