ves package
Submodules
ves.basis module
Classes defining basis expansions.
- class ves.basis.LegendreBasis1D(degree, min, max, axis='x', weights=None)
Bases:
torch.nn.modules.module.Module
Legendre polynomial basis expansion along x- or y- coordinate.
The coordinates x (or y) are scaled to lie within [-1, 1] using the [min, max] attributes, as
\[x' = \frac{x - (min + max) / 2}{(max - min) / 2}\]A legendre basis expansion is defined over \(x'\), as
\[B(x') = \sum_{i=0}^{d} w_i P_i(x')\]where \(P_i\) is the legendre polynomial of degree \(i\), \(w_i\) is its coefficient in the expansion, and \(d\) is the degree of the expansion.
Notes
Weights \(w_i\) are learnable as self.weights is a torch Parameter.
This basis expansion can directly be used a TorchForce module.
- degree
Degree of expansion.
- Type
int
- min
Min x-/y-value for scaling.
- Type
float
- max
Max x-/y-value for scaling.
- Type
float
- axis
‘x’ or ‘y’ (default=’x’).
- Type
str
- weights
Legendre polynomial coefficients (array len = degree).
- Type
torch.nn.Parameter
- forward(positions)
The forward method returns the energy computed from positions.
- Parameters
positions – torch.Tensor with shape (1, 3) positions[0, k] is the position (in nanometers) of spatial dimension k of particle 0
- Returns
- torch.Scalar
The potential energy (in kJ/mol)
- Return type
potential
- classmethod legendre_polynomial(x, degree: int) torch.Tensor
Computes a legendre polynomial of degree \(n\) using dynamic programming and the Bonnet’s recursion formula:
\[(n + 1) P_{n+1}(x) = (2n + 1) x P_n(x) - nP_{n-1}(x)\]
- training: bool
- class ves.basis.LegendreBasis2DPathCV(degree, x_i, y_i, gamma_i, lam, weights=None, kappa=0)
Bases:
torch.nn.modules.module.Module
Legendre polynomial basis expansion along a path collective variable (CV) defined using the points \(x_i, y_i\), which represent a path in 2D space.
For a point \((x, y)\) the CVs \(s\) and \(z\), which measure the location of the point parallel to the path and perpendicular to the path respectively, are computed as
\[s = \frac{1}{N - 1} \frac{\sum_{i=0}^{N-1} i\ e^{-\lambda [(x - x_i) ^ 2 + (y - y_i) ^ 2]}}{\sum_{i=0}^{N-1} e^{-\lambda [(x - x_i) ^ 2 + (y - y_i) ^ 2]}}\]\[z = -\frac{1}{\lambda} \ln (\sum_{i=0}^{N-1} e^{-\lambda [(x - x_i) ^ 2 + (y - y_i) ^ 2]})\]The CV \(s\) is scaled to lie within [-1, 1] as
\[s' = (s - 1/2) / (1/2)\]A legendre basis expansion is defined over \(s'\), as
\[B(s') = \sum_{i=0}^{d} w_i P_i(s')\]where \(P_i\) is the legendre polynomial of order \(i\), \(w_i\) is its coefficient in the expansion, and \(d\) is the degree of the expansion.
A restraining harmonic bias is defined to restrict sampling to regions near the path, as
\[U(z) = \frac{\kappa}{2} z^2\]This basis expansion can directly be used a TorchForce module.
- degree
Degree of basis.
- Type
int
- x_i
x-coordinates of images along the path.
- Type
torch.DoubleTensor
- y_i
y-coordinates of images along the path.
- Type
torch.DoubleTensor
- gamma_i
- Type
torch.DoubleTensor
- lam
Value of \(\lambda\).
- Type
float
- weights
Legendre polynomial coefficients of expansion along \(s\) (array len = degree).
- Type
torch.nn.Parameter
- kappa
Strength of restraining harmonic potential along \(z\).
- Type
torch.DoubleTensor
- forward(positions)
The forward method returns the energy computed from positions.
- Parameters
positions – torch.Tensor with shape (1, 3) positions[0, k] is the position (in nanometers) of spatial dimension k of particle 0
- Returns
- torch.Scalar
The potential energy (in kJ/mol)
- Return type
potential
- classmethod legendre_polynomial(x, degree: int) torch.Tensor
Computes a legendre polynomial of degree \(n\) using dynamic programming and the Bonnet’s recursion formula:
\[(n + 1) P_{n+1}(x) = (2n + 1) x P_n(x) - nP_{n-1}(x)\]
- training: bool
- class ves.basis.LegendreBasis2DRadialCV(degree, x_min, y_min, x_max, y_max, weights=None)
Bases:
torch.nn.modules.module.Module
Legendre polynomial basis expansion along a radial collective variable (CV) defined in the neighborhood of points \((x_min, y_min)\) and \((x_max, y_max)\).
A point \((x, y)\) is mapped to a CV \(s\) as
\[s = \frac{(x - x_{min})^2 + (y - y_{min})^2}{(x_{max} - x_{min})^2 + (y_{max} - y_{min})^2}\]and then scaled to lie within [-1, 1] as
\[s' = (s - 1/2) / (1/2)\]A legendre basis expansion is defined over \(s'\), as
\[B(s') = \sum_{i=0}^{d} w_i P_i(s')\]where \(P_i\) is the legendre polynomial of order \(i\), \(w_i\) is its coefficient in the expansion, and \(d\) is the degree of the expansion.
Notes
Weights \(w_i\) are learnable as self.weights is a torch Parameter.
This basis expansion can directly be used a TorchForce module.
- degree
Degree of basis.
- Type
int
- x_min
Start x-coordinate.
- Type
float
- y_min
Start y-coordinate.
- Type
float
- x_max
End x-coordinate.
- Type
float
- y_max
End y-coordinate.
- Type
float
- weights
Legendre polynomial coefficients (array len = degree).
- Type
torch.nn.Parameter
- forward(positions)
The forward method returns the energy computed from positions.
- Parameters
positions – torch.Tensor with shape (1, 3) positions[0, k] is the position (in nanometers) of spatial dimension k of particle 0
- Returns
- torch.Scalar
The potential energy (in kJ/mol)
- Return type
potential
- classmethod legendre_polynomial(x, degree: int) torch.Tensor
Computes a legendre polynomial of degree \(n\) using dynamic programming and the Bonnet’s recursion formula:
\[(n + 1) P_{n+1}(x) = (2n + 1) x P_n(x) - nP_{n-1}(x)\]
- training: bool
- class ves.basis.NNBasis1D(min, max, axis='x', hidden_layer_sizes=[], act=<class 'torch.nn.modules.activation.ReLU'>)
Bases:
torch.nn.modules.module.Module
Neural network basis expansion along x- or y- coordinate.
The coordinates x (or y) are scaled to lie within [0, 1] using the [min, max] attributes, as
\[x' = (x - min) / (max - min)\]A nonlinear basis expansion is defined over \(x'\) as
\[B(x') = N(x')\]where \(N\) is a neural network.
Notes
Neural network weights are learnable.
This basis expansion can directly be used a TorchForce module.
- min
Min x-/y-value for scaling.
- Type
float
- max
Max x-/y-value for scaling.
- Type
float
- axis
‘x’ or ‘y’ (default=’x’).
- Type
str
- layers
List of neural layers.
- Type
torch.nn.ModuleList
- forward(positions)
The forward method returns the energy computed from positions.
- Parameters
positions – torch.Tensor with shape (1, 3) positions[0, k] is the position (in nanometers) of spatial dimension k of particle 0
- Returns
- torch.Scalar
The potential energy (in kJ/mol)
- Return type
potential
- training: bool
ves.bias module
Classes defining bias potentials and update mechanisms.
- class ves.bias.HarmonicBias_SingleParticle_1D(k, s0, axis='x', model_loc='model.pt')
Bases:
ves.bias.NNBias
Applies a harmonic bias potential \(k/2 (s - s0)^2\) along the x- or y- coordinate of a single particle.
- k
Strength of harmonic bias.
- Type
float
- s0
Location (x-/y- value) of harmonic bias’ center.
- Type
float
- axis
Axis (x/y) along which bias is applied.
- Type
str
- model_loc
Path to save module to.
- class ves.bias.HarmonicBias_SingleParticle_1D_ForceModule(k, s0, axis='x')
Bases:
torch.nn.modules.module.Module
A harmonic potential \(k/2 (s - s0)^2\) as a static compute graph.
The axis attribute specifies whether the potential is applied to the x- or the y- coordinate of the particle.
- k
Strength of harmonic bias.
- Type
float
- s0
Location (x-/y- value) of harmonic bias’ center.
- Type
float
- axis
Axis (x/y) along which bias is applied.
- Type
str
- forward(positions)
The forward method returns the energy computed from positions.
- Parameters
positions – torch.Tensor with shape (1, 3) positions[0, k] is the position (in nanometers) of spatial dimension k of particle 0
- Returns
- torch.Scalar
The potential energy (in kJ/mol)
- Return type
potential
- training: bool
- class ves.bias.NNBias(model_loc='model.pt')
Bases:
ves.bias.Bias
Base class defining neural network biases.
Child classes must set the self.model parameter to a TorchScript compatible neural network instance before calling super().__init__().
Warning
Failing to set the self.model parameter will throw an exception.
- model
TorchScript-compatible neural network instance which returns bias energy from particle coordinates.
- model_loc
Path to save module to.
- property force
- class ves.bias.StaticBias_SingleParticle(V_module, model_loc)
Bases:
ves.bias.NNBias
Applies a static basis set expanded potential to a single particle.
- The potential (V_module) can be an instance of:
ves.basis.LegendreBasis1D: an expansion over a basis set (Valsson and Parrinello 2014) along the x- or the y- coordinate.
ves.basis.LegendreString2D: an expansion over a basis set (Valsson and Parrinello 2014) along a string in the x- and y- coordinates.
ves.basis.NNBasis1D: a neural network (Valsson and Parrinello 2014) along the x- or the y- coordinate.
- V_module
TorchScript compatible neural network instance.
- model_loc
Path to save module to.
- class ves.bias.VESBias_SingleParticle_1D(V_module, target, beta, optimizer_type, optimizer_params, model_loc, update_steps=5000)
Bases:
ves.bias.NNBias
Applies a dynamic basis set expanded potential along the x- or y- coordinate of a single particle. The parameters of the basis set can be updated by calling the update function with a trajectory.
- The potential (V_module) can be an instance of:
ves.basis.LegendreBasis1D: an expansion over a basis set (Valsson and Parrinello 2014).
ves.basis.NNBasis1D: a neural network (Valsson and Parrinello 2014).
The axis attribute of the V_module instance specifies whether the potential acts along the x- or y- axis.
- V_module
- target
- beta
- optimizer_type
- optimizer_params
- model_loc
- update_steps
- update(traj)
ves.config_creation module
Functions to help create initial configurations.
- ves.config_creation.singleParticle1D_init_coord(potential: ves.potentials.Potential1D, temp: float, ntrials: int = 100, xmin: float = 0, xmax: float = 1)
Uses monte carlo trials to place a particle on a free energy surface defined by potential.
Attempts are made within bounds (xmin, xmax).
- ves.config_creation.singleParticle2D_init_coord(potential: ves.potentials.Potential2D, temp: float, ntrials: int = 100, xmin: float = 0, xmax: float = 1, ymin: float = 0, ymax: float = 1)
Uses monte carlo trials to place a particle on a free energy surface defined by potential.
Attempts are made within bounds (xmin, xmax) and (ymin, ymax).
ves.langevin_dynamics module
Classes for performing biased langevin dynamics simulations using OpenMM.
- class ves.langevin_dynamics.SingleParticleSimulation(potential: openmm.openmm.CustomExternalForce, mass: int = 1, temp: float = 300, friction: float = 100, timestep: float = 10, init_state: Optional[openmm.openmm.State] = None, init_coord: numpy.ndarray = array([[0, 0, 0]]), gpu: bool = False, cpu_threads: Optional[int] = None, seed: Optional[int] = None, traj_in_mem: bool = True)
Bases:
object
Performs langevin dynamics simulation of a particle on a 2D potential energy surface with an optional VES bias.
- init_ves(bias: ves.bias.Bias, static: bool = False, startafter: int = 2, learnevery: int = 50)
Initialize simulation with VES bias.
ves.potentials module
Classes defining potential energy surfaces.
- class ves.potentials.CatchBondPotential2D(force_x=0, force_y=0, y_0=1, y_scale=5, y_shift=4, xy_0=0, xy_scale=2, gx_0=2, gx_scale=0.5, gy_0=- 2.5, gy_scale=0.25)
Bases:
ves.potentials.Potential2D
3-basin catch bond potential (slip bond potential with an extra harmonic basin).
\[U(x, y) = -\ln \left[ e^{-\left( \left(\frac{(y - y_0)^2}{y_{scale}} - y_{shift} \right)^2 + \frac{(x - y - xy_{shift})^2}{2} \right) } + e^-{(x - gx_0)^2 / gx_scale + (y - gy_0)^2 / gy_scale} \right]\]- potential(x, y)
Computes the catch bond potential at a given point (x, y).
- class ves.potentials.DoubleWellPotential1D
Bases:
ves.potentials.Potential1D
1D double well potential.
- potential(x)
Computes the potential at a given point x.
- Parameters
x (float) – Point to compute potential at.
- Returns
Value of potential at x.
- Return type
V (float)
- class ves.potentials.DoubleWellPotential2D
Bases:
ves.potentials.Potential2D
Double well potential.
\[U(x, y) = 10(x^4 - 4 x^2 + 0.7 x) + 5 y^2\]- potential(x, y)
Computes the double well potential at a given point (x, y).
- class ves.potentials.MullerBrownPotential
Bases:
ves.potentials.Potential2D
2D Muller-Brown potential.
- A = [-200, -100, -170, 15]
- a = [-1, -1, -6.5, 0.7]
- b = [0, 0, 11, 0.6]
- c = [-10, -10, -6.5, 0.7]
- potential(x, y)
Compute the Muller-Brown potential at a given point (x, y).
- x_bar = [1, 0, -0.5, -1]
- y_bar = [0, 0.5, 1.5, 1]
- class ves.potentials.Potential1D
Bases:
openmm.openmm.CustomExternalForce
Abstract class defining basic 1D potential behavior.
A harmonic restraining potential of magnitude 1000 kJ/mol is applied on the y and z coordinates about y=0 and z=0.
Note
Child classes must call super.__init__() only after initializing the force attribute in the x variable.
- force
OpenMM-compatible custom force expression.
- Type
str
- potential(x: float)
Computes the potential at a given point x.
- Parameters
x (float) – Point to compute potential at.
- Returns
Value of potential at x.
- Return type
V (float)
- class ves.potentials.Potential2D
Bases:
openmm.openmm.CustomExternalForce
Abstract class defining basic 2D potential behavior.
A harmonic restraining potential of magnitude 1000 kJ/mol is applied on the z coordinates about z=0.
Note
Child classes must call super.__init__() only after initializing the force attribute in x and y variables.
- force
OpenMM-compatible custom force expression.
- Type
str
- potential(x: float, y: float)
Computes the potential at a given point (x, y).
- Parameters
x (float) – x-coordinate of the point to compute potential at.
y (float) – y-coordinate of the point to compute potential at.
- Returns
Value of potential at (x, y).
- Return type
V (float)
- class ves.potentials.SingleWellPotential1D
Bases:
ves.potentials.Potential1D
1D single well potential.
- potential(x)
Computes the potential at a given point x.
- Parameters
x (float) – Point to compute potential at.
- Returns
Value of potential at x.
- Return type
V (float)
- class ves.potentials.SingleWellPotential2D
Bases:
ves.potentials.Potential2D
Single well potential.
\[U(x, y) = x^2 + y^2\]- potential(x, y)
Computes the single well potential at a given point (x, y).
- class ves.potentials.SlipBondPotential2D(force_x=0, force_y=0, y_0=1, y_scale=5, y_shift=4, xy_0=0, xy_scale=2)
Bases:
ves.potentials.Potential2D
2-basin slip bond potential.
\[U(x, y) = \left( \left(\frac{(y - y\_0)^2}{y\_scale} - y\_shift \right)^2 + \frac{(x - y - xy\_0)^2}{xy\_scale} \right)\]- potential(x, y)
Computes the slip bond potential at a given point (x, y).
ves.target module
Class defining VES target free energy landscape.
- class ves.target.Target_Uniform_HardSwitch_x(mesh)
Bases:
ves.target.Target_x
Distribution is uniform over interval [-1, 1]. Outside of this, it is zero.
ves.utils module
Miscellanious utility classes and functions.
- class ves.utils.TrajectoryReader(traj_file, comment_char='#', format='txyz')
Bases:
object
Utility class for reading large trajectories.
- Parameters
traj_file (str) – Path to trajectory file.
comment_char (str) – Character marking the beginning of a comment line.
format (str) – Format of each line (options = ‘txyz’ or ‘xyz’; default = ‘txyz’)
- read_traj(skip=1)
Reads trajectory.
- Parameters
skip (int) – Number of frames to skip between reads (default = 1).
- Returns
tuple(T, traj) if self.format == ‘txyz’ traj if self.format == ‘xyz’
ves.visualization module
Classes for visualizing trajectory data.
- class ves.visualization.VisualizePotential1D(potential1D: ves.potentials.Potential1D, temp: float, xrange: tuple, mesh: int = 200)
Bases:
object
Class defining functions to generate scatter plots and animated trajectories of a particle on a 1D potential.
- Parameters
potential1D –
temp –
xrange –
mesh –
- animate_traj(traj, outdir, every=1, s=3, c='black', call_ffmpeg: bool = True)
Plots positions at timesteps defined by interval every on potential energy surface and stitches together plots using ffmpeg to make a movie.
- Parameters
traj (iterable) –
outdir (str) –
every (int) –
s (int) –
c (str) –
call_ffmpeg (bool) –
- plot_potential()
Plots the potential within (xrange[0], xrange[1]).
- scatter_traj(traj, outimg, every=1, s=1, c='black')
Scatters entire trajectory onto potential energy surface.
- Parameters
traj –
outimg –
every –
s –
c –
- class ves.visualization.VisualizePotential2D(potential2D: ves.potentials.Potential2D, temp: float, xrange: tuple, yrange: tuple, contourvals=None, clip=None, mesh: int = 200, cmap: str = 'jet')
Bases:
object
Class defining functions to generate scatter plots and animated trajectories of a particle on a 2D potential surface.
- Parameters
potential2D (pyib.md.potentials.Potential2D) – 2D potential energy surface.
temp (float) – Temperature (required, as free energies are plotted in kT).
xrange (tuple of length 2) – Range of x-values to plot.
yrange (tuple of length 2) – Range of y-values to plot.
contourvals (int or array-like) – Determines the number and positions of the contour lines / regions. Refer to the matplotlib documentation for details.
clip (float) – Value of free energy (in kT) to clip contour plot at.
mesh – Number of mesh points in each dimension for contour plot.
cmap – Matplotlib colormap.
- animate_traj(traj, outdir, every=1, s=3, c='black', call_ffmpeg: bool = True)
Plots positions at timesteps defined by interval every on potential energy surface and stitches together plots using ffmpeg to make a movie.
- Parameters
traj (iterable) –
outdir (str) –
every (int) –
s (int) –
c (str) –
call_ffmpeg (bool) –
- animate_traj_projection_x(traj, outdir, every=1, s=3, c='black', call_ffmpeg: bool = True)
Plots positions at timesteps defined by interval every on the x-projection of the potential energy surface and stitches together plots using ffmpeg to make a movie.
- Parameters
traj (iterable) –
outdir (str) –
every (int) –
s (int) –
c (str) –
call_ffmpeg (bool) –
- animate_traj_projection_y(traj, outdir, every=1, s=3, c='black', call_ffmpeg: bool = True)
Plots positions at timesteps defined by interval every on the x-projection of the potential energy surface and stitches together plots using ffmpeg to make a movie.
- Parameters
traj (iterable) –
outdir (str) –
every (int) –
s (int) –
c (str) –
call_ffmpeg (bool) –
- plot_potential()
Plots the potential within (xrange[0], xrange[1]) and (yrange[0], yrange[1]).
- plot_projection_x()
Plots the x-projection of potential within (xrange[0], xrange[1]) and (yrange[0], yrange[1]).
- plot_projection_y()
Plots the y-projection of potential within (xrange[0], xrange[1]) and (yrange[0], yrange[1]).
- scatter_traj(traj, outimg, every=1, s=1, c='black')
Scatters entire trajectory onto potential energy surface.
- Parameters
traj –
outimg –
every –
s –
c –
- scatter_traj_projection_x(traj, outimg, every=1, s=1, c='black')
Scatters x-projection of entire trajectory onto potential energy surface.
- Parameters
traj –
outimg –
every –
s –
c –
- scatter_traj_projection_y(traj, outimg, every=1, s=1, c='black')
Scatters y-projection of entire trajectory onto potential energy surface.
- Parameters
traj –
outimg –
every –
s –
c –
- ves.visualization.visualize_free_energy_2D(xvals, yvals, xrange, yrange, nbins_x=100, nbins_y=100, contourvals=None, clip=None, cmap='jet', dpi=150)
Plots 2D free energy profile from 2D trajectory data. :param xvals: Array of x coordinates of points to bin. :type xvals: numpy.ndarray :param yvals: Array of y coordinates of points to bin. :type yvals: numpy.ndarray :param xrange: Range of x-values to plot. :type xrange: tuple of length 2 :param yrange: Range of y-values to plot. :type yrange: tuple of length 2 :param nbins_x: Number of bins along the x-axis (default=100). :type nbins_x: int :param nbins_y: Number of bins along the y-axis (default=100). :type nbins_y: int :param contourvals: Determines the number and positions of the contour lines / regions. Refer to the matplotlib documentation for details. :type contourvals: int or array-like :param clip: Value of free energy (in kT) to clip contour plot at. :type clip: float :param cmap: Matplotlib colormap (default=jet). :param dpi: Output DPI (default=150).