Yank¶
Interface for automated free energy calculations.
-
class
yank.yank.
Topography
(topology, ligand_atoms=None, solvent_atoms='auto')[source]¶ A class mapping and labelling the different components of a system.
The object holds the topology of a system and offers convenience functions to identify its various parts such as solvent, receptor, ions and ligand atoms.
A molecule should be labeled as a ligand, only if there is also a receptor. If there is only a single molecule its atom indices can be obtained from solute_atoms instead. In ligand-receptor system, solute_atoms provides the atom indices for both molecules.
Parameters: topology : mdtraj.Topology or simtk.openmm.app.Topology
The topology object specifying the system.
ligand_atoms : iterable of int or str, optional
The atom indices of the ligand. A string is interpreted as an mdtraj DSL specification of the ligand atoms.
solvent_atoms : iterable of int or str, optional
The atom indices of the solvent. A string is interpreted as an mdtraj DSL specification of the solvent atoms. If ‘auto’, a list of common solvent residue names will be used to automatically detect solvent atoms (default is ‘auto’).
Attributes
ligand_atoms
The atom indices of the ligand as list receptor_atoms
The atom indices of the receptor as list (read-only). solute_atoms
The atom indices of the non-solvent molecule(s) (read-only). solvent_atoms
The atom indices of the solvent molecules. ions_atoms
The indices of all ions atoms in the solvent (read-only). -
topology
¶ mdtraj.Topology: A copy of the topology (read-only).
-
ligand_atoms
¶ The atom indices of the ligand as list
This can be empty if this
Topography
doesn’t represent a receptor-ligand system. Use solute_atoms to obtain the atom indices of the molecule if this is the case.If assigned to a string, it will be interpreted as an mdtraj DSL specification of the atom indices.
-
receptor_atoms
¶ The atom indices of the receptor as list (read-only).
This can be empty if this Topography doesn’t represent a receptor-ligand system. Use solute_atoms to obtain the atom indices of the molecule if this is the case.
-
solute_atoms
¶ The atom indices of the non-solvent molecule(s) (read-only).
Practically, this are all the indices of the atoms that are not considered solvent. In a receptor-ligand system, this includes the atom indices of both the receptor and the ligand.
-
solvent_atoms
¶ The atom indices of the solvent molecules.
This includes eventual ions. If assigned to a string, it will be interpreted as an mdtraj DSL specification of the atom indices. If assigned to ‘auto’, a set of solvent auto indices is automatically built from common solvent residue names.
-
ions_atoms
¶ The indices of all ions atoms in the solvent (read-only).
-
-
class
yank.yank.
IMultiStateSampler
[source]¶ A sampler for multiple thermodynamic states.
This is the interface documents the properties and methods that need to be exposed by the sampler object to be compatible with the class
AlchemicalPhase
.Attributes
number_of_iterations
int: the total number of iterations to run. iteration
int: the current iteration. metadata
dict: a copy of the metadata dictionary passed on creation. sampler_states
list of SamplerState: the sampler states at the current iteration. -
number_of_iterations
¶ int: the total number of iterations to run.
-
iteration
¶ int: the current iteration.
-
metadata
¶ dict: a copy of the metadata dictionary passed on creation.
-
sampler_states
¶ list of SamplerState: the sampler states at the current iteration.
-
create
(thermodynamic_state, sampler_states, storage, unsampled_thermodynamic_states, metadata)[source]¶ Create new simulation and initialize the storage.
Parameters: thermodynamic_state : list of openmmtools.states.ThermodynamicState
The thermodynamic states for the simulation.
sampler_states : openmmtools.states.SamplerState or list
One or more sets of initial sampler states. If a list of SamplerStates, they will be assigned to thermodynamic states in a round-robin fashion.
storage : str or Reporter
If str: The path to the storage file. Reads defaults from the
yank.repex.Reporter
classIf
yank.repex.Reporter
: Reads the reporter settings for files and optionsIn the future this will be able to take a Storage class as well.
unsampled_thermodynamic_states : list of openmmtools.states.ThermodynamicState
These are ThermodynamicStates that are not propagated, but their reduced potential is computed at each iteration for each replica. These energy can be used as data for reweighting schemes.
metadata : dict
Simulation metadata to be stored in the file.
-
minimize
(tolerance, max_iterations)[source]¶ Minimize all states.
Parameters: tolerance : simtk.unit.Quantity
Minimization tolerance (units of energy/mole/length, default is
1.0 * unit.kilojoules_per_mole / unit.nanometers
).max_iterations : int
Maximum number of iterations for minimization. If 0, minimization continues until converged.
-
equilibrate
(n_iterations, mcmc_moves=None)[source]¶ Equilibrate all states.
Parameters: n_iterations : int
Number of equilibration iterations.
mcmc_moves : MCMCMove or list of MCMCMove, optional
Optionally, the MCMCMoves to use for equilibration can be different from the ones used in production (default is None).
-
run
(n_iterations=None)[source]¶ Run the simulation.
This runs at most
number_of_iterations
iterations. Useextend()
to pass the limit.Parameters: n_iterations : int, optional
If specified, only at most the specified number of iterations will be run (default is None).
-
-
class
yank.yank.
AlchemicalPhase
(sampler)[source]¶ A single thermodynamic leg (phase) of an alchemical free energy calculation.
This class wraps around a general MultiStateSampler and handle the creation of an alchemical free energy calculation.
Parameters: sampler : MultiStateSampler
The sampler instance implementing the
IMultiStateSampler
interface.Attributes
iteration
int: the current iteration (read-only). number_of_iterations
int: the total number of iterations to run. is_completed
Boolean check if if the sampler has been completed by its own determination or if we have exceeded number of -
static
from_storage
(storage)[source]¶ Static constructor from an existing storage file.
Parameters: storage : str or Reporter
If str: The path to the primary storage file. Default checkpointing options are stored in this case
If
yank.repex.Reporter
: loads from the reporter class, including checkpointing informationIn the future this will be able to take a Storage class as well.
Returns: alchemical_phase : AlchemicalPhase
A new instance of
AlchemicalPhase
in the same state of the last stored iteration.
-
iteration
¶ int: the current iteration (read-only).
-
number_of_iterations
¶ int: the total number of iterations to run.
-
is_completed
¶ Boolean check if if the sampler has been completed by its own determination or if we have exceeded number of iterations
-
create
(thermodynamic_state, sampler_states, topography, protocol, storage, restraint=None, anisotropic_dispersion_cutoff=None, alchemical_regions=None, alchemical_factory=None, metadata=None)[source]¶ Create a new AlchemicalPhase calculation for a specified protocol.
If
anisotropic_dispersion_cutoff
is different thanNone
. The end states of the phase will be reweighted. The fully interacting state accounts for:- The truncation of nonbonded interactions.
2. The reciprocal space which is not modeled in alchemical states if an Ewald method is used for long-range interactions.
Parameters: thermodynamic_state : openmmtools.states.ThermodynamicState
Thermodynamic state holding the reference system, temperature and pressure.
sampler_states : openmmtools.states.SamplerState or list
One or more sets of initial sampler states. If a list of SamplerStates, they will be assigned to replicas in a round-robin fashion.
topography : Topography
The object holding the topology and labelling the different components of the system. This is used to discriminate between ligand-receptor and solvation systems.
protocol : dict
The dictionary
{parameter_name: list_of_parameter_values}
defining the protocol. All the parameter values list must have the same number of elements.storage : str or initialized Reporter class
If str: Path to the storage file. The default checkpointing options (see the
yank.repex.Reporter
class) will be used in this caseIf
yank.repex.Reporter
: Uses files and checkpointing options of the reporter class passed inrestraint : ReceptorLigandRestraint, optional
Restraint to add between protein and ligand. This must be specified for ligand-receptor systems in non-periodic boxes.
anisotropic_dispersion_cutoff : simtk.openmm.Quantity, ‘auto’, or None, optional, default None
If specified, this is the cutoff at which to reweight long range interactions of the end states to correct for anisotropic dispersions.
If ‘auto’, then the distance is automatically chosen based on the minimum possible size it can be given the box volume, then behaves as if a Quantity was passed in.
If None, the correction won’t be applied (units of length, default is None).
alchemical_regions : openmmtools.alchemy.AlchemicalRegion, optional
If specified, this is the
AlchemicalRegion
that will be passed to theAbsoluteAlchemicalFactory
, otherwise the ligand will be alchemically modified according to the given protocol.alchemical_factory : openmmtools.alchemy.AbsoluteAlchemicalFactory, optional
If specified, this
AbsoluteAlchemicalFactory
will be used instead of the one created with default options.metadata : dict, optional
Simulation metadata to be stored in the file.
-
minimize
(tolerance=Quantity(value=1.0, unit=kilojoule/(nanometer*mole)), max_iterations=0)[source]¶ Minimize all the states.
The minimization is performed in two steps. In the first one, the positions are minimized in the reference thermodynamic state (i.e. non alchemically-modified). Only then, the positions are minimized in their alchemically softened state.
Parameters: tolerance : simtk.unit.Quantity, optional
Minimization tolerance (units of energy/mole/length, default is
1.0 * unit.kilojoules_per_mole / unit.nanometers
).max_iterations : int, optional
Maximum number of iterations for minimization. If 0, minimization continues until converged.
-
randomize_ligand
(sigma_multiplier=2.0, close_cutoff=Quantity(value=1.5, unit=angstrom))[source]¶ Randomize the ligand positions in every state.
The position and orientation of the ligand in each state will be randomized. This works only if the system is a ligand-receptor system.
If you call this before minimizing, each positions will be minimized separately in the reference state, so you may want to call it afterwards to speed up minimization.
Parameters: sigma_multiplier : float, optional
The ligand will be placed close to a random receptor atom at a distance that is normally distributed with standard deviation
sigma_multiplier * receptor_radius_of_gyration
(default is 2.0).close_cutoff : simtk.unit.Quantity, optional
Each random placement proposal will be rejected if the ligand ends up being closer to the receptor than this cutoff (units of length, default is
1.5*unit.angstrom
).
-
equilibrate
(n_iterations, mcmc_moves=None)[source]¶ Equilibrate all states.
Parameters: n_iterations : int
Number of equilibration iterations.
mcmc_moves : MCMCMove or list of MCMCMove, optional
Optionally, the MCMCMoves to use for equilibration can be different from the ones used in production.
-
static