Public API¶
- class sire.morph.Perturbation(mol, map=None)[source]¶
This class provides lots of convenience functions that make it easier to work with an visualise perturbations
- extract_perturbed(properties: list[str] = None, remove_ghosts: bool = True)[source]¶
Extract the perturbed state properties of the molecule and remove the reference state.
If a list of properties is passed then only those properties will be extracted to the reference molecule
- Parameters:
properties (list[str]) – The list of properties to extract to the perturbed molecule
remove_ghosts (bool) – If True then any ghost atoms will be removed from the molecule
- extract_reference(properties: list[str] = None, remove_ghosts: bool = True)[source]¶
Extract the reference state properties of the molecule and remove the perturbed state.
If a list of properties is passed then only those properties will be extracted to the reference molecule
- Parameters:
properties (list[str]) – The list of properties to extract to the reference molecule
remove_ghosts (bool) – If True then any ghost atoms will be removed from the molecule
- link_to_perturbed(properties: list[str] = None, auto_commit: bool = True)[source]¶
Link all of the properties of the molecule to their values in the perturbed molecule (lambda=1).
If a list of properties is passed then only those properties will be linked to the perturbed molecule
- Parameters:
properties (list[str]) – The list of properties to link to the perturbed molecule
auto_commit (bool) – If True then the molecule will be committed and returned
- link_to_reference(properties: list[str] = None, auto_commit: bool = True)[source]¶
Link all of the properties of the molecule to their values in the reference molecule (lambda=0).
If a list of properties is passed then only those properties will be linked to the reference molecule
- Parameters:
properties (list[str]) – The list of properties to link to the reference molecule
auto_commit (bool) – If True then the molecule will be committed and returned
- to_openmm(constraint: str = None, perturbable_constraint: str = None, swap_end_states: bool = None, include_constrained_energies: bool = None, map: dict = None)[source]¶
Return the PerturbableOpenMMMolecule that corresponds to this perturbation in the OpenMM simulation. The arguments to this function have the same meaning as the equivalents in the dynamics() and minimisation() functions of a molecule.
- Parameters:
constraint (str) – The constraint algorithm to use
perturbable_constraint (str) – The constraint algorithm to use for perturbable atoms
swap_end_states (bool) – If True then the end states will be swapped
include_constrained_energies (bool) – If True then the constrained energies will be included
map (dict) – The property map to use
- Returns:
The perturbable OpenMM molecule
- Return type:
- zero_ghost_angles(auto_commit: bool = True)[source]¶
Zero the angles in the reference and perturbed states where any of the atoms in those states is a ghost atom
- sire.morph.annihilate(mol, as_new_molecule: bool = True, map=None)[source]¶
Return a merged molecule that represents the perturbation that completely annihilates the molecule. The returned merged molecule will be suitable for using in a double-annihilation free energy simulation, e.g. to calculate absolute binding free energies. Note that this perturbation will remove all intramolecular interactions, not just nonbonded intramolecular interactions. You should add positional restraints to all atoms in the molecule to prevent to prevent it drifting apart.
- Parameters:
mol (Molecule view) – The molecule (or part of molecule) to annihilate. This will only annihilate the atoms in this molecule view. Normally, you would want to pass in the entire molecule.
as_new_molecule (bool, optional) – Whether to return the merged molecule as a new molecule, or to assign a new molecule number to the result. Default is True.
map (dict, optional) – Property map to assign properties in the returned, merged molecule, plus to find the properties that will be annihilated.
- Returns:
The merged molecule representing the annihilation perturbation
- Return type:
- sire.morph.create_from_pertfile(mol, pertfile, map=None)[source]¶
Create a merged molecule from the passed molecule and pertfile. This will create and return a merged molecule with an initial and reference state that follows the instructions in the pertfile.
- Parameters:
mol (sire.mol.Molecule) – The molecule to be merged
pertfile (str) – The path to the pertfile
- Returns:
The merged molecule
- Return type:
- sire.morph.decouple(mol, as_new_molecule: bool = True, map=None)[source]¶
Return a merged molecule that represents the perturbation that completely decouples the molecule. The returned merged molecule will be suitable for using in a double-decoupling free energy simulation, e.g. to calculate absolute binding free energies.
- Parameters:
mol (Molecule view) – The molecule (or part of molecule) to decouple. This will only decouple the atoms in this molecule view. Normally, you would want to pass in the entire molecule.
as_new_molecule (bool, optional) – Whether to return the merged molecule as a new molecule, or to assign a new molecule number to the result. Default is True.
map (dict, optional) – Property map to assign properties in the returned, merged molecule, plus to find the properties that will be decoupled.
- Returns:
The merged molecule representing the decoupling perturbation
- Return type:
- sire.morph.evaluate_xml_force(mols, xml, force)[source]¶
Evaluate the custom force defined in the passed XML file. The passed molecules must be the ones used to create the OpenMM context associated with the XML file.
- Parameters:
mols (sire.system.System, sire.mol.Molecule) – The perturbable molecular system or molecule to evaluate the force on. This should have already been linked to the appropriate end state.
xml (str) – The path to the XML file containing the custom force.
force (str) – The name of the custom force to evaluate. Options are: “ghost-ghost”, “ghost-nonghost”, “ghost-14”.
- Returns:
pairs ([(sire.mol.Atom, sire.mol.Atom)]) – The atom pairs that interacted.
nrg_coul ([sr.units.GeneralUnit]) – The Coulomb energies for each atom pair.
nrg_lj ([sr.units.GeneralUnit]) – The Lennard-Jones energies for each atom pair.
- sire.morph.extract_perturbed(mols, remove_ghosts: bool = True, map=None)[source]¶
Return the passed molecule(s) where any perturbable molecules have been extracted to their perturbed (lambda=1) state (i.e. deleting their reference state)
- sire.morph.extract_reference(mols, remove_ghosts: bool = True, map=None)[source]¶
Return the passed molecule(s) where any perturbable molecules have been extracted to their reference (lambda=0) state (i.e. deleting their perturbed state)
- sire.morph.link_to_perturbed(mols, map=None)[source]¶
Return the passed molecule(s), where any perturbable molecules are linked to their perturbed (lambda=1) states
- sire.morph.link_to_reference(mols, map=None)[source]¶
Return the passed molecule(s), where any perturbable molecules are linked to their reference (lambda=0) states
- sire.morph.match(mol0, mol1, match=None, prematch=None, match_light_atoms=False, map0=None, map1=None)¶
Perform a simple match that tries to identify the mapping from atoms in ‘mol0’ to the atoms in ‘mol1’. This uses the AtomMCSMatcher to match the atoms, using the passed prematch argument.
However, if the
match
argument is provided, this will be used as the atom mapping directly (it can either be a dictionary mapping atom identifiers, or an AtomMatcher object).- Parameters:
mol0 (Molecule view) – The reference state molecule (or part of molecule)
mol1 (Molecule view) – The perturbed state molecule (or part of molecule)
match (dict, AtomMatcher, optional) – The atom matcher to use to match atoms. If this is a dictionary of atom identifiers, then this will be passed to a AtomIDMatcher object. If this is an AtomMatcher object, then this will be used directly.
prematch (dict, AtomMatcher, optional) – The atom matcher to use to prematch atoms. If
match
is not supplied, then this will be used as the prematch argument to the AtomMCSMatcher used to find the maximum common subgraph match.match_light_atoms (bool, optional) – Whether to match light atoms (i.e. hydrogen atoms) if using the default AtomMCSMatcher. Default is False.
map0 (dict, optional) – Property map to find properties in mol0
map1 (dict, optional) – Property map to find properties in mol1
- Returns:
The atom mapping between the two molecules (or parts of molecules)
- Return type:
- sire.morph.merge(mol0, mol1, match=None, prematch=None, map=None, map0=None, map1=None)[source]¶
Merge together the atoms in ‘mol0’ and ‘mol1’ and return as a single merged (perturbable) molecule. This function will conduct a merge and return a perturbable molecule such that it is equivalent to mol0 at the reference state (lambda=0) and equivalent to mol1 at the perturbed state (lambda=1).
The sr.morph.match_atoms function will be called with the passed
match
and prematch arguments to determine the atom mapping between the two molecules.- Parameters:
mol0 (Molecule view) – The reference state molecule (or part of molecule)
mol1 (Molecule view) – The perturbed state molecule (or part of molecule)
match (dict, AtomMapping, optional) – If provided, this will be passed as the
match
argument to sr.morph.match_atoms, to aid in the atom mapping.prematch (dict, AtomMapping, optional) – If provided, this will be passed as the prematch argument to sr.morph.match_atoms, to aid in the atom mapping.
map (dict, optional) – Property map to assign properties in the returned, merged molecule.
map0 (dict, optional) – Property map to find properties in mol0
map1 (dict, optional) – Property map to find properties in mol1
- Returns:
The merged molecule
- Return type:
- sire.morph.mutate(mol0, mol1, match=None, prematch=None, map=None, map0=None, map1=None)[source]¶
Mutate mol0 to mol1, returning the mutated (new) molecule. This is equivalent to calling
merge
on the two molecules (or parts of molecules) and then extracting the perturbed state.This function is most useful for mutating parts of molecules, e.g. passing in two residues as mol0 and mol1 would mutate that residue to the other within the larger molecule containing mol0. This can be used for mutating residues in proteins, or for copying and pasting parts of one molecule into another.
- Parameters:
mol0 (Molecule view) – The molecule (or part of molecule) that will be mutated.
mol1 (Molecule view) – The molecule (or part of molecule) that will be mutated to. This will replace the atoms in mol0.
match (dict, AtomMatcher, optional) – If provided, this will be passed as the
match
argument to sr.morph.match_atoms, to aid in the atom mapping.prematch (dict, AtomMatcher, optional) – If provided, this will be passed as the prematch argument to sr.morph.match_atoms, to aid in the atom mapping.
map (dict, optional) – Property map to assign properties in the returned, mutated molecule.
map0 (dict, optional) – Property map to find properties in mol0
map1 (dict, optional) – Property map to find properties in mol1
- Returns:
The mutated molecule
- Return type:
- sire.morph.repartition_hydrogen_masses(mols, mass_factor=1.5, ignore_water: bool = False, ignore_non_water: bool = False, include_end_states: bool = True, map=None)[source]¶
Increase the mass of hydrogen atoms to hmass times * amu, and subtract the mass increase from the heavy atom the hydrogen is bonded to.
- (based heavily on the repartitionMasses function in
Sire.Tools.OpenMMMD)
- Parameters:
mol (sire.mol.Molecule or list of molecules, or System) – The molecule(s) whose hydrogen masses should be repartitioned
mass_factor (float) – The factor to multiply the mass of hydrogen atoms by. Using the default of 1.5 is known to be a good value to use to achieve a 4 fs timestep with the (default) LangevinMiddle integrator
ignore_water (bool) – Whether or not to ignore water molecules (default False)
ignore_non_water (bool) – Whether or not to ignore non-water molecules (default False)
include_end_states (bool) – Whether or not to repartition the masses of the end states of perturbable molecules (default True) (i.e. this will automatically repartition mass0 and mass1 in addition to mass)
map (dict) – The property map used to identify molecular properties
- Returns:
The repartitioned molecule(s)
- Return type:
sire.mol.Molecule, list of molecules or System
- sire.morph.replica_exchange(replica0, replica1, rangenerator=None, rescale_velocities: bool = True)[source]¶
Perform a replica exchange move between the passed two replicas. This will be a Hamiltonian Replica Exchange, where the energies at the two λ-values of the replicas will be evaluated for both, and a Monte Carlo test applied to the difference of the difference of those energies to decide if the exchange should be accepted.
- Parameters:
replica0 (sire.mol.Dynamics) – The dynamics object holding the first replica to exchange
replica1 (sire.mol.Dynamics) – The dynamics object holding the second replica to exchange
rangenerator (sire.maths.RanGenerator, optional) – A random number generator to use for the Monte Carlo test. This should return uniformly distributed random numbers between 0 and 1. If this is not specified, then the global random number generator will be used.
rescale_velocities (bool, optional) – If True, then the velocities of the replicas will be rescaled to the new temperature after the exchange.
- Returns:
The replicas, either swapped if the move was accepted, or in the same order if the move was rejected. Plus, a boolean to say if the move was or was not accepted.
- Return type:
(replica0, replica1, bool)
Examples
>>> from sire.morph import replica_exchange >>> replica0 = mols.dynamics(lambda_value=0.0, ...) >>> replica1 = mols.dynamics(lambda_value=0.2, ...) >>> replica0, replica1, accepted = replica_exchange(replica0, replica1)
- sire.morph.shrink_ghost_atoms(mols, length=None, map=None)[source]¶
Update all of the molecules (or single molecule) in ‘mols’ so that the bond lengths involving ghost atoms at one or other end state (but not both) are shrunk to length. This is 0.6 A by default (this seems to be a value that causes fewest NaNs).
- sire.morph.to_alchemlyb(energy_trajectories, temperature=None, energy_unit='kcal/mol')[source]¶
Convert the passed list of energy trajectories into a single alchemlyb-compatible DataFrame, ready for free energy analysis via alchemlyb.
- Parameters:
energy_trajectories (list of sire.maths.EnergyTrajectory objects,) –
- list of filenames of s3 files,
globbed expression for list of filenames etc.
A list of EnergyTrajectory objects, each containing the energy trajectory for a single simulation at a single lambda value.
temperature (temperature, optional) – The temperature of the simulation. If not provided, the temperature will be taken from the values in each EnergyTrajectory
energy_unit (str) – Whichever of the alchemlyb energy units you want the output DataFrame to use. This is in alchemlyb format, e.g. kcal/mol, kJ/mol, or kT
- Returns:
A single DataFrame containing the energy trajectories from all simulations, ready for free energy analysis via alchemlyb.
- Return type:
pandas.DataFrame
- sire.morph.zero_ghost_angles(mols, map=None)[source]¶
Zero any angles in the reference or perturbed states where any of the atoms in those states is a ghost (dummy) atom