Metadynamics
- class mlptrain.sampling.metadynamics.Metadynamics(cvs: Sequence['_PlumedCV'] | '_PlumedCV', bias: 'mlptrain.PlumedBias' | None = None, temp: float | None = None)
Bases:
objectMetadynamics class for running biased molecular dynamics using metadynamics bias and analysing the results
- __init__(cvs: Sequence['_PlumedCV'] | '_PlumedCV', bias: 'mlptrain.PlumedBias' | None = None, temp: float | None = None)
Molecular dynamics using metadynamics bias. Used for calculating free energies (by using well-tempered metadynamics bias) and sampling configurations for active learning.
variables are required, e.g. Metadynamics(cvs=cv1, bias=PlumedBias(cvs=(cv1, cv2)), where cv1 will be biased using metadynamics, and cv2 might be an additional CV with WALLS to constrain the system
- block_analysis(start_time: float, idx: int = 1, energy_units: str = 'kcal mol-1', n_bins: int = 300, min_n_blocks: int = 10, min_blocksize: int = 10, blocksize_interval: int = 10, bandwidth: float = 0.02, cvs_bounds: Sequence | None = None, temp: float | None = None, dt: float | None = None, interval: int | None = None) None
Perform block averaging analysis on the sliced trajectory of the most recent metadynamics run. Plot the block analysis and save mean FES grids with a range of block sizes, which, if the block analysis converged, can be used for plotting the FES error using plot_fes() method.
(in ps)
- Parameters:
idx – (int) Index specifying which metadynamics run (from n_runs) to use for block analysis
energy_units – (str) Energy units to be used in plotting, available units: ‘eV’, ‘kcal mol-1’, ‘kJ mol-1’
n_bins – (int) Number of bins to use when dumping histograms
min_n_blocks – (int) Minimum number of blocks
min_blocksize – (int) Minimum size of a block
blocksize_interval – (int) Block size interval describing the difference in the number of frames per block for adjacent block sizes
bandwidth – (float) Bandwidth parameter used in the histogram estimations
cvs_bounds – Specifies the range between which to compute the free energy for each collective variable, e.g. [(0.8, 1.5), (80, 120)]
temp – (float) Temperature in K during the analysed simulation
dt – (float) Time-step in fs during the analysed simulation
interval – (int) Interval between saving the geometry during the analysed simulation
- compute_fes(energy_units: str = 'kcal mol-1', n_bins: int = 300, cvs_bounds: Sequence | None = None, via_reweighting: bool = False, start_time: float = 0.0, bandwidth: float = 0.02, temp: float | None = None, dt: float | None = None, interval: int | None = None) ndarray
Compute a grid containing collective variable grids and free energy surface grids, which is saved in the current directory as .npy file.
- Parameters:
n_bins – (int) Number of bins to use in every dimension for fes file generation from HILLS
cvs_bounds – Specifies the range between which to compute the free energy for each collective variable e.g. [(0.8, 1.5), (80, 120)]
via_reweighting – (bool) If true the free energy is computed using the reweighting scheme
start_time – (float) Start time of the sliced trajectory which is going to be used in reweighting (in ps)
bandwidth – (float) Bandwidth parameter used in the histogram estimations
temp – (float) Temperature in K during the analysed simulation, only needed when reweighting
dt – (float) Time-step in fs during the analysed simulation, only needed when reweighting
interval – (int) Interval between saving the geometry during the analysed simulation, only needed when reweighting
- Returns:
The total grid containing CVs and FESs
- Return type:
(np.ndarray)
- estimate_width(configurations: 'mlptrain.Configuration' | 'mlptrain.ConfigurationSet', mlp: MLPotential, temp: float = 300, interval: int = 10, dt: float = 1, plot: bool = True, **kwargs) List
Estimate optimal widths (σ) to be used in metadynamics.
- Parameters:
mlp – Machine learnt potential
temp – (float) Temperature in K to initialise velocities and to run NVT MD. Must be positive
interval – (int) Interval between saving the geometry
dt – (float) Time-step in fs
plot – (bool) If True plots trajectories of collective variables as a function of time
- List of optimal width values for each CV,
e.g. [0.02, 0.03] for two CVs
- Return type:
(List)
- property kbt: float
Value of k_B*T in ASE units
- property n_cvs: int
Number of collective variables used in metadynamics
- plot_fes(energy_units: str = 'kcal mol-1', confidence_level: float = 0.95, n_bins: int = 300, cvs_bounds: Sequence | None = None, fes_npy: str | None = None, blocksize: int | None = None) None
Plot the free energy surface with a confidence interval. If the .npy file is not supplied, the file is computed (if metadynamics has been run). If blocksize is supplied, the FES error is plotted using block_analysis.npz file with the given blocksize generated during block averaging analysis.
‘kJ mol-1’. If fes_npy and/or blocksize options are used, the units of fes_npy and/or blocksize files must match energy units here
- Parameters:
confidence_level – (float) Specifies what confidence level to use in plots (probability for FES to lie in the plotted range)
n_bins – (int) Number of bins to use in every dimension for fes file generation from HILLS
cvs_bounds – Specifies the range between which to compute the free energy for each collective variable, e.g. [(0.8, 1.5), (80, 120)]
fes_npy – (str) File name of the .npy file used for plotting. FES obtained through block analysis should be specified using blocksize
blocksize – (int) If block analysis has been performed, the integer specifies which block size to use for plotting the FES standard error
- plot_fes_convergence(stride: int, n_surfaces: int = 5, time_units: str = 'ps', energy_units: str = 'kcal mol-1', n_bins: int = 300, cvs_bounds: Sequence | None = None, idx: int = 1) None
Compute multiple fes.dat files from a HILLS_idx.dat file by summing the deposited gaussians using a stride. Use the computed files to plot multiple FESs as a function of simulation time.
- Parameters:
n_surfaces – (int) Number of surfaces to be plotted (counting from the last computed surface), must not exceed the number of computed surfaces
time_units – (str) Time units to be used in plotting, available units: ‘fs’, ‘ps’, ‘ns’
energy_units – (str) Energy units to be used in plotting, available units: ‘eV’, ‘kcal mol-1’, ‘kJ mol-1’
n_bins – (int) Number of bins to use in every dimension for fes file generation from HILLS
cvs_bounds – Specifies the range between which to compute the free energy for each collective variable, e.g. [(0.8, 1.5), (80, 120)]
idx – (int) Integer which specifies which metadynamics run to use for plotting the FES convergence
- plot_gaussian_heights(energy_units: str = 'kcal mol-1', time_units: str = 'ps', path: str = 'plumed_files/metadynamics') None
Plot the height of deposited gaussians as a function of time (using HILLS_{idx}.dat files).
- Parameters:
time_units – (str) Time units to be used in plotting, available units: ‘fs’, ‘ps’, ‘ns’
path – (str) Directory where HILLS_{idx}.dat files are located
- run_metadynamics(configuration: mlptrain.Configuration, mlp: MLPotential, temp: float, interval: int, dt: float, pace: int = 100, height: float | None = None, width: list[float] | float | None = None, biasfactor: float | None = None, al_iter: int | None = None, n_runs: int = 1, save_sep: bool = True, all_to_xyz: bool = False, restart: bool = False, **kwargs) None
Perform multiple metadynamics runs in parallel, generate .xyz and .traj files containing trajectories of the runs, generate PLUMED files containing deposited gaussians and trajectories in terms of the CVs.
- Parameters:
temp – (float) Temperature in K to initialise velocities and to run NVT MD. Must be positive
interval – (int) Interval between saving the geometry
dt – (float) Time-step in fs
pace – (int) τ_G/dt, interval at which a new gaussian is placed
height – (float) ω, initial height of placed gaussians (in eV)
width – (List[float] | float | None) σ, standard deviation (parameter describing the width) of placed gaussians, if not supplied it is estimated automatically
biasfactor – (float | None) γ, describes how quickly gaussians shrink, larger values make gaussians to be placed less sensitive to the bias potential
al_iter – (int | None) If not None, can be used to load metadynamics bias generated during inherited-bias active learning. The index specifies from which AL iteration to load the bias
n_runs – (int) Number of times to run metadynamics
save_sep – (bool) If True saves trajectories of each simulation separately
all_to_xyz – (bool) If True all .traj trajectory files are saved as .xyz files (when using save_fs, save_ps, save_ns)
restart – (bool) If True restarts the most recent metadynamics simulation
in some units
- Keyword Arguments:
grid_bin} ({grid_min, grid_max,) – Grid parameters: if defined, the deposited gaussians are dumped to a grid that is used to compute the biased forces
constraints – (List) List of ASE constraints to use in the dynamics e.g. [ase.constraints.Hookean(a1, a2, k, rt)]
- try_multiple_biasfactors(configuration: mlptrain.Configuration, mlp: MLPotential, temp: float, interval: int, dt: float, biasfactors: Sequence[float], pace: int = 500, height: float | None = None, width: list[float] | float | None = None, plotted_cvs: Sequence[_PlumedCV] | None = None, **kwargs) None
Execute multiple well-tempered metadynamics runs in parallel with a provided sequence of biasfactors and plot the resulting trajectories, useful for estimating the optimal biasfactor value.
- Parameters:
temp – (float) Temperature in K to initialise velocities and to run NVT MD. Must be positive
interval – (int) Interval between saving the geometry
dt – (float) Time-step in fs
biasfactors – Sequence of γ, describes how quickly gaussians shrink, larger values make gaussians to be placed less sensitive to the bias potential
pace – (int) τ_G/dt, interval at which a new gaussian is placed
height – (float) ω, initial height of placed gaussians (in eV)
width – (List[float] | float | None) σ, standard deviation (parameter describing the width) of placed gaussians, if not supplied it is estimated automatically
plotted_cvs – Sequence of one or two PlumedCV objects which are going to be plotted, must be a subset for CVs used to define the Metadynamics object
deposited gaussians are dumped to a grid that is used to compute the biased forces
- Keyword Arguments:
constraints – (List) List of ASE constraints to use in the dynamics e.g. [ase.constraints.Hookean(a1, a2, k, rt)]