chiLife Functions

General Functions

distance_distribution(*args: SpinLabel, r: ArrayLike | None = None, sigma: float = 1.0, use_spin_centers: bool = True, dependent: bool = False, uq: bool = False) ndarray

Calculates total distribution of spin-spin distances among an arbitrary number of spin labels, using the distance range r (in angstrom). The distance distribution is obtained by summing over all label pair distance distributions. These in turn are calculated by convolving a weighted histogram of pairwise distances between lables with a normal distribution. Distance between labels are calculated either between label spin_center if spin_populations=False, or between all pairs of spn-bearing atoms if spin_populations=True.

Parameters:
  • *args (SpinLabel) – Any number of spin label objects.

  • r (ArrayLike) – Evenly spaced array of distances (in angstrom) to calculate the distance distribution over.

  • sigma (float) – The standard deviation of the normal distribution used in convolution with the distance histogram, in angstrom. Default is 1.

  • use_spin_centers (bool) – If False, distances are computed between the weighted centroids of the spin atoms. If True, distances are computed by summing over the distributed spin density on spin-bearing atoms on the labels.

  • dependent (bool) – Consider the (clash) effects of spin label rotamers on each other.

  • uq (bool) – Perform uncertainty analysis (experimental)

Returns:

P – Predicted distance distribution, in 1/angstrom

Return type:

np.ndarray

create_library(libname: str, pdb: str, dihedral_atoms: List[List[str]] | None = None, site: int = 1, resname: str | None = None, dihedrals: ArrayLike | None = None, weights: ArrayLike | None = None, sigmas: ArrayLike | None = None, permanent: bool = False, default: bool = False, force: bool = False, spin_atoms: List[str] | None = None, aln_atoms: List[str] | None = None, backbone_atoms: List[str] | None = None, description: str | None = None, comment: str | None = None, reference: str | None = None) None

Create a chiLife rotamer library from a pdb file and a user defined set of parameters.

Parameters:
  • libname (str) – Name for the new rotamer library.

  • pdb (str) – Name of (and path to) pdb file containing the user defined spin label structure. This pdb file should contain only the desired spin label and no additional residues.

  • site (int) – The residue number of the side chain in the pdb file you would like to add.

  • resname (str) – 3-letter residue code.

  • dihedral_atoms (list) –

    a list of rotatable dihedrals. List should contain sub-lists of 4 atom names. Atom names must be the same as defined in the pdb file eg:

    [['CA', 'CB', 'SG', 'SD'],
    ['CB', 'SG', 'SD', 'CE']...]
    

  • dihedrals (ArrayLike, optional) – Array of dihedral angles. If provided the new label object will be stored as a rotamer library with the dihedrals provided. Array should be n x m where n is the number of rotamers and m is the number of dihedrals.

  • weights (ArrayLike, optional) – Weights associated with the dihedral angles provided by the dihedrals keyword argument. There should be one weight per rotamer and the rotamer weights should sum to 1.

  • sigmas (ArrayLike, optional) – Sigma parameter for distributions of dihedral angles. Should be n x m matrix where n is the number of rotamers and m is the number of dihedrals. This feature will be used when performing off rotamer samplings.

  • permanent (bool) – If set to True the library will be stored in the chilife user_rotlibs directory in addition to the current working directory.

  • default (bool) – If set to true and permanent is also set to true then this rotamer library will become the default rotamer library for the given resname

  • force (bool = False,) – If set to True and permanent is also set to true this library will overwrite any existing library with the same name.

  • spin_atoms (list) – A list of atom names on which the spin density is localized.

  • aln_atoms (List[str]) – The three atoms defining the alignment coordinate frame of the target residue. Usually the N-CA-C atoms of the protein backbone or similar for nucleotide backbones.

  • backbone_atoms (List[str]) – Names of the atoms that correspond to backbone atoms. This includes aln_atoms and any other atoms that may need to be adjusted to fit the target residue backbone, e.g. the hydrogen atoms of carboxyl oxygen atoms. Any atoms not part of the backbone will be considered side chain atoms and may be subject to sampling.

  • description (str) – A short (< 60 characters) description of the side chain library being created

  • comment (str) – Additional information about the rotamer library that may not fit in description.

  • reference (str) – Any relevant citations associated with the rotamer library.

create_dlibrary(libname: str, pdb: str, sites: Tuple, dihedral_atoms: List[List[List[str]]], resname: str | None = None, dihedrals: ArrayLike | None = None, weights: ArrayLike | None = None, sort_atoms: List[str] | Dict | str = True, permanent: bool = False, default: bool = False, force: bool = False, spin_atoms: List[str] | None = None, description: str | None = None, comment: str | None = None, reference: str | None = None) None

Create a chiLife bifunctional rotamer library from a pdb file and a user defined set of parameters.

Parameters:
  • libname (str,) – Name for the user defined label.

  • pdb (str) – Name of (and path to) pdb file containing the user defined spin label structure. This pdb file should contain only the desired spin label and no additional residues.

  • sites (Tuple[int]) – The residue numbers that the bifunctional label is attached to.

  • dihedral_atoms (List[List[List[str]]]) –

    A list defining the mobile dihedral atoms of both halves of the bifunctional label. There should be two entries, one for each site that the bifunctional label attaches to. Each entry should contain a sublist of dihedral definitions. Dihedral definitions are a list of four strings that are the names of the four atoms defining the dihedral angle.

    [
    # Side chain 1
    [['CA', 'CB', 'SG', 'SD'],
    ['CB', 'SG', 'SD', 'CE']...],
    
    # Side chain 2
    [['CA', 'CB', 'SG', 'SD'],
    ['CB', 'SG', 'SD', 'CE']...]
    ]
    

  • resname (str) – Residue type 3-letter code.

  • dihedrals (ArrayLike, optional) – Array of dihedral angles. If provided the new label object will be stored as a rotamer library with the dihedrals provided.

  • weights (ArrayLike, optional) – Weights associated with the dihedral angles provided by the dihedrals keyword argument permanent: bool If set to True the library will be stored in the chilife user_rotlibs directory in addition to the current working directory.

  • sort_atoms (bool) – Rotamer library atoms are sorted to achieve consistent internal coordinates by default. Sometimes this sorting is not desired or the PDB file is already sorted. Setting this option to False turns off atom sorting.

  • permanent (bool) – If set to True, a copy of the rotamer library will be stored in the chiLife rotamer libraries directory

  • default (bool) – If set to true and permanent is also set to true then this rotamer library will become the default rotamer library for the given resname

  • force (bool = False,) – If set to True and permanent is also set to true this library will overwrite any existing library with the same name.

  • spin_atoms (list, dict) – List of atom names on which the spin density is localized, e.g [‘N’, ‘O’], or dictionary with spin atom names (key ‘spin_atoms’) and spin atom weights (key ‘spin_weights’).

  • description (str) – A short (< 60 characters) description of the side chain library being created

  • comment (str) – Additional information about the rotamer library that may not fit in description.

  • reference (str) – Any relevant citations associated with the rotamer library.

add_rotlib_dir(directory: Path | str) None

Add a directory to search for rotlibs when none are found in te working directory. This directory will be searched before the chiLife default rotlib directory.

Parameters:

directory (Path, str) – String or Path object of the desired directory.

remove_rotlib_dir(directory: Path | str) None

Remove a user added directory to search for rotlibs.

Parameters:

directory (Path, str) – String or Path object of the desired directory.

add_library(filename: str | Path, libname: str | None = None, default: bool = False, force: bool = False)

Add the provided rotamer library to the chilife rotamer library directory so that it does not need to be in the working directory when utilizing.

Note

It is generally preferred to keep a custom directory of rotlibs and use add_rotlib_dir() to force chiLife to search that directory for rotlibs since this function adds rotlibs to the chilife install directory that can be overwritten with updates and is not preserved when updating python.

Parameters:
  • filename (str, Path) – Name or Path object oof the rotamer library npz file.

  • libname (str) – Unique name of the rotamer library

  • default (bool) – If True, chilife will make this the default rotamer library for this residue type.

  • force (bool) – If True, chilife will overwrite any existing rotamer library with the same name if it exists.

remove_library(name: str, prompt: bool = True)

Removes a library from the chilife rotamer library directory and from the chilife dihedral definitions

Parameters:
  • name (str) – Name of the rotamer library

  • prompt (bool) – If set to False al warnings will be ignored.

list_available_rotlibs()

Lists residue types and rotamer libraries that are currently available. More information on any individual rotamer library by using the rotlib_info function.

repack(protein: ~MDAnalysis.core.universe.Universe | ~MDAnalysis.core.groups.AtomGroup, *spin_labels: ~chilife.RotamerEnsemble.RotamerEnsemble, repetitions: int = 200, temp: float = 1, energy_func: ~typing.Callable = <function get_lj_rep>, off_rotamer=False, **kwargs) Tuple[Universe, ArrayLike]

Markov chain Monte Carlo repack a protein around any number of SpinLabel or RotamerEnsemble objects.

Parameters:
  • protein (MDAnalysis.Universe, MDAnalysis.AtomGroup) – MolSys to be repacked

  • spin_labels (RotamerEnsemble, SpinLabel) – RotamerEnsemble or SpinLabel object placed at site of interest

  • repetitions (int) – Number of successful MC samples to perform before terminating the MC sampling loop

  • temp (float, ArrayLike) – Temperature (Kelvin) for both clash evaluation and metropolis-hastings acceptance criteria. Accepts a list or array like object of temperatures if a temperature schedule is desired

  • energy_func (Callable) – Energy function to be used for clash evaluation. Must accept a protein and RotamerEnsemble object and return an array of potentials in kcal/mol, with one energy per rotamer in the rotamer ensemble

  • off_rotamer (bool) – Boolean argument that decides whether off rotamer sampling is used when repacking the provided residues.

  • kwargs (dict) – Additional keyword arguments to be passed to mutate .

Returns:

  • protein (MDAnalysis.Universe) – MCMC trajectory of local repack

  • deltaEs (np.ndarray:) – Change in energy_func score at each accept of MCMC trajectory

IO Functions

read_distance_distribution(file_name: str) Tuple[ndarray, ndarray]

Reads a DEER distance distribution file in the DeerAnalysis or similar format.

Parameters:

file_name (str) – File name of distance distribution

Returns:

  • r (np.ndarray) – Domain of distance distribution in the same units as the file.

  • p (np.ndarray) – Probability density over r normalized such that the integral of p over r is 1.

read_library(rotlib: str, Phi: float | None = None, Psi: float | None = None) Dict

Function to read rotamer libraries and bifunctional rotamer libraries as dictionaries.

Parameters:
  • rotlib (str, Path) – 3 letter residue code or path the rotamer library file.

  • Phi (float) – Backbone Phi dihedral angle for the provided residue. Only applicable for canonical amino acids.

  • Psi (float) – Backbone Psi dihedral angle for the provided residue. Only applicable for canonical amino acids.

Returns:

lib – Dictionary of arrays containing rotamer library information in cartesian and dihedral space.

Return type:

dict

save(file_name: str, *molecules: RotamerEnsemble | MolecularSystemBase | Universe | AtomGroup | str, protein_path: str | Path | None = None, mode: str = 'w', conect=False, **kwargs) None

Save a pdb file of the provided labels and proteins

Parameters:
  • file_name (str) – Desired file name for output file. Will be automatically made based off of protein name and labels if not provided.

  • *molecules (RotmerLibrary, chiLife.MolSys, mda.Universe, mda.AtomGroup, str) – Object(s) to save. Includes RotamerEnsemble, SpinLabels, dRotamerEnsembles, proteins, or path to protein pdb. Can add as many as desired, except for path, for which only one can be given.

  • protein_path (str, Path) – The Path to a pdb file to use as the protein object.

  • mode (str) – Which mode to open the file in. Accepts ‘w’ or ‘a’ to overwrite or append.

  • **kwargs

    Additional arguments to pass to write_labels

    write_spin_centersbool

    Write spin centers (atoms named NEN) as a separate object with weights mapped to q-factor.

write_atoms(file: TextIO, atoms: Universe | AtomGroup | MolecularSystemBase, coords: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] | None = None) None

Write a set of atoms to a file in PDB format without any prefix or suffix, i.e. Not wrapped in MODEL or ENDMDL

Parameters:
  • file (TextIO) – A writable file (TextIO) object. The file that the atoms will be written to.

  • atoms (MDAnalysis.Universe, MDAnalysis.AtomGroup, MolecularSystemBase) – Object containing atoms to be written. May also work with duck-typed set of atoms that are iterable and each atom has properties index, name,

  • coords (ArrayLike) – Array of atom coordinates corresponding to atoms