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 expanded ensembles for the target distribution.
This method is similar to the OPES method (OPES) with expanded ensembles target distribution (replica-exchange-like) [60].
An expanded ensemble is obtained by summing a set of ensembles at slightly different termodynamic conditions, or with slightly different Hamiltonians. Such ensembles can be sampled via methods like replica exchange, or this OPES_EXPANDED bias action. A typical example is a multicanonical simulation, in which a whole range of temperatures is sampled instead of a single one.
In oreder to define an expanded target ensemble we use expansion collective variables (ECVs), \(\Delta u_\lambda(\mathbf{x})\). The bias at step \(n\) is
\[ V_n(\mathbf{x})=-\frac{1}{\beta}\log \left(\frac{1}{N_{\{\lambda\}}}\sum_\lambda e^{-\Delta u_\lambda(\mathbf{x})+\beta\Delta F_n(\lambda)}\right)\, . \]
See Ref.[60] for more details on the method.
Notice that the estimates in the DELTAFS file are expressed in energy units, and should be multiplied by \(\beta\) to be dimensionless as in Ref.[60]. The DELTAFS file also contains an estimate of \(c(t)=\frac{1}{\beta} \log \langle e^{\beta V}\rangle\). Similarly to OPES_METAD, it is printed only for reference, since \(c(t)\) reaches a fixed value as the bias converges, and should NOT be used for reweighting. Its value is also needed for restarting a simulation.
You can store the instantaneous \(\Delta F_n(\lambda)\) estimates also in a more readable format using STATE_WFILE and STATE_WSTRIDE. Restart can be done either from a DELTAFS file or from a STATE_RFILE, it is equivalent.
Contrary to OPES_METAD, OPES_EXPANDED does not use kernel density estimation.
# simulate multiple temperatures, as in parallel tempering ene: ENERGY ecv: ECV_MULTITHERMALARG=enecompulsory keyword the label of the internal energy of the system.TEMP_MAX=1000 opes: OPES_EXPANDEDthe maximum of the temperature rangeARG=ecv.*compulsory keyword the label of the ECVs that define the expansion.PACE=500compulsory keyword how often the bias is updatedYou can easily combine multiple ECVs. The OPES_EXPANDED bias will create a multidimensional target grid to sample all the combinations.
Click on the labels of the actions for more information on what each action computes# simulate multiple temperatures while biasing a CV ene: ENERGY dst: DISTANCEATOMS=1,2 ecv1: ECV_MULTITHERMALthe pair of atom that we are calculating the distance between.ARG=enecompulsory keyword the label of the internal energy of the system.TEMP_SET_ALL=200,300,500,1000 ecv2: ECV_UMBRELLAS_LINEmanually set all the temperaturesARG=dstthe input for this action is the scalar output from one or more other actions.CV_MIN=1.2compulsory keyword the minimum of the CV range to be exploredCV_MAX=4.3compulsory keyword the maximum of the CV range to be exploredSIGMA=0.5 opes: OPES_EXPANDEDcompulsory keyword sigma of the umbrella GaussiansARG=ecv1.*,ecv2.*compulsory keyword the label of the ECVs that define the expansion.PACE=500compulsory keyword how often the bias is updatedOBSERVATION_STEPS=1compulsory keyword ( default=100 ) number of unbiased initial PACE steps to collect statistics for initializationIf an ECV is based on more than one CV you must provide all the output component, in the proper order. You can use Regular Expressions for that, like in the following example.
Click on the labels of the actions for more information on what each action computes# simulate multiple temperatures and pressures while biasing a two-CVs linear path ene: ENERGY vol: VOLUME ecv_mtp: ECV_MULTITHERMAL_MULTIBARIC ...ARG=ene,volcompulsory keyword the labels of the potential energy and of the volume of the system.TEMP=300compulsory keyword ( default=-1 ) temperature.TEMP_MIN=200the minimum of the temperature rangeTEMP_MAX=800the maximum of the temperature rangePRESSURE=0.06022140857*1000 #1 kbarcompulsory keyword pressure.PRESSURE_MIN=0the minimum of the pressure rangePRESSURE_MAX=0.06022140857*2000 #2 kbar ... cv1: DISTANCEthe maximum of the pressure rangeATOMS=1,2 cv2: DISTANCEthe pair of atom that we are calculating the distance between.ATOMS=3,4 ecv_umb: ECV_UMBRELLAS_LINEthe pair of atom that we are calculating the distance between.ARG=cv1,cv2the input for this action is the scalar output from one or more other actions.TEMP=300compulsory keyword ( default=-1 ) temperature.CV_MIN=0.1,0.1compulsory keyword the minimum of the CV range to be exploredCV_MAX=1.5,1.5compulsory keyword the maximum of the CV range to be exploredSIGMA=0.2compulsory keyword sigma of the umbrella GaussiansBARRIER=70 opes: OPES_EXPANDEDa guess of the free energy barrier to be overcome (better to stay higher than lower)ARG=(ecv_.*)compulsory keyword the label of the ECVs that define the expansion.PACE=500compulsory keyword how often the bias is updatedWALKERS_MPI( default=off ) switch on MPI version of multiple walkersPRINT_STRIDE=1000compulsory keyword ( default=100 ) stride for printing to DELTAFS file, in units of PACE
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 |
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 |
ARG | the label of the ECVs that define the expansion. You can use an * to make sure all the output components of the ECVs are used, as in the examples above |
PACE | how often the bias is updated |
OBSERVATION_STEPS | ( default=100 ) number of unbiased initial PACE steps to collect statistics for initialization |
FILE | ( default=DELTAFS ) a file with the estimate of the relative Delta F for each component of the target and of the global c(t) |
PRINT_STRIDE | ( default=100 ) stride for printing to DELTAFS file, in units of PACE |
NUMERICAL_DERIVATIVES | ( default=off ) calculate the derivatives for these quantities numerically |
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 |
FMT | specify format for DELTAFS file |
STATE_RFILE | read from this file the Delta F estimates and all the info needed to RESTART the simulation |
STATE_WFILE | write to this file the Delta F estimates 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) |
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 |