OPES_METAD_EXPLORE
This is part of the opes module
It is only available if you configure PLUMED with ./configure –enable-modules=opes . Furthermore, this feature is still being developed so take care when using it and report any problems on the mailing list.

On-the-fly probability enhanced sampling with well-tempered target distribution in exploreation mode.

On-the-fly probability enhanced sampling with well-tempered target distribution (OPES) with well-tempered target distribution, exploration mode [59] .

This OPES_METAD_EXPLORE action samples the well-tempered target distribution, that is defined via its marginal \(p^{\text{WT}}(\mathbf{s})\propto [P(\mathbf{s})]^{1/\gamma}\) over some collective variables (CVs), \(\mathbf{s}=\mathbf{s}(\mathbf{x})\). While OPES_METAD does so by estimating the unbiased distribution \(P(\mathbf{s})\), OPES_METAD_EXPLORE instead estimates on-the-fly the target \(p^{\text{WT}}(\mathbf{s})\) and uses it to define the bias. The bias at step \(n\) is

\[ V_n(\mathbf{s}) = (\gamma-1)\frac{1}{\beta}\log\left(\frac{p^{\text{WT}}_n(\mathbf{s})}{Z_n}+\epsilon\right)\, . \]

See Ref.[59] for a complete description of the method.

Intuitively, while OPES_METAD aims at quickly converging the reweighted free energy, OPES_METAD_EXPLORE aims at quickly sampling the target well-tempered distribution. Given enough simulation time, both will converge to the same bias potential but they do so in a qualitatively different way. Compared to OPES_METAD, OPES_METAD_EXPLORE is more similar to METAD, because it allows the bias to vary significantly, thus enhancing exploration. This goes at the expenses of a typically slower convergence of the reweight estimate. OPES_METAD_EXPLORE can be useful e.g.~for simulating a new system with an unknown BARRIER, or for quickly testing the effectiveness of a new CV that might be degenerate.

Similarly to OPES_METAD, also OPES_METAD_EXPLORE uses a kernel density estimation with the same on-the-fly compression algorithm. The only difference is that the kernels are not weighted, since it estimates the sampled distribution and not the reweighted unbiased one.

All the options of OPES_METAD are also available in OPES_METAD_EXPLORE, except for those that modify the target distribution, since only a well-tempered target is allowed in this case.

Examples

The following is a minimal working example:

Click on the labels of the actions for more information on what each action computes
tested on master
cv: DISTANCE 
ATOMS
the pair of atom that we are calculating the distance between.
=1,2 opes: OPES_METAD_EXPLORE
ARG
compulsory keyword the labels of the scalars on which the bias will act
=cv
PACE
compulsory keyword the frequency for kernel deposition
=500
BARRIER
compulsory keyword the free energy barrier to be overcome.
=40 PRINT
STRIDE
compulsory keyword ( default=1 ) the frequency with which the quantities of interest should be output
=100
FILE
the name of the file on which to output these quantities
=COLVAR
ARG
compulsory keyword the labels of the values that you would like to print to the file
=cv,opes.*
Glossary of keywords and components
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
rct estimate of c(t). log(exp(beta V)/beta, should become flat as the simulation converges. Do NOT use for reweighting
zed estimate of Z_n. should become flat once no new CV-space region is explored
neff effective sample size
nker total number of compressed kernels used to represent the bias

In addition the following quantities can be calculated by employing the keywords listed below

Quantity Keyword Description
work CALC_WORK total accumulated work done by the bias
nlker NLIST number of kernels in the neighbor list
nlsteps NLIST number of steps from last neighbor list update
Compulsory keywords
ARG the labels of the scalars on which the bias will act
TEMP ( default=-1 ) temperature. If not set, it is taken from MD engine, but not all MD codes provide it
PACE the frequency for kernel deposition
SIGMA ( default=ADAPTIVE ) the initial widths of the kernels, divided by the square root of gamma. If not set, an adaptive sigma will be used with the given ADAPTIVE_SIGMA_STRIDE
BARRIER the free energy barrier to be overcome. It is used to set BIASFACTOR, EPSILON, and KERNEL_CUTOFF to reasonable values
COMPRESSION_THRESHOLD ( default=1 ) merge kernels if closer than this threshold, in units of sigma
FILE ( default=KERNELS ) a file in which the list of all deposited kernels is stored
Options
NUMERICAL_DERIVATIVES ( default=off ) calculate the derivatives for these quantities numerically
NLIST ( default=off ) use neighbor list for kernels summation, faster but experimental
NLIST_PACE_RESET ( default=off ) force the reset of the neighbor list at each PACE. Can be useful with WALKERS_MPI
FIXED_SIGMA ( default=off ) do not decrease sigma as the simulation proceeds. Can be added in a RESTART, to keep in check the number of compressed kernels
RECURSIVE_MERGE_OFF ( default=off ) do not recursively attempt kernel merging when a new one is added
NO_ZED ( default=off ) do not normalize over the explored CV space, Z_n=1
STORE_STATES ( default=off ) append to STATE_WFILE instead of ovewriting it each time
CALC_WORK ( default=off ) calculate the total accumulated work done by the bias since last restart
WALKERS_MPI ( default=off ) switch on MPI version of multiple walkers
SERIAL

( default=off ) perform calculations in serial

ADAPTIVE_SIGMA_STRIDE number of steps for measuring adaptive sigma. Default is 10xPACE
SIGMA_MIN never reduce SIGMA below this value
BIASFACTOR the gamma bias factor used for the well-tempered target distribution. Cannot be 'inf'
EPSILON the value of the regularization constant for the probability
KERNEL_CUTOFF truncate kernels at this distance, in units of sigma
NLIST_PARAMETERS ( default=3.0,0.5 ) the two cutoff parameters for the kernels neighbor list
FMT specify format for KERNELS file
STATE_RFILE read from this file the compressed kernels and all the info needed to RESTART the simulation
STATE_WFILE write to this file the compressed kernels and all the info needed to RESTART the simulation
STATE_WSTRIDE number of MD steps between writing the STATE_WFILE. Default is only on CPT events (but not all MD codes set them)
EXCLUDED_REGION kernels are not deposited when the action provided here has a nonzero value, see example above
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