metadynamics
|
Implements a metadynamics update scheme.
This class implements an integration scheme for metadynamics.
For a review of metadynamics, a free-energy measurement technique, see Barducci et al., Metadynamics, Wiley Interdiscipl. Rev.: Comput. Mol. Sci. 5, pp. 826-843 (2011).
Some of the features of this class are loosely inspired by the PLUMED plugin for Metadynamics, http://www.plumed-code.org/.
The IntegratorMetaDynamics class takes a number of CollectiveVariables, and uses their values during the course of the simulation to update the bias potential, a sum of Gaussians, which disfavors revisiting visited states. The derivative of the bias potential (with respect to the collective coordinate) in turn multiplies the force on the particles as defined by the collective variable.
The bias potential is updated every stride steps.
Well-tempered metadynamics is supported by default, a value of T_shift defines the temperature shift entering the scaling factor for the bias potential (see above ref. for details). Very large values of T_shift correspond to standard metadynamics, T_shift->0 corresponds to standard molecular dynamics.
Two modes of operation are supported: in the first mode, the Gaussians are resummed every time step. With increasing simulation time, this process will take longer, and slow down the simulation. In the second mode, a grid mode, a current grid of values of the collective variables is maintained and updated whenever a new Gaussian is deposited. The value of the bias potential and its derivatives of the potential are calculated by using numerical interpolation. In this mode, the time per step is O(1), so this mode is preferrable for long simulations.
It is possible to output the grid after the simulation, and to restart from the grid file. It is also possible to restart from the grid file and turn off the deposition of new Gaussians, e.g. to equilibrate the system in the bias potential landscape and measure the histogram of the collective variable, to correct for errors.
Public Member Functions | |
IntegratorMetaDynamics (boost::shared_ptr< SystemDefinition > sysdef, Scalar deltaT, Scalar W, Scalar T_shift, unsigned int stride, bool add_hills=true, const std::string &filename="", bool overwrite=false, bool use_grid=false) | |
virtual void | update (unsigned int timestep) |
virtual void | prepRun (unsigned int timestep) |
virtual void | printStats () |
void | registerCollectiveVariable (boost::shared_ptr< CollectiveVariable > cv, Scalar sigma, Scalar cv_min=Scalar(0.0), Scalar cv_max=Scalar(0.0), int num_points=0) |
void | removeAllVariables () |
std::vector< std::string > | getProvidedLogQuantities () |
Scalar | getLogValue (const std::string &quantity, unsigned int timestep) |
void | setGrid (bool use_grid) |
bool | isInitialized () |
void | dumpGrid (const std::string &filename) |
void | restartFromGridFile (const std::string &filename) |
void | setAddHills (bool add_hills) |
Private Member Functions | |
void | updateBiasPotential (unsigned int timestep) |
Internal helper function to update the bias potential. | |
void | openOutputFile () |
Helper function to open output file for logging. | |
void | writeFileHeader () |
Helper function to write file header. | |
void | setupGrid () |
Helper function to initialize the grid. | |
Scalar | interpolateBiasPotential (const std::vector< Scalar > &val) |
Helper function to get value of bias potential by multilinear interpolation. | |
void | readGrid (const std::string &filename) |
Helper function to read in data from grid file. |
Private Attributes | |
Scalar | m_W |
Height of Gaussians. | |
Scalar | m_T_shift |
Temperature shift. | |
unsigned int | m_stride |
Number of timesteps between Gaussian depositions. | |
std::vector < CollectiveVariableItem > | m_variables |
The list of collective variables. | |
std::vector< std::vector < Scalar > > | m_cv_values |
History of CV values. | |
unsigned int | m_num_update_steps |
Number of update steps performed in this run thus far. | |
std::vector< std::string > | m_log_names |
Names of logging quantities. | |
std::string | m_suffix |
Suffix for unbiased variables. | |
Scalar | m_curr_bias_potential |
The sum of Gaussians. | |
std::vector< Scalar > | m_bias_potential |
List of values of the bias potential. | |
bool | m_is_initialized |
True if history-dependent potential has been initialized. | |
const std::string | m_filename |
Name of output file. | |
bool | m_overwrite |
True if the file should be overwritten. | |
bool | m_is_appending |
True if we are appending to the file. | |
std::ofstream | m_file |
Output log file. | |
std::string | m_delimiter |
Delimiting string. | |
bool | m_use_grid |
True if we are using a grid. | |
std::vector< Scalar > | m_grid |
d-dimensional grid to store values of bias potential | |
IndexGrid | m_grid_index |
Indexer for the d-dimensional grid. | |
bool | m_add_hills |
True if hills should be added during the simulation. | |
std::string | m_restart_filename |
Name of file to read restart information from. |
IntegratorMetaDynamics | ( | boost::shared_ptr< SystemDefinition > | sysdef, |
Scalar | deltaT, | ||
Scalar | W, | ||
Scalar | T_shift, | ||
unsigned int | stride, | ||
bool | add_hills = true , |
||
const std::string & | filename = "" , |
||
bool | overwrite = false , |
||
bool | use_grid = false |
||
) |
Constructor
sysdef | System definition |
deltaT | Time step |
W | Weight of Gaussians (multiplicative factor) |
T_shift | Temperature shift for well-tempered metadynamics |
stride | Number of time steps between deposition of Gaussians |
add_hills | True if new Gaussians should be added during the simulation |
filename | Name of file to output hill information to |
overwrite | True f the file should be overwritten when it exists already |
use_grid | True if grid should be used |
void dumpGrid | ( | const std::string & | filename | ) |
Output the grid to a file
filename | Name of file to dump grid to |
|
inline |
Obtain the value of a specific log quantity
quantity | The quantity to obtain the value of |
timestep | The current value of the time step |
|
inline |
Returns the names of log quantitites provided by this integrator
|
inline |
Returns true if the integration has been initialized (i.e. the simulation has been run at least once)
|
virtual |
Prepare for the run
timestep | The current value of the timestep |
|
inlinevirtual |
Output statistics at end of run (unimplemented as of now)
|
inline |
Register a new collective variable
cv | The collective variable |
sigma | The standard deviation of Gaussians for this collective variable |
cv_min | Minimum value, if using grid |
cv_max | Maximum value, if using grid |
num_points | Number of grid points to use for interpolation |
|
inline |
Remove all collective variables
|
inline |
Restart from a grid file. Upon running the simulation, the information will be read in.
filename | The name of the file that contains the grid information |
|
inline |
Set a flag that controls deposition of new Gaussian hills
add_hills | True if hills should be generated during the simulation |
void setGrid | ( | bool | use_grid | ) |
Set up grid interpolation
use_grid | True if the grid is to be used |
|
virtual |
Sample collective variables, update bias potential and derivatives
timestep | The current value of the timestep |