Action: RESCALE
| Module | isdb |
|---|---|
| Description | Usage |
| Scales the value of an another action, being a Collective Variable or a Bias. |
Details and examples
Scales the value of an another action, being a Collective Variable or a Bias.
The rescaling factor is determined by a parameter defined on a logarithmic grid of dimension NBIN in the range from 1 to MAX_RESCALE. The current value of the rescaling parameter is stored and shared across other actions using a SELECTOR. A Monte Carlo procedure is used to update the value of the rescaling factor every MC_STRIDE steps of molecular dynamics. Well-tempered metadynamics, defined by the parameters W0 and BIASFACTOR, is used to enhance the sampling in the space of the rescaling factor. The well-tempered metadynamics bias potential is written to the file BFILE every BSTRIDE steps and read when restarting the simulation using the directive RESTART.
scaling arguments
Additional arguments not to be scaled, one for each bin in the rescaling parameter ladder, can be provided at the end of the ARG list. The number of such arguments is specified by the option NOT_RESCALED. These arguments will be not be scaled, but they will be considered as bias potentials and used in the computation of the Metropolis acceptance probability when proposing a move in the rescaling parameter. See example below.
multiple replicas
If PLUMED is running in a multiple-replica framework (for example using the -multi option in GROMACS), the arguments will be summed across replicas, unless the NOT_SHARED option is used. Also, the value of the SELECTOR will be shared and thus will be the same in all replicas.
Examples
In this example we use RESCALE to implement a simulated-tempering like approach. The total potential energy of the system is scaled by a parameter defined on a logarithmic grid of 5 bins in the range from 1 to 1.5. A well-tempered metadynamics bias potential is used to ensure diffusion in the space of the rescaling parameter.
ene: ENERGYCalculate the total potential energy of the simulation box. More details
SELECTORDefines a variable (of the type double) inside the PLUMED code that can be used and modified by other actions. More details NAMEname of the SELECTOR=GAMMA VALUEset (initial) value of the SELECTOR=0 RESCALEScales the value of an another action, being a Collective Variable or a Bias. More details ... LABELa label for the action so that its output can be referenced in the input to other actions=res ARGthe labels of the scalars on which the bias will act=ene TEMPtemperature=300 SELECTORname of the SELECTOR used for rescaling=GAMMA MAX_RESCALEmaximum values for rescaling=1.5 NBINnumber of bins for gamma grid=5 W0initial bias height=1000 BIASFACTORbias factor=100.0 BSTRIDEstride for writing bias=2000 BFILEfile name for bias=bias.dat ...
PRINTPrint quantities to a file. More details FILEthe name of the file on which to output these quantities=COLVAR ARGthe labels of the values that you would like to print to the file=* STRIDE the frequency with which the quantities of interest should be output=100
In this second example, we add to the simulated-tempering approach introduced above one Parallel Bias metadynamics simulation (see PBMETAD) for each value of the rescaling parameter. At each moment of the simulation, only one of the PBMETAD actions is activated, based on the current value of the associated SELECTOR. The PBMETAD bias potentials are not scaled, but just used in the calculation of the Metropolis acceptance probability when proposing a move in the rescaling parameter.
ene: ENERGYCalculate the total potential energy of the simulation box. More details d: DISTANCECalculate the distance/s between pairs of atoms. More details ATOMSthe pair of atom that we are calculating the distance between=1,2 SELECTORDefines a variable (of the type double) inside the PLUMED code that can be used and modified by other actions. More details NAMEname of the SELECTOR=GAMMA VALUEset (initial) value of the SELECTOR=0 pbmetad0: PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d SELECTORadd forces and do update based on the value of SELECTOR=GAMMA SELECTOR_IDvalue of SELECTOR=0 SIGMAthe widths of the Gaussian hills=0.1 PACEthe frequency for hill addition, one for all biases=500 HEIGHTthe height of the Gaussian hills, one for all biases=1 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS.0 pbmetad1: PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d SELECTORadd forces and do update based on the value of SELECTOR=GAMMA SELECTOR_IDvalue of SELECTOR=1 SIGMAthe widths of the Gaussian hills=0.1 PACEthe frequency for hill addition, one for all biases=500 HEIGHTthe height of the Gaussian hills, one for all biases=1 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS.1 pbmetad2: PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d SELECTORadd forces and do update based on the value of SELECTOR=GAMMA SELECTOR_IDvalue of SELECTOR=2 SIGMAthe widths of the Gaussian hills=0.1 PACEthe frequency for hill addition, one for all biases=500 HEIGHTthe height of the Gaussian hills, one for all biases=1 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS.2 pbmetad3: PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d SELECTORadd forces and do update based on the value of SELECTOR=GAMMA SELECTOR_IDvalue of SELECTOR=3 SIGMAthe widths of the Gaussian hills=0.1 PACEthe frequency for hill addition, one for all biases=500 HEIGHTthe height of the Gaussian hills, one for all biases=1 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS.3 pbmetad4: PBMETADUsed to performed Parallel Bias metadynamics. More details ARGthe labels of the scalars on which the bias will act=d SELECTORadd forces and do update based on the value of SELECTOR=GAMMA SELECTOR_IDvalue of SELECTOR=4 SIGMAthe widths of the Gaussian hills=0.1 PACEthe frequency for hill addition, one for all biases=500 HEIGHTthe height of the Gaussian hills, one for all biases=1 BIASFACTORuse well tempered metadynamics with this bias factor, one for all biases=8 FILEfiles in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found=HILLS.4 RESCALEScales the value of an another action, being a Collective Variable or a Bias. More details ... LABELa label for the action so that its output can be referenced in the input to other actions=res TEMPtemperature=300 ARGthe labels of the scalars on which the bias will act=ene,pbmetad0.bias,pbmetad1.bias,pbmetad2.bias,pbmetad3.bias,pbmetad4.bias SELECTORname of the SELECTOR used for rescaling=GAMMA MAX_RESCALEmaximum values for rescaling=1.5 NOT_RESCALEDthese last N arguments will not be scaled=5 NBINnumber of bins for gamma grid=5 W0initial bias height=1000 BIASFACTORbias factor=100.0 BSTRIDEstride for writing bias=2000 BFILEfile name for bias=bias.dat ...
PRINTPrint quantities to a file. More details FILEthe name of the file on which to output these quantities=COLVAR ARGthe labels of the values that you would like to print to the file=* STRIDE the frequency with which the quantities of interest should be output=100
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 |
Output components
This action calculates the values in the following table. These values can be referenced elsewhere in the input by using this Action's label followed by a dot and the name of the value required from the list below.
| Name | Type | Description |
|---|---|---|
| bias | scalar | the instantaneous value of the bias potential |
| igamma | scalar | gamma parameter |
| accgamma | scalar | MC acceptance for gamma |
| wtbias | scalar | well-tempered bias |
Full list of keywords
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 |
| TEMP | compulsory | none | temperature |
| SELECTOR | compulsory | none | name of the SELECTOR used for rescaling |
| MAX_RESCALE | compulsory | none | maximum values for rescaling |
| NBIN | compulsory | none | number of bins for gamma grid |
| W0 | compulsory | none | initial bias height |
| BIASFACTOR | compulsory | none | bias factor |
| BSTRIDE | compulsory | none | stride for writing bias |
| BFILE | compulsory | none | file name for bias |
| NUMERICAL_DERIVATIVESThis keyword do not have examples | optional | false | calculate the derivatives for these quantities numerically |
| NOT_SHAREDThis keyword do not have examples | optional | not used | list of arguments (from 1 to N) not summed across replicas |
| NOT_RESCALED | optional | not used | these last N arguments will not be scaled |
| MC_STEPSThis keyword do not have examples | optional | not used | number of MC steps |
| MC_STRIDEThis keyword do not have examples | optional | not used | MC stride |
| PACEThis keyword do not have examples | optional | not used | Pace for adding bias, in MC stride unit |