Plumed
- class mlptrain.sampling.plumed.PlumedAverageCV(name: str, atom_groups: Sequence | None = None)
Bases:
_PlumedCVClass used to initialise a PLUMED collective variable as an average between multiple degrees of freedom
- __init__(name: str, atom_groups: Sequence | None = None)
PLUMED collective variable as an average between multiple degrees of freedom (distances, angles, torsions),
e.g. [(0, 1), (2, 3), (4, 5)] gives ζ = 1/3 * (r_01 + r_23 + r_45)
which are used to generate DOFs
- class mlptrain.sampling.plumed.PlumedBias(cvs: Sequence[_PlumedCV] | _PlumedCV | None = None, filename: str | None = None)
Bases:
ASEConstraint@DynamicAttrs Class for storing collective variables and parameters used in biased simulations
- __init__(cvs: Sequence[_PlumedCV] | _PlumedCV | None = None, filename: str | None = None)
Class for storing collective variables and parameters used in biased simulations, parameters are not initialised with the object and have to be defined seperately using PlumedBias methods.
Can be initialised from a complete PLUMED input file as well.
- adjust_forces(atoms, forces) None
Adjust the forces of the system by adding PLUMED forces
- adjust_positions(atoms, newpositions) None
Method required for ASE but not used in mlp-train
- adjust_potential_energy(atoms) float
Adjust the energy of the system by adding PLUMED bias
- property biasfactor_setup: str
String specifying biasfactor in PLUMED input file
- property cv_sequence: str
String containing names of all collective variables separated by commas
- property from_file: bool
Whether the bias is initialised using PLUMED input file
- initialise_for_metad_al(width: Sequence[float] | float, pace: int = 20, height: float | None = None, biasfactor: float | None = None, cvs: Sequence[_PlumedCV] | _PlumedCV | None = None, grid_min: Sequence[float] | float | None = None, grid_max: Sequence[float] | float | None = None, grid_bin: Sequence[float] | float | None = None) None
Initialise PlumedBias for metadynamics active learning by setting the required parameters.
- Parameters:
pace – (int) τ_G/dt, interval at which a new gaussian is placed
height – (float) ω, initial height of placed gaussians (in eV). If not supplied will be set to 5*k_B*T, where T is the temperature at which metadynamics active learning is performed
biasfactor – (float) γ, describes how quickly gaussians shrink, larger values make gaussians to be placed less sensitive to the bias potential
cvs – (mlptrain._PlumedCV) Sequence of PLUMED collective variables which will be biased. If this variable is not set, all CVs attached to self will be biased
grid_min – (float) Lower bound of the grid
grid_max – (float) Upper bound of the grid
grid_bin – (float) Number of bins to use for each collective variable, if not specified PLUMED automatically sets bin width to 1/5 of the gaussian width (σ) value
- property metad_cv_sequence: str
String containing names of collective variables used in metadynamics separated by commas
- property metad_grid_setup: str
String specifying metadynamics grid parameters in PLUMED input file
- property metadynamics: bool
True if any parameters required for metadynamics are set
- property n_cvs: int
Number of collective variables attached to the bias
- property n_metad_cvs: int
Number of collective variables attached to the bias that will be used in metadynamics
- strip() None
Change the bias such that it would only contain the definitions of collective variables and their associated upper and lower walls
- property width_sequence: str
String containing width values separated by commas
- write_cv_files() None
Write files attached to the CVs to the current directory
- class mlptrain.sampling.plumed.PlumedCNCV(name: str, r_ref: float, d_ref: float = 0, n: int = 6, m: int = 12, atom_groups: Sequence | None = None)
Bases:
_PlumedCV- __init__(name: str, r_ref: float, d_ref: float = 0, n: int = 6, m: int = 12, atom_groups: Sequence | None = None)
PLUMED collective variable as a coordination number (CN) between two atoms or groups of atoms
e.g. [(0, 2)] gives ζ =[1-((r_02 - d_ref)/r_ref)^n]/ [1-((r_02 - d_ref)/r_ref)^m]
which corresponds to the CN between atoms 0 and 2. To ensure that the CN has continuous derivatives, we use a rational switching function consistent with PLUMED. More information: https://www.plumed.org/doc-v2.9/user-doc/html/_c_o_o_r_d_i_n_a_t_i_o_n.html
- Parameters:
d_ref – (float) Parmeter d_0 in the original switching function. Default = 0.
n – (int) Parameters used to adjust the shape of the function
m – (int) Parameters used to adjust the shape of the function
atom_groups – (Sequence[Sequence[int]]) List of atom index sequences which are used to generate DOFs
- class mlptrain.sampling.plumed.PlumedCalculator(calc, input, timestep, atoms=None, kT=1.0, log='', restart=False, use_charge=False, update_charge=False)
Bases:
PlumedModified ASE calculator, instead of returning biased energies and forces, this calculator computes unbiased energies and forces, and computes PLUMED energy and force biases separately.
- calculate(atoms=None, properties=['energy', 'forces', 'energy_bias', 'forces_bias'], system_changes=['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms']) None
Compute the properties and attach them to the results
- compute_energy_and_forces(pos, istep) Tuple
Compute unbiased energies and forces, and PLUMED energy and force biases separately
- implemented_properties: List[str] = ['energy', 'forces', 'energy_bias', 'forces_bias']
Properties calculator can handle (energy, forces, …)
- class mlptrain.sampling.plumed.PlumedCustomCV(filename: str, component: str | None = None, units: str | None = None)
Bases:
_PlumedCVClass used to initialise a PLUMED collective variable from a file
- __init__(filename: str, component: str | None = None, units: str | None = None)
PLUMED collective variable from a file. The file must be written in the style of a PLUMED input file, but only contain input used in the definition of the CV, and the CV itself. There are no constraints on the CVs that can be defined.
- e.g. dof0: COM ATOMS=1-5
dof1: CENTER ATOMS=3-7 cv0: CUSTOM ARG=dof0,dof1 FUNC=x^2*y^2 PERIODIC=NO
- Parameters:
component – (str) Name of a component of the last CV in the supplied PLUMED input file to use as a collective variable
units – (str) Units of the collective variable, used in plots
- class mlptrain.sampling.plumed.PlumedDifferenceCV(name: str, atom_groups: Sequence | None = None)
Bases:
_PlumedCVClass used to initialise a PLUMED collective variable as a difference between two degrees of freedom
- __init__(name: str, atom_groups: Sequence | None = None)
PLUMED collective variable as a difference between two degrees of freedom (distances, angles, torsions),
e.g. [(0, 2), (0, 1)] gives ζ = r_02 - r_01
which are used to generate DOFs
- mlptrain.sampling.plumed.get_colvar_filename(cv: _PlumedCV, **kwargs) str
Return the name of the file where the trajectory in terms of collective variable values would be written
- mlptrain.sampling.plumed.get_hills_filename(**kwargs) str
Return the name of the file where a list of deposited gaussians would be written
- mlptrain.sampling.plumed.plot_cv1_and_cv2(filenames: Sequence[str], style: str = 'scatter', cvs_units: Sequence[str] | None = None, cvs_limits: Sequence[Sequence[float]] | None = None, label: str | None = None) None
Plot the trajectory of the system by tracking two collective variables using two colvar files. The function only works for two collective variables.
‘histogram’
- Parameters:
cvs_units – Units of the CVs to be plotted
cvs_limits – Min and max limits of the CVs in the plot
label – (str) Label attached to the name of the plot, useful when multiple plots of the same CVs are generated in the same directory
- mlptrain.sampling.plumed.plot_cv_versus_time(filename: str, style: str = 'trajectory', time_units: str = 'ps', cv_units: str | None = None, cv_limits: Sequence[float] | None = None, label: str | None = None) None
Plot a collective variable as a function of time from a given colvar file. Only plot the first collective variable in the colvar file.
- Parameters:
time_units – (str) Units of time
cv_units – Units of the CV to be plotted
cv_limits – Min and max limits of the CV in the plot
label – (str) Label attached to the name of the plot, useful when multiple plots of the same CVs are generated in the same directory
- mlptrain.sampling.plumed.plumed_setup(bias: PlumedBias, temp: float, interval: int, **kwargs) List[str]
Generate a list which represents the PLUMED input file
- Parameters:
interval – (int) Interval between saving the geometry