Side Chain Ensembles

RotamerEnsemble

The RotamerEnsemble is the base class for (monofunctional) side chain ensemble objects. It contains all the methods and attributes that are shared between standard canonical amino acids and specialized and non-canonical amino acids.

class RotamerEnsemble(res, site=None, protein=None, chain=None, rotlib=None, **kwargs)

Create new RotamerEnsemble object.

Parameters:
  • res (string) – 3-character name of desired residue, e.g. R1A.

  • site (int) – MolSys residue number to attach library to.

  • protein (MDAnalysis.Universe, MDAnalysis.AtomGroup) – Object containing all protein information (coords, atom types, etc.)

  • chain (str) – MolSys chain identifier to attach spin label to.

  • rotlib (str) – Rotamer library to use for constructing the RotamerEnsemble

  • **kwargs (dict) –

    minimize: bool

    Switch to turn on/off minimization. During minimization each rotamer is optimized in dihedral space with respect to both internal and external clashes, and deviation from the starting conformation in dihedral space.

    min_method: str

    Name of the minimization algorithm to use. All scipy.optimize.minimize algorithms are available and include: ‘Nelder-Mead’, ‘Powell’, ‘CG’, ‘BFGS’, ‘Newton-CG’, ‘L-BFGS-B’, ‘TNC’, ‘COBYLA’, ‘SLSQP’, ‘trust-constr’, ‘dogleg’, ‘trust-ncg’, ‘trust-exact’, ‘trust-krylov’, and custom.

    exclude_nb_interactions: int:

    When calculating internal clashes, ignore 1-exclude_nb_interactions interactions and below. Defaults to ignore 1-3 interactions, i.e. atoms that are connected by 2 bonds or fewer will not have a steric effect on each other.

    eval_clashbool

    Switch to turn clash evaluation on (True) and off (False).

    forcefield: str

    Name of the forcefield you wish to use to parameterize atoms for the energy function. Currently, supports charmm and uff

    energy_funccallable

    Python function or callable object that takes a protein and a RotamerEnsemble object as input and returns an energy value (kcal/mol) for each atom of each rotamer in the ensemble. See also Scoring . Defaults to chiLife.get_lj_rep .

    forgivefloat

    Softening factor to be passed to energy_func. Only used if energy_func uses a softening factor. Defaults to 1.0. See Scoring .

    tempfloat

    Temperature to use when running energy_func. Only used if energy_func accepts a temperature argument Defaults to 298 K.

    clash_radiusfloat

    Cutoff distance (angstroms) for inclusion of atoms in clash evaluations. This distance is measured from clash_ori Defaults to the longest distance between any two atoms in the rotamer ensemble plus 5 angstroms.

    clash_oristr

    Atom selection to use as the origin when finding atoms within the clash_radius. Defaults to ‘cen’, the centroid of the rotamer ensemble heavy atoms.

    protein_treeScipy.spatial.cKDTree

    KDTree of atom positions for fast distance calculations and neighbor detection. Defaults to None

    trim: bool

    When true, the lowest trim_tol fraction of rotamers in the ensemble will be removed.

    trim_tol: float

    Tolerance for trimming rotamers from the ensemble. trim_tol=0.005 means the bottom 0.5% of rotamers will be removed.

    alignment_methodstr

    Method to use when attaching or aligning the rotamer ensemble backbone with the protein backbone. Defaults to 'bisect' which aligns the CA atom, the vectors bisecting the N-CA-C angle and the N-CA-C plane.

    sampleint, bool

    Argument to use the off-rotamer sampling method. If False or 0 the off-rotamer sampling method will not be used. If int the ensemble will be generated with that many off-rotamer samples.

    dihedral_sigmasfloat, numpy.ndarray

    Standard deviations of dihedral angles (degrees) for off rotamer sampling. Can be a single number for isotropic sampling, a vector to define each dihedral individually or a matrix to define a value for each rotamer and each dihedral. Setting this value to np.inf will force uniform (accessible volume) sampling. Defaults to 35 degrees.

    weighted_samplingbool

    Determines whether the rotamer ensemble is sampled uniformly or based off of their intrinsic weights. Defaults to False.

    use_Hbool

    Determines if hydrogen atoms are used or not. Defaults to False.

classmethod from_pdb(pdb_file, res=None, site=None, protein=None, chain=None)
Parameters:
  • pdb_file

  • res – (Default value = None)

  • site – (Default value = None)

  • protein – (Default value = None)

  • chain – (Default value = None)

classmethod from_mda(residue, **kwargs)

Create RotamerEnsemble from the MDAnalysis.Residue for chiLife.Residue

Parameters:
  • residue (MDAnalysis.Residue, chiLife.Residue) – A Residue object from which to create the RotamerEnsemble.

  • **kwargs (dict) – Additional keyword arguments to use for creating the RotamerEnsemble.

Returns:

cls – A RotamerEnsemble object created from the provided Residue object.

Return type:

RotamerEnsemble

classmethod from_trajectory(traj, site, chain=None, energy=None, burn_in=0, **kwargs)

Create a RotamerEnsemble object from a trajectory.

Parameters:
  • traj (MDAnalysis.Universe, MDAnalysis.AtomGroup, chiLife.MolecularSystemBase) – The trajectory over which the RotamerEnsemble is to be created.

  • site (res) – Residue number of the residue to create a rotamer ensemble for.

  • chain (str) – Chain identifier of the residue to create a RotamerEnsemble for.

  • energy (ArrayLike) – Potential energy of the residue at each frame of the trajectory in kcal/mol

  • burn_in (int) – Number of frames to skip from the trajectory. Used to skip burn-in period when performing MCMC sampling.

  • **kwargs (dict) –

    Additional keyword arguments to use for creating the RotamerEnsemble.

    dihedral_atomslist

    A list of atoms defining the mobile dihedral angles of the residue. Primarily used when creating an ensemble of a residue type chiLife is unaware of (see Defining Mobile Dihedrals

    tempfloat

    Temperature to use when evaluating conformational populations from user supplied energy argument.

    spin_atomslist, dict

    A list or dict defining the spin atoms of the label and their relative populations (see Defining Spin-atoms and Their Weights).

Returns:

cls – A RotamerEnsemble object created from the provided trajectory.

Return type:

RotamerEnsemble

to_rotlib(libname: str | None = None, description: str | None = None, comment: str | None = None, reference: str | None = None) None

Saves the current RotamerEnsemble as a rotamer library file.

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

  • description (str) – Description of the rotamer library.

  • comment (str) – Additional comments about the rotamer library.

  • reference (str) – References associated with the rotamer library.

set_site(site)

Assign RotamerLibrary to a site.

Parameters:

site (int) – Residue number to assign SpinLabel to.

copy(site: int | None = None, chain: str | None = None, rotlib: dict | None = None)

Create a deep copy of the spin label. Assign new site, chain or rotlib if desired. Useful when labeling homo-oligomers or proteins with repeating units/domains.

Parameters:
  • site (int) – New site number for the copy.

  • chain (str) – New chain identifier for the copy.

  • rotlib (dict) – New (base) rotamer library for the dict. Used primarily when copying dRotamerEnsembels.

Returns:

new_copy – A deep copy of the original RotamerLibrary

Return type:

RotamerEnsemble

to_site(site_pos: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] | None = None) None

Move spin label to new site

Parameters:

site_pos (ArrayLike) – 3x3 array of ordered backbone atom coordinates of new site (N CA C) (Default value = None)

ICs_to_site(cori, cmx)

Modify the internal coordinates to be aligned with the site that the RotamerEnsemble is attached to

backbone_to_site()

Modify backbone atoms to match the site that the RotamerEnsemble is attached to

sample(n=1, off_rotamer=False, **kwargs)

Randomly sample a rotamer in the library.

Parameters:
  • n (int) – number of rotamers to sample.

  • off_rotamer (bool | List[bool]) – If False the rotamer library used to construct the RotamerEnsemble will be sampled using the exact dihedrals. If True, all mobile dihedrals will undergo minor perturbations. A list indicates a per-dihedral switch to turn on/off off-rotamer sampling of a particular dihedrals, corresponding to self.dihedral_atoms e.g. off_+rotamer = [False, False, False, True, True] for R1 will only sample off-rotamers for χ4 and χ5.

  • **kwargs (dict) –

    return_dihedralsbool

    If True, sample will return a MolSysIC object of the sampled rotamer

Returns:

  • coords (ArrayLike) – The 3D coordinates of the sampled rotamer(s)

  • new_weights (ArrayLike) – New weights (relative) of the sampled rotamer(s)

  • ICs (List[chilife.MolSysIC] (Optional)) – Internal coordinate (MolSysIC) objects of the rotamer(s).

update_weight(weight: float) None

Function to assign self.current_weight, which is the estimated weight of the rotamer currently occupying the attachment site on self.protein. This is only relevant if the residue type on self.protein is the same as the RotamerLibrary.

Parameters:

weight (float) – New weight for the current residue.

minimize(callback=None)

Minimize the energy of the rotamers in dihedral space and re-weight with respect to the new energies.

Parameters:

callback (callable) – A callable function to be passed as the scipy.optimize.minimize function. See the scipy documentation for details.

trim_rotamers(keep_idx: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] | None = None) None

Remove rotamers with small weights from ensemble

Parameters:

keep_idx (ArrayLike[int]) – Indices of rotamers to keep. If None rotamers will be trimmed based off of self.trim_tol.

property centroid

Get the centroid of the whole rotamer ensemble.

evaluate()

Place rotamer ensemble on protein site and recalculate rotamer weights.

save_pdb(name: str | None = None) None

Save a pdb file containing the rotamer ensemble.

Parameters:

name (str | None) – Name of the file. If None self.name will be used.

get_lib(rotlib)

Parse backbone information from protein and fetch the appropriate rotamer library.

Parameters:

rotlib (str | dict) – The name of the rotlib or a dictionary containing all the information required for a rotamer library.

Returns:

lib – Dictionary containing all underlying attributes of a rotamer library.

Return type:

dict

property bonds

Set of atom indices corresponding to the atoms that form covalent bonds

property non_bonded

Set of atom indices corresponding to the atoms that are not covalently bonded. Also excludes atoms that have 1-n non-bonded interactions where n=self._exclude_nb_interactions . By default, 1-3 interactions are excluded

property clash_ori

Location to use as ‘center’ when searching for atoms with potential clashes

property coords

Cartesian coordinates of all rotamers in the ensemble.

property dihedrals

Values of the (mobile) dihedral angles of all rotamers in the ensemble

property mx

The rotation matrix to rotate a residue from the local coordinate frame to the current residue. The local coordinate frame is defined by the alignment method

property origin

Origin of the local coordinate frame

property CB

The coordinates of the β-carbon atom of the RotamerEnsemble

intra_fit()

Align all rotamers of the ensemble to the ensemble backbone

get_sasa()

Calculate the solvent accessible surface area (SASA) of each rotamer in the protein environment.

set_dihedral_sampling_sigmas(value)

Helper function to assign new sigmas for off-rotamer sampling.

Parameters:

value (float | ArrayLike[float]) – The sigma value(s) to assign. A single float will assign all dihedrals of all rotamers to the same sigma (isotropic). An array with as many elements as dihedrals will assign different sigmas to each dihedral (Anisotropic). An array of dimensions n_rots by n_dihedrals will assign different anisotropic sigmas to different rotamers.

SpinLabel

The chilife.SpinLabel object is the primary feature of χLife. chilife.SpinLabel inherits from the chilife.RotamerEnsemble object and therefore has all the same properties and methods . Additionally chilife.SpinLabel have several other features unique to spin labels and useful for protein and spin label modeling.

class SpinLabel(label, site=None, protein=None, chain=None, rotlib=None, **kwargs)

A RotamerEnsemble made from a side chain with one or more unpaired electrons.

Parameters:
  • label (str) – Name of desired spin label, e.g. R1A.

  • site (int) – MolSys residue number to attach spin label to.

  • chain (str) – MolSys chain identifier to attach spin label to.

  • protein (MDAnalysis.Universe, MDAnalysis.AtomGroup) – Object containing all protein information (coords, atom types, etc.)

  • protein_tree (cKDtree) – k-dimensional tree object associated with the protein coordinates.

property spin_coords

get the spin coordinates of the rotamer ensemble

property spin_centers

get the spin center of the rotamers in the ensemble

property spin_centroid

Average location of all the label’s spin_coords weighted based off of the rotamer weights

classmethod from_mmm(label, site=None, protein=None, chain=None, **kwargs)

Create a SpinLabel object using the default MMM protocol with any modifications passed via kwargs

classmethod from_wizard(label, site=None, protein=None, chain=None, to_find=200, to_try=10000, vdw=2.5, clashes=5, **kwargs)

Create a SpinLabel object using the default MTSSLWizard protocol with any modifications passed via kwargs

dRotamerEnsemble

The dRotamerEnsemble is the base class for bifunctional side chain ensemble objects. Like chilife.RotamerEnsemble, it contains all the methods and attributes that are shared between all bifunctional amino acids whether they are spin labels or other bifunctional non-canonical amino acids.

class dRotamerEnsemble(res, sites, protein=None, chain=None, rotlib=None, **kwargs)

Create new dRotamerEnsemble object.

Parameters:
  • res (string) – 3-character name of desired residue, e.g. RXA.

  • site (tuple[int, int]) – MolSys residue numbers to attach the bifunctional library to.

  • protein (MDAnalysis.Universe, MDAnalysis.AtomGroup) – Object containing all protein information (coords, atom types, etc.)

  • chain (str) – MolSys chain identifier to attach the bifunctional ensemble to.

  • rotlib (str) – Rotamer library to use for constructing the RotamerEnsemble

  • **kwargs (dict) –

    restraint_weight: float

    Force constant (kcal/mol/A^2) for calculating energetic penalty of restraint satisfaction, i.e. the alignment of the overlapping atoms of the two mono-functional subunits of the bifunctional label.

    torsion_weight: float

    Force constant (kcal/mol/radian^2) for calculating energetic penalty of the deviation from rotamer starting dihedral angles.

    minimize: bool

    Switch to turn on/off minimization. During minimization each rotamer is optimized in dihedral space with respect to alignment of the “cap” atoms of the two mono-functional subunits, internal clashes and deviation from the starting conformation in dihedral space.

    min_method: str

    Name of the minimization algorithm to use. All scipy.optimize.minimize algorithms are available and include: ‘Nelder-Mead’, ‘Powell’, ‘CG’, ‘BFGS’, ‘Newton-CG’, ‘L-BFGS-B’, ‘TNC’, ‘COBYLA’, ‘SLSQP’, ‘trust-constr’, ‘dogleg’, ‘trust-ncg’, ‘trust-exact’, ‘trust-krylov’, and custom.

    exclude_nb_interactions: int:

    When calculating internal clashes, ignore 1-exclude_nb_interactions interactions and below. Defaults to ignore 1-3 interactions, i.e. atoms that are connected by 2 bonds or fewer will not have a steric effect on each other.

    eval_clashbool

    Switch to turn clash evaluation on (True) and off (False).

    forcefield: str

    Name of the forcefield you wish to use to parameterize atoms for the energy function. Currently, supports charmm and uff

    energy_funccallable

    Python function or callable object that takes a protein and a RotamerEnsemble object as input and returns an energy value (kcal/mol) for each atom of each rotamer in the ensemble. See also Scoring . Defaults to chiLife.get_lj_energy .

    forgivefloat

    Softening factor to be passed to energy_func. Only used if energy_func uses a softening factor. Defaults to 0.95. See Scoring .

    tempfloat

    Temperature to use when running energy_func. Only used if energy_func accepts a temperature argument Defaults to 298 K.

    clash_radiusfloat

    Cutoff distance (angstroms) for inclusion of atoms in clash evaluations. This distance is measured from clash_ori Defaults to the longest distance between any two atoms in the rotamer ensemble plus 5 angstroms.

    clash_oristr

    Atom selection to use as the origin when finding atoms within the clash_radius. Defaults to ‘cen’, the centroid of the rotamer ensemble heavy atoms.

    protein_treeScipy.spatial.cKDTree

    KDTree of atom positions for fast distance calculations and neighbor detection. Defaults to None

    trim: bool

    When true, the lowest trim_tol fraction of rotamers in the ensemble will be removed.

    trim_tol: float

    Tolerance for trimming rotamers from the ensemble. trim_tol=0.005 means the bottom 0.5% of rotamers will be removed.

    alignment_methodstr

    Method to use when attaching or aligning the rotamer ensemble backbone with the protein backbone. Defaults to 'bisect' which aligns the CA atom, the vectors bisecting the N-CA-C angle and the N-CA-C plane.

    dihedral_sigmasfloat, numpy.ndarray

    Standard deviations of dihedral angles (degrees) for off rotamer sampling. Can be a single number for isotropic sampling, a vector to define each dihedral individually or a matrix to define a value for each rotamer and each dihedral. Setting this value to np.inf will force uniform (accessible volume) sampling. Defaults to 35 degrees.

    weighted_samplingbool

    Determines whether the rotamer ensemble is sampled uniformly or based off of their intrinsic weights. Defaults to False.

    use_Hbool

    Determines if hydrogen atoms are used or not. Defaults to False.

res

3-character name of desired residue, e.g. RXA.

Type:

string

site1

MolSys residue number of the first attachment site.

Type:

int

site2

MolSys residue number of the second attachment site.

Type:

int

increment

Number of residues between two attachment sites.

Type:

int

protein

Object containing all protein information (coords, atom types, etc.)

Type:

MDAnalysis.Universe, MDAnalysis.AtomGroup, MolSys

chain

Chain identifier of the site the bifunctional ensemble is attached to.

Type:

str

name

Name of the ensemble. Usually include the native amino acid, the site number and the label that was attached changing the name will change the object name when saving in a PDB.

Type:

str

RE1

Monofunctional ensemble subunit attached to the first site.

Type:

RotamerEnsemble

RE2

Monofunctional ensemble subunit attached to the second site.

Type:

RotamerEnsemble

property weights

Array of the fractions of rotamer populations for each rotamer in the library.

property coords

The 3D cartesian coordinates of each atom of each rotamer in the library.

property atom_names

The names of each atom in the rotamer

property atom_types

The element or atom type of each atom in the rotamer.

property dihedral_atoms

Four atom sets defining each flexible dihedral of the side chain

property dihedrals

attr::~dihedral_atoms for each rotamer in the library

Type:

Dihedral angle values of each dihedral defined in

Type:

py

property centroid

The geometric center of all atoms of all rotamers in the rotamer library

property clash_ori

The origin used to determine if an external atom will be considered for clashes using the clash_radius property of the ensemble

property side_chain_idx

Indices of the atoms that correspond to the side chain atoms (e.g. CB, CG, etc. and not N, CA, C)

property bonds

Array of intra-label atom pairs indices that are covalently bonded.

property non_bonded

Array of indices of intra-label non-bonded atom pairs. Primarily used for internal clash evaluation when sampling the dihedral space

guess_chain()

Function to guess the chain based off of the attached protein and residue number.

Returns:

chain – The chain name

Return type:

str

get_lib(rotlib)

Specil function to get a rotamer library for a dRotamerEnsemble and apply all attributes to the object instance. rotlib can be a residue name, a file name, or a dictionary containing all the standard entries of a rotamer library.

Parameters:

rotlib (Union[Path, str, dict]) – The name of the residue, the Path to the rotamer library file or a dictionary containing all entries of a rotamer library

create_ensembles()

Creates monofunctional components of the bifunctional rotamer ensemble.

save_pdb(name: str | None = None)

Save a PDB file of the ensemble

Parameters:

name (str) – Name of the PDB file

minimize(callback=None)

Minimize rotamers in dihedral space in the current context. Performed by default unless the minimize=False keyword argument is used during construction. Note that the minimization method is controlled by the

Parameters:

callback (Callable) –

A callable function to be passed as the scipy.optimize.minimize function. See the scipy documentation for details.

trim_rotamers()

Remove low probability rotamers from the ensemble. All rotamers accounting for the population less than self.trim_tol will be removed.

evaluate()

Place rotamer ensemble on protein site and recalculate rotamer weights.

copy()

Create a deep copy of the dRotamerEnsemble object

dSpinLabel

class dSpinLabel(label, sites, protein=None, chain=None, rotlib=None, **kwargs)

The dSpinLabel constructor has all the same arguments and keyword arguments as the dRotamerEnsemble class in addition to a few attributes relating to the unpaired electron.

SL1

Monofunctional ensemble subunit attached to the first site.

Type:

SpinLabel

SL2

Monofunctional ensemble subunit attached to the first site.

Type:

SpinLabel

property spin_atoms

Names of the atoms where the unpaired electron density is localized.

property spin_idx

Indices of the atoms where the unpaired electron density is localized.

property spin_coords

coordinates of the atoms where the unpaired electron density is localized for each rotamer.

property spin_centers

Weighted average location of the unpaired electron density for each rotamer.

property spin_weights

Relative unpaired electron density of the spin_atoms.

property spin_centroid

Weighted average location of the spin coordinates for reach rotamer of the ensemble.

create_ensembles()

Overrides parent create_ensemble method to create monofunctional components of the bifunctional rotamer ensemble. Major difference is that it creates self.SL1/2 properties in addition to self.RE1/2

IntrinsicLabel

class IntrinsicLabel(res, atom_selection, spin_atoms=None, name='IntrinsicLabel')

A helper object to assign spin density to part of a protein that is already present, e.g. metal ligands. IntrinsicLabels can be used with the distance_distribution() function

Parameters:
  • res (str) – 3-character identifier of desired residue, e.g. NC2 (Native Cu(II)).

  • atom_selection (MDAnalysis.AtomGroup, chiLife.System) – Group of atoms constituting the intrinsic label. These selections can consist of multiple residues, which can be useful in the case of ions with multiple coordinating residues

  • spin_atoms (str, list, tuple, array, dict, MDAnalysis.AtomGroup) – Atoms of the intrinsic label that host the unpaired electron density. Can be a single atom name, a list/tuple/array of atom names or a dictionary mapping atom names to their relative populations. spin_atoms `` can also be an ``MDAnalysis.AtomGroup object derived from the same MDAnalysis.Universe as the atom_selection keyword argument. Use of an AtomGroup is particularly useful when there is spin density is distributed on several atoms with the same name on different residues within the IntrinsicLabel.

property spin_coords

get the spin coordinates of the rotamer ensemble

property spin_centers

get the spin center of the rotamers in the ensemble

property spin_centroid

Average location of all the label’s spin_coords weighted based off of the rotamer weights