ampal package

Submodules

ampal.align module

A module containing classes for aligning structure.

class ampal.align.MMCAlign(eval_fn, eval_args: list | None, polypeptide: Polypeptide)[source]

Bases: object

A alignment protocol that uses Metropolis Monte Carlo.

Notes

THIS IS CURRENTLY SUPER INEFFICIENT DUE TO THE DEEPCOPIES. I plan to improve this by aligning arrays of atoms and only recording the rotation and translation that led to that alignment.

Parameters:
  • eval_fn (Polypeptide -> float) – A function to evaluate the quality of your fit.

  • eval_args (list) – A list of static args to be used in the eval_fn, these will be unpacked into the evaluation function _i.e._ `eval_fn(polypeptide, *eval_args).

  • polypeptide (Polypeptide) – An ampal polypeptide containing the model to be aligned.

static check_move(new, old, t)[source]

Determines if a model will be accepted.

static float_f(f)[source]

Formats a float for printing to std out.

start_optimisation(rounds: int, max_angle: float, max_distance: float, temp: float = 298.15, stop_when=None, verbose=None)[source]

Starts the loop fitting protocol.

Parameters:
  • rounds (int) – The number of Monte Carlo moves to be evaluated.

  • max_angle (float) – The maximum variation in rotation that can moved per step.

  • max_distance (float) – The maximum distance the can be moved per step.

  • temp (float, optional) – Temperature used during fitting process.

  • stop_when (float, optional) – Stops fitting when energy is less than or equal to this value.

ampal.align.align_backbones(reference, mobile, stop_when=None, verbose=False)[source]

ampal.amino_acids module

ampal.amino_acids.add_amino_acid_to_json(code, description, letter='X', modified=None, force_add=False)[source]

Add an amino acid to the amino_acids.json file used to populate the amino_acid table.

Parameters:
  • code (str) – New code to be added to amino acid table.

  • description (str) – Description of the amino acid, e.g. ‘amidated terminal carboxy group’.

  • letter (str, optional) – One letter code for the amino acid. Defaults to ‘X’

  • modified (str or None, optional) – Code of modified amino acid, e.g. ‘ALA’, or None. Defaults to None

  • force_add (bool, optional) – If True, will over-write existing dictionary value for code if already in amino_acids.json. If False, then an IOError is raised if code is already in amino_acids.json.

Raises:

IOError – If code is already in amino_acids.json and force_add is False.

Return type:

None

ampal.amino_acids.get_aa_code(aa_letter)[source]

Get three-letter aa code if possible. If not, return None.

If three-letter code is None, will have to find this later from the filesystem.

Parameters:

aa_letter (str) – One-letter amino acid code.

Returns:

aa_code – Three-letter aa code.

Return type:

str, or None

ampal.amino_acids.get_aa_info(code)[source]

Get dictionary of information relating to a new amino acid code not currently in the database.

Notes

Use this function to get a dictionary that is then to be sent to the function add_amino_acid_to_json(). use to fill in rows of amino_acid table for new amino acid code.

Parameters:

code (str) – Three-letter amino acid code.

Raises:

IOError – If unable to locate the page associated with the amino acid name on the PDBE site.

Returns:

aa_dict – Keys are AminoAcidDB field names. Values are the str values for the new amino acid, scraped from the PDBE if possible. None if not found.

Return type:

dict

ampal.amino_acids.get_aa_letter(aa_code)[source]

Get one-letter version of aa_code if possible. If not, return ‘X’.

Parameters:

aa_code (str) – Three-letter amino acid code.

Returns:

aa_letter – One-letter aa code. Default value is ‘X’.

Return type:

str

ampal.ampal_warnings module

exception ampal.ampal_warnings.DependencyNotFoundWarning[source]

Bases: RuntimeWarning

exception ampal.ampal_warnings.MalformedPDBWarning[source]

Bases: RuntimeWarning

exception ampal.ampal_warnings.NoncanonicalWarning[source]

Bases: RuntimeWarning

exception ampal.ampal_warnings.NotParameterisedWarning[source]

Bases: RuntimeWarning

ampal.ampal_warnings.check_availability(program, test_func, global_settings)[source]

ampal.analyse_protein module

Contains various tools for analysing protein structure.

ampal.analyse_protein.alpha_angles(p, reference_axis, tag=True, reference_axis_name='ref_axis')[source]

Alpha angle calculated using points on the primitive of helix and axis.

Notes

The final value is None, since the angle calculation requires pairs of points along the primitive and axis. This is a generalisation of the calculation used to measure the tilt of a helix in a coiled-coil with respect to the central axis of the coiled coil.

Parameters:
  • p (ampal.Polymer) – Reference Polymer.

  • reference_axis (list(numpy.array or tuple or list)) – Length of reference_axis must equal length of the Polymer. Each element of reference_axis represents a point in R^3.

  • tag (bool, optional) – If True, tags the Chain with the reference axis coordinates and each Residue with its alpha angle. Alpha angles are stored at the Residue level, but are calculated using the CA atom.

  • reference_axis_name (str, optional) – Used to name the keys in tags at Chain and Residue level.

Returns:

alphas – The alpha angle for the Polymer at each point of its primitive, in degrees.

Return type:

list of float

Raises:

ValueError – If the Polymer and the reference_axis have unequal length.

ampal.analyse_protein.cc_to_local_params(pitch, radius, oligo)[source]

Returns local parameters for an oligomeric assembly.

Parameters:
  • pitch (float) – Pitch of assembly

  • radius (float) – Radius of assembly

  • oligo (int) – Oligomeric state of assembly

Returns:

  • pitchloc (float) – Local pitch of assembly (between 2 adjacent component helices)

  • rloc (float) – Local radius of assembly

  • alphaloc (float) – Local pitch-angle of assembly

ampal.analyse_protein.charge_series(seq, granularity=0.1)[source]

Calculates the charge for pH 1-13.

Parameters:
  • seq (str) – Sequence of amino acids.

  • granularity (float, optional) – Granularity of pH values i.e. if 0.1 pH = [1.0, 1.1, 1.2…]

ampal.analyse_protein.classify_angle_as_rotamer(angle: float) int | None[source]

Classifies a chi angle into a rotamer.

Notes

This function uses the convention 1:g+ 2:t 3:g- as from 10.1006/jmbi.1993.1170 R L Dunbrack Jr, M Karplus (1993) Backbone-dependent rotamer library for proteins. Application to side-chain prediction.

Parameters:

angle (float) – Chi angle to classify.

ampal.analyse_protein.crick_angles(p, reference_axis, tag=True, reference_axis_name='ref_axis')[source]

Returns the Crick angle for each CA atom in the Polymer.

Notes

The final value is in the returned list is None, since the angle calculation requires pairs of points on both the primitive and reference_axis.

Parameters:
  • p (ampal.Polymer) – Reference Polymer.

  • reference_axis (list(numpy.array or tuple or list)) – Length of reference_axis must equal length of the Polymer. Each element of reference_axis represents a point in R^3.

  • tag (bool, optional) – If True, tags the Polymer with the reference axis coordinates and each Residue with its Crick angle. Crick angles are stored at the Residue level, but are calculated using the CA atom.

  • reference_axis_name (str, optional) – Used to name the keys in tags at Chain and Residue level.

Returns:

cr_angles – The crick angles in degrees for each CA atom of the Polymer.

Return type:

list(float)

Raises:

ValueError – If the Polymer and the reference_axis have unequal length.

ampal.analyse_protein.flip_reference_axis_if_antiparallel(p, reference_axis, start_index=0, end_index=-1)[source]

Flips reference axis if direction opposes the direction of the Polymer.

Notes

If the angle between the vector for the Polymer and the vector for the reference_axis is > 90 degrees, then the reference axis is reversed. This is useful to run before running polymer_to_reference_axis_distances, crick_angles, or alpha_angles. For more information on the start and end indices, see chain_vector.

Parameters:
  • p (ampal.Polymer) – Reference Polymer.

  • reference_axis (list(numpy.array or tuple or list)) – Length of reference_axis must equal length of the Polymer. Each element of reference_axis represents a point in R^3.

  • start_index (int, optional) – Default is 0 (start at the N-terminus of the Polymer)

  • end_index (int, optional) – Default is -1 (start at the C-terminus of the Polymer)

Returns:

reference_axis

Return type:

list(numpy.array or tuple or list)

ampal.analyse_protein.make_primitive(cas_coords, window_length=3)[source]

Calculates running average of cas_coords with a fixed averaging window_length.

Parameters:
  • cas_coords (list(numpy.array or float or tuple)) – Each element of the list must have length 3.

  • window_length (int, optional) – The number of coordinate sets to average each time.

Returns:

s_primitive – Each array has length 3.

Return type:

list(numpy.array)

Raises:

ValueError – If the length of cas_coords is smaller than the window_length.

ampal.analyse_protein.make_primitive_extrapolate_ends(cas_coords, smoothing_level=2)[source]

Generates smoothed helix primitives and extrapolates lost ends.

Notes

From an input list of CA coordinates, the running average is calculated to form a primitive. The smoothing_level dictates how many times to calculate the running average. A higher smoothing_level generates a ‘smoother’ primitive - i.e. the points on the primitive more closely fit a smooth curve in R^3. Each time the smoothing level is increased by 1, a point is lost from either end of the primitive. To correct for this, the primitive is extrapolated at the ends to approximate the lost values. There is a trade-off then between the smoothness of the primitive and its accuracy at the ends.

Parameters:
  • cas_coords (list(numpy.array or float or tuple)) – Each element of the list must have length 3.

  • smoothing_level (int) – Number of times to run the averaging.

Returns:

final_primitive – Each array has length 3.

Return type:

list(numpy.array)

ampal.analyse_protein.make_primitive_smoothed(cas_coords, smoothing_level=2)[source]

Generates smoothed primitive from a list of coordinates.

Parameters:
  • cas_coords (list(numpy.array or float or tuple)) – Each element of the list must have length 3.

  • smoothing_level (int, optional) – Number of times to run the averaging.

Returns:

s_primitive – Each array has length 3.

Return type:

list(numpy.array)

Raises:

ValueError – If the smoothing level is too great compared to the length of cas_coords.

ampal.analyse_protein.measure_sidechain_torsion_angles(residue, verbose=True) Tuple[List[float | None], List[int | None]][source]

Calculates sidechain dihedral angles for a residue

Parameters:
  • residue ([ampal.Residue]) – Residue object.

  • verbose (bool, optional) – If true, tells you when a residue does not have any known dihedral angles to measure.

Returns:

  • chi_angles ([float]) – Length depends on residue type, in range [-pi, pi]

    [0] = chi1 [if applicable] [1] = chi2 [if applicable] [2] = chi3 [if applicable] [3] = chi4 [if applicable]

  • rotamers ([int]) – Rotamer class according to Dunbrack 1993 {1, 2, 3}

    [0] = chi1 [if applicable] [1] = chi2 [if applicable] [2] = chi3 [if applicable] [3] = chi4 [if applicable]

ampal.analyse_protein.measure_torsion_angles(residues)[source]

Calculates the dihedral angles for a list of backbone atoms.

Parameters:

residues ([ampal.Residue]) – List of Residue objects.

Returns:

torsion_angles – One triple for each residue, containing torsion angles in the range [-pi, pi].

[0] omega [1] phi [2] psi

For the first residue, omega and phi are not defined. For the final residue, psi is not defined.

Return type:

(float, float, float)

Raises:

ValueError – If the number of input residues is less than 2.

ampal.analyse_protein.partial_charge(aa, pH)[source]

Calculates the partial charge of the amino acid.

Parameters:
  • aa (str) – Amino acid single-letter code.

  • pH (float) – pH of interest.

ampal.analyse_protein.polymer_to_reference_axis_distances(p, reference_axis, tag=True, reference_axis_name='ref_axis')[source]

Returns distances between the primitive of a Polymer and a reference_axis.

Notes

Distances are calculated between each point of the Polymer primitive and the corresponding point in reference_axis. In the special case of the helical barrel, if the Polymer is a helix and the reference_axis represents the centre of the barrel, then this function returns the radius of the barrel at each point on the helix primitive. The points of the primitive and the reference_axis are run through in the same order, so take care with the relative orientation of the reference axis when defining it.

Parameters:
  • p (ampal.Polymer) –

  • reference_axis (list(numpy.array or tuple or list)) – Length of reference_axis must equal length of the Polymer. Each element of reference_axis represents a point in R^3.

  • tag (bool, optional) – If True, tags the Chain with the reference axis coordinates and each Residue with its distance to the ref axis. Distances are stored at the Residue level, but refer to distances from the CA atom.

  • reference_axis_name (str, optional) – Used to name the keys in tags at Chain and Residue level.

Returns:

distances – Distance values between corresponding points on the reference axis and the Polymer Primitive.

Return type:

list(float)

Raises:

ValueError – If the Polymer and the reference_axis have unequal length.

ampal.analyse_protein.polypeptide_vector(p, start_index=0, end_index=-1, unit=True)[source]

Vector along the Chain primitive (default is from N-terminus to C-terminus).

Notes

start_index and end_index can be changed to examine smaller sections of the Chain, or reversed to change the direction of the vector.

Parameters:
  • p (ampal.Polymer) – Reference Polymer.

  • start_index (int, optional) – Default is 0 (start at the N-terminus of the Chain)

  • end_index (int, optional) – Default is -1 (start at the C-terminus of the Chain)

  • unit (bool) – If True, the vector returned has a magnitude of 1.

Returns:

vector – vector has shape (1, 3)

Return type:

a numpy.array

ampal.analyse_protein.reference_axis_from_chains(chains)[source]

Average coordinates from a set of primitives calculated from Chains.

Parameters:

chains (list(Chain)) –

Returns:

reference_axis – The averaged (x, y, z) coordinates of the primitives for the list of Chains. In the case of a coiled coil barrel, this would give the central axis for calculating e.g. Crick angles.

Return type:

numpy.array

Raises:

ValueError : – If the Chains are not all of the same length.

ampal.analyse_protein.residues_per_turn(p)[source]

The number of residues per turn at each Monomer in the Polymer.

Notes

Each element of the returned list is the number of residues per turn, at a point on the Polymer primitive. Calculated using the relative positions of the CA atoms and the primitive of the Polymer. Element i is the calculated from the dihedral angle using the CA atoms of the Monomers with indices [i, i+1] and the corresponding atoms of the primitive. The final value is None.

Parameters:

p (ampal.Polypeptide) – Polypeptide from which residues per turn will be calculated.

Returns:

rpts – Residue per turn values.

Return type:

[float]

ampal.analyse_protein.sequence_charge(seq, pH=7.4)[source]

Calculates the total charge of the input polypeptide sequence.

Parameters:
  • seq (str) – Sequence of amino acids.

  • pH (float) – pH of interest.

ampal.analyse_protein.sequence_isoelectric_point(seq, granularity=0.1)[source]

Calculates the isoelectric point of the sequence for ph 1-13.

Parameters:
  • seq (str) – Sequence of amino acids.

  • granularity (float, optional) – Granularity of pH values i.e. if 0.1 pH = [1.0, 1.1, 1.2…]

ampal.analyse_protein.sequence_molar_extinction_280(seq)[source]

Returns the molar extinction coefficient of the sequence at 280 nm.

Notes

Units = M/cm

Parameters:

seq (str) – Sequence of amino acids.

ampal.analyse_protein.sequence_molecular_weight(seq)[source]

Returns the molecular weight of the polypeptide sequence.

Notes

Units = Daltons

Parameters:

seq (str) – Sequence of amino acids.

ampal.assembly module

Defines various containers for AMPAL objects.

class ampal.assembly.AmpalContainer(ampal_objects=None, id=None)[source]

Bases: object

Custom list type class that holds multiple model states.

Notes

In this case, a state is defined as a set of coordinates that represents a protein model and an associated score or set of scores.

Parameters:
  • ampal_objects ([AMPAL], optional) – A list of AMPAL objects with which to initialise the AMPAL container. This can be an Assembly, Polymer or Monomer.

  • id (str, optional) – Identifier for the AMPAL container.

id

Identifier for the AMPAL container.

Type:

str

append(item)[source]

Adds an AMPAL object to the AmpalContainer.

extend(ampal_container)[source]

Extends an AmpalContainer with another AmpalContainer.

property pdb

Compiles the PDB strings for each state into a single file.

sort_by_tag(tag)[source]

Sorts the AmpalContainer by a tag on the component objects.

Parameters:

tag (str) – Key of tag used for sorting.

class ampal.assembly.Assembly(molecules=None, assembly_id='')[source]

Bases: BaseAmpal

A container that holds Polymer type objects.

Notes

Has a simple hierarchy: Assembly contains one or more Polymer, which in turn contains one or more Monomer.

Parameters:
  • molecules (Polymer or [Polymer], optional) – Polymer or list containing Polymer objects to be assembled.

  • assembly_id (str, optional) – An ID that the user can use to identify the Assembly. This is used when generating a pdb file using Assembly().pdb.

Raises:

TypeErrorAssembly objects can only be initialised empty, using a Polymer or a list of Polymers.

append(item)[source]

Adds a Polymer to the Assembly.

Raises:

TypeError – Raised if other is any type other than Polymer.

property backbone

Generates a new Assembly containing only the backbone atoms.

Notes

Metadata is not currently preserved from the parent object. Sequence data is retained, but only the main chain atoms are retained.

Returns:

bb_assemblyAssembly containing only the backbone atoms of the original Assembly.

Return type:

ampal.Protein

extend(assembly)[source]

Extends the Assembly with the contents of another Assembly.

Raises:

TypeError – Raised if other is any type other than Assembly.

property fasta

Generates a FASTA string for the Assembly.

Notes

Explanation of FASTA format: https://en.wikipedia.org/wiki/FASTA_format Recommendation that all lines of text be shorter than 80 characters is adhered to. Format of PDBID|CHAIN|SEQUENCE is consistent with files downloaded from the PDB. Uppercase PDBID used for consistency with files downloaded from the PDB. Useful for feeding into cdhit and then running sequence clustering.

Returns:

fasta_str – String of the fasta file for the Assembly.

Return type:

str

filter_mol_types(mol_types)[source]
get_atoms(ligands=True, pseudo_group=False, inc_alt_states=False)[source]

Flat list of all the Atoms in the Assembly.

Parameters:
  • ligands (bool, optional) – Include ligand Atoms.

  • pseudo_group (bool, optional) – Include pseudo_group Atoms.

  • inc_alt_states (bool, optional) – Include alternate sidechain conformations.

Returns:

atoms – All the Atoms as a iterator.

Return type:

itertools.chain

get_ligands(solvent=True)[source]

Retrieves all ligands from the Assembly.

Parameters:

solvent (bool, optional) – If True, solvent molecules will be included.

get_monomers(ligands=True, pseudo_group=False)[source]

Retrieves all the Monomers from the Assembly object.

Parameters:
  • ligands (bool, optional) – If true, will include ligand Monomers.

  • pseudo_group (bool, optional) – If True, will include pseudo atoms.

is_within(cutoff_dist, point, ligands=True)[source]

Returns all atoms in AMPAL object within cut-off distance from the point.

property isoelectric_point

Returns the isoelectric point of the Assembly.

make_pdb(ligands=True, alt_states=False, pseudo_group=False, header=True, footer=True)[source]

Generates a PDB string for the Assembly.

Parameters:
  • ligands (bool, optional) – If True, will include ligands in the output.

  • alt_states (bool, optional) – If True, will include alternate conformations in the output.

  • pseudo_group (bool, optional) – If True, will include pseudo atoms in the output.

  • header (bool, optional) – If True will write a header for output.

  • footer (bool, optional) – If True will write a footer for output.

Returns:

pdb_str – String of the pdb for the Assembly. Generated by collating Polymer().pdb calls for the component Polymers.

Return type:

str

property molar_extinction_280

Returns the extinction co-efficient of the Assembly at 280 nm.

property molecular_weight

Returns the molecular weight of the Assembly in Daltons.

property pdb

Runs make_pdb in default mode.

property primitives

Generates a new Assembly containing the primitives of each Polymer.

Notes

Metadata is not currently preserved from the parent object.

Returns:

prim_assemblyAssembly containing only the primitives of the Polymers in the original Assembly.

Return type:

ampal.Protein

relabel_all()[source]

Relabels all Polymers, Monomers and Atoms with default labeling.

relabel_atoms(start=1)[source]

Relabels all Atoms in numerical order, offset by the start parameter.

Parameters:

start (int, optional) – Defines an offset for the labelling.

relabel_monomers()[source]

Relabels all Monomers in the component Polymers in numerical order.

relabel_polymers(labels=None)[source]
Relabels the component Polymers either in alphabetical order or

using a list of labels.

Parameters:

labels (list, optional) – A list of new labels.

Raises:

ValueError – Raised if the number of labels does not match the number of component Polymer objects.

property sequences

Returns the sequence of each Polymer in the Assembly as a list.

Returns:

sequences – List of sequences.

Return type:

[str]

tag_ca_geometry(force=False, reference_axis=None, reference_axis_name='ref_axis')[source]

Tags each Monomer in the Assembly with its helical geometry.

Parameters:
  • force (bool, optional) – If True the tag will be run even if Monomers are already tagged.

  • reference_axis (list(numpy.array or tuple or list), optional) – Coordinates to feed to geometry functions that depend on having a reference axis.

  • reference_axis_name (str, optional) – Used to name the keys in tags at Chain and Residue level.

tag_torsion_angles(force=False)[source]

Tags each Monomer in the Assembly with its torsion angles.

Parameters:

force (bool, optional) – If True, the tag will be run even if Monomers are already tagged.

ampal.base_ampal module

Contains the base and common classes for all AMPAL objects.

class ampal.base_ampal.Atom(coordinates, element, atom_id=' ', res_label=None, occupancy=1.0, bfactor=1.0, charge=' ', state='A', parent=None)[source]

Bases: object

Object containing atomic coordinates and element information.

Notes

Atom is an AMPAL object, but it does not inherit from BaseAmpal.

Parameters:
  • coordinates (3D Vector (tuple, list, numpy.array)) – Position of Atom in 3D space.

  • element (str) – Element of Atom.

  • atom_id (str) – Identifier for Atom, usually a number.

  • res_label (str, optional) – Label used in Monomer to refer to the Atom type i.e. “CA” or “OD1”.

  • occupancy (float, optional) – The occupancy of the Atom.

  • bfactor (float, optional) – The bfactor of the Atom.

  • charge (str, optional) – The point charge of the Atom.

  • state (str) – The state of this Atom. Used to identify Atoms with a number of conformations.

  • parent (ampal.Monomer, optional) – A reference to the Monomer containing this Atom.

id

Identifier for Atom, usually a number.

Type:

str

res_label

Label used in Monomer to refer to the Atom type i.e. “CA” or “OD1”.

Type:

str

element

Element of Atom.

Type:

str

parent
A reference to the Monomer containing this Atom.

number of conformations.

Type:

ampal.Monomer

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

property array

Converts the atomic coordinate to a numpy.array.

rotate(angle, axis, point=None, radians=False)[source]

Rotates Atom by angle.

Parameters:
  • angle (float) – Angle that Atom will be rotated.

  • axis (3D Vector (tuple, list, numpy.array)) – Axis about which the Atom will be rotated.

  • point (3D Vector (tuple, list, numpy.array), optional) – Point that the axis lies upon. If None then the origin is used.

  • radians (bool, optional) – True is angle is define in radians, False is degrees.

translate(vector)[source]

Translates Atom.

Parameters:
  • vector (3D Vector (tuple, list, numpy.array)) – Vector used for translation.

  • inc_alt_states (bool, optional) – If true, will rotate atoms in all states i.e. includes alternate conformations for sidechains.

property unique_id

Creates a unique ID for the Atom based on its parents.

Returns:

unique_id – (polymer.id, residue.id, atom.id)

Return type:

(str, str, str)

property x

The x coordinate.

property y

The y coordinate.

property z

The z coordinate.

class ampal.base_ampal.BaseAmpal[source]

Bases: object

Base class for all AMPAL objects except ampal.atom.

Raises:

NotImplementedErrorBaseAmpal is an abstract base class and is not intended to be instanciated. A NotImplementedError is raised if a method is called that is required on a child class but is not implemented in BaseAmpal.

property centre_of_mass

Returns the centre of mass of AMPAL object.

Notes

All atoms are included in calculation, call centre_of_mass manually if another selection is require.

Returns:

centre_of_mass – 3D coordinate for the centre of mass.

Return type:

numpy.array

get_atoms(ligands=True, inc_alt_states=False)[source]
is_within(cutoff_dist, point)[source]

Returns all atoms in ampal object within cut-off distance from the point.

make_pdb()[source]
property pdb

Runs make_pdb in default mode.

rmsd(other, backbone=False)[source]

Calculates the RMSD between two AMPAL objects.

Notes

No fitting operation is performs and both AMPAL objects must have the same number of atoms.

Parameters:
  • other (AMPAL Object) – Any AMPAL object with get_atoms method.

  • backbone (bool, optional) – Calculates RMSD of backbone only.

rotate(angle, axis, point=None, radians=False, inc_alt_states=True)[source]

Rotates every atom in the AMPAL object.

Parameters:
  • angle (float) – Angle that AMPAL object will be rotated.

  • axis (3D Vector (tuple, list, numpy.array)) – Axis about which the AMPAL object will be rotated.

  • point (3D Vector (tuple, list, numpy.array), optional) – Point that the axis lies upon. If None then the origin is used.

  • radians (bool, optional) – True is angle is define in radians, False is degrees.

  • inc_alt_states (bool, optional) – If true, will rotate atoms in all states i.e. includes alternate conformations for sidechains.

translate(vector, inc_alt_states=True)[source]

Translates every atom in the AMPAL object.

Parameters:
  • vector (3D Vector (tuple, list, numpy.array)) – Vector used for translation.

  • inc_alt_states (bool, optional) – If true, will rotate atoms in all states i.e. includes alternate conformations for sidechains.

class ampal.base_ampal.Monomer(atoms=None, monomer_id=' ', parent=None)[source]

Bases: BaseAmpal

Groups of Atoms that form Polymers.

Parameters:
  • atoms (OrderedDict or {OrderedDict}, optional) – OrderedDict containing Atoms for the Monomer. OrderedDict is used to maintain the order items were added to the dictionary.

  • monomer_id (str, optional) – String used to identify the residue.

  • parent (Polymer, optional) – A reference to the Polymer containing this Monomer.

states

Contains an OrderedDicts containing atom information for each state available for the Monomer.

Type:

dict

id

String used to identify the residue.

Type:

str

parent

A reference to the Polymer containing this Monomer.

Type:

Polymer or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

property active_state

Defines which state dictionary should be used currently.

property atoms

Atoms in the currently active state.

close_monomers(group, cutoff=4.0)[source]

Returns a list of Monomers from within a cut off distance of the Monomer

Parameters:
  • group (BaseAmpal or Subclass) – Group to be search for Monomers that are close to this Monomer.

  • cutoff (float) – Distance cut off.

Returns:

nearby_residues – List of Monomers within cut off distance.

Return type:

[Monomers]

get_atoms(inc_alt_states=False)[source]

Returns all atoms in the Monomer.

Parameters:

inc_alt_states (bool, optional) – If True, will return Atoms for alternate states.

get_monomers()[source]

Returns the this Monomer.

Notes

This function is only present for consistency in the interface.

make_pdb()[source]

Generates a PDB string for the Monomer.

class ampal.base_ampal.Polymer(monomers=None, ligands=None, polymer_id=' ', molecule_type='', parent=None, sl=2)[source]

Bases: BaseAmpal

A container that holds Monomer type objects.

Notes

Polymer has a simple hierarchy: A Polymer contains one or more Monomer.

Parameters:
  • monomers (Monomer or [Monomer], optional) – Monomer or list containing Monomer objects to form the Polymer().

  • ligands (LigandGroup, optional) – Ligands associated with the Polymer.

  • polymer_id (str, optional) – An ID that the user can use to identify the Polymer. This is used when generating a pdb file using Polymer().pdb.

  • molecule_type (str, optional) – A description of the type of Polymer i.e. Protein, DNA etc.

  • parent (ampal.Assembly, optional) – Reference to Assembly containing the Polymer.

  • sl (int, optional) – The default smoothing level used when calculating the backbone primitive.

id

Polymer ID

Type:

str

parent

Reference to Assembly containing the Polymer.

Type:

ampal.Assembly or None

molecule_type

A description of the type of Polymer i.e. Protein, DNA etc.

Type:

str

ligands

A LigandGroup containing all the Ligands associated with this Polymer chain.

Type:

ampal.LigandGroup

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

sl

The default smoothing level used when calculating the backbone primitive.

Type:

int

Raises:

TypeError – Polymer objects can only be initialised empty, using a Monomer or a list of Monomers.

append(item)[source]

Appends a Monomer to the `Polymer.

Notes

Does not update labelling.

extend(polymer)[source]

Extends the Polymer with the contents of another Polymer.

Notes

Does not update labelling.

get_atoms(ligands=True, inc_alt_states=False)[source]

Flat list of all the Atoms in the Polymer.

Parameters:

inc_alt_states (bool) – If true atoms from alternate conformations are included rather than only the “active” states.

Returns:

atoms – Returns an iterator of all the atoms. Convert to list if you require indexing.

Return type:

itertools.chain

get_monomers(ligands=True)[source]

Retrieves all the Monomers from the AMPAL object.

Parameters:

ligands (bool, optional) – If true, will include ligand Monomers.

get_reference_coords()[source]

Gets list of coordinates of all reference atoms in the Polymer.

Returns:

ref_coords – List has the same length as the Polymer. The first, second and third elements of array i contain the x, y and z coordinates of the i-th reference atom.

Return type:

[numpy.array]

make_pdb(alt_states=False, inc_ligands=True)[source]

Generates a PDB string for the Polymer.

Parameters:
  • alt_states (bool, optional) – Include alternate conformations for Monomers in PDB.

  • inc_ligands (bool, optional) – Includes Ligands in PDB.

Returns:

pdb_str – String of the pdb for the Polymer. Generated using information from the component Monomers.

Return type:

str

relabel_all()[source]

Relabels all Monomers and Atoms using default labeling.

relabel_atoms(start=1)[source]

Relabels all Atoms in numerical order.

Parameters:

start (int, optional) – Offset the labelling by start residues.

relabel_monomers(labels=None)[source]

Relabels the either in numerically or using a list of labels.

Parameters:

labels (list, optional) – A list of new labels.

Raises:

ValueError – Raised if the number of labels does not match the number of component Monoer objects.

ampal.base_ampal.cap(v, l)[source]

Shortens string is above certain length.

ampal.base_ampal.centre_of_atoms(atoms, mass_weighted=True)[source]

Returns centre point of any list of atoms.

Parameters:
  • atoms (list) – List of AMPAL atom objects.

  • mass_weighted (bool, optional) – If True returns centre of mass, otherwise just geometric centre of points.

Returns:

centre_of_mass – 3D coordinate for the centre of mass.

Return type:

numpy.array

ampal.base_ampal.find_atoms_within_distance(atoms, cutoff_distance, point)[source]

Returns atoms within the distance from the point.

Parameters:
  • atoms ([ampal.atom]) – A list of ampal.atoms.

  • cutoff_distance (float) – Maximum distance from point.

  • point ((float, float, float)) – Reference point, 3D coordinate.

Returns:

filtered_atomsatoms list filtered by distance.

Return type:

[ampal.atoms]

ampal.base_ampal.write_pdb(residues, chain_id=' ', alt_states=False, strip_states=False)[source]

Writes a pdb file for a list of residues.

Parameters:
  • residues (list) – List of Residue objects.

  • chain_id (str) – String of the chain id, defaults to ‘ ‘.

  • alt_states (bool, optional) – If true, include all occupancy states of residues, else outputs primary state only.

  • strip_states (bool, optional) – If true, remove all state labels from residues. Only use with alt_states false.

Returns:

pdb_str – String of the PDB file.

Return type:

str

ampal.data module

Loads data required for running ISAMBARD.

ampal.dssp module

This module provides an interface to the program DSSP.

For more information on DSSP see [4].

References

ampal.dssp.dssp_available()[source]

True if mkdssp is available on the path.

ampal.dssp.extract_all_ss_dssp(in_dssp, path=True)[source]

Uses DSSP to extract secondary structure information on every residue.

Parameters:
  • in_dssp (str) – Path to DSSP file.

  • path (bool, optional) – Indicates if pdb is a path or a string.

Returns:

dssp_residues

Each internal list contains:

[0] int Residue number [1] str Secondary structure type [2] str Chain identifier [3] str Residue type [4] float Phi torsion angle [5] float Psi torsion angle [6] int dssp solvent accessibility

Return type:

[tuple]

ampal.dssp.find_ss_regions(dssp_residues, loop_assignments=(' ', 'B', 'S', 'T'))[source]

Separates parsed DSSP data into groups of secondary structure.

Notes

Example: all residues in a single helix/loop/strand will be gathered into a list, then the next secondary structure element will be gathered into a separate list, and so on.

Parameters:

dssp_residues ([tuple]) –

Each internal list contains:

[0] int Residue number [1] str Secondary structure type [2] str Chain identifier [3] str Residue type [4] float Phi torsion angle [5] float Psi torsion angle [6] int dssp solvent accessibility

Returns:

fragments – Lists grouped in continuous regions of secondary structure. Innermost list has the same format as above.

Return type:

[[list]]

ampal.dssp.get_ss_regions(assembly, ss_types)[source]

Returns an Assembly containing Polymers for each region of structure.

Parameters:
  • assembly (ampal.Assembly) – Assembly object to be searched secondary structure regions.

  • ss_types (list) – List of secondary structure tags to be separate i.e. [‘H’] would return helices, [‘H’, ‘E’] would return helices and strands.

Returns:

fragmentsAssembly containing a Polymer for each region of specified secondary structure.

Return type:

Assembly

ampal.dssp.run_dssp(pdb, path=True)[source]

Uses DSSP to find helices and extracts helices from a pdb file or string. :param pdb: Path to pdb file or string. :type pdb: str :param path: Indicates if pdb is a path or a string. :type path: bool, optional

Returns:

dssp_out – Std out from DSSP.

Return type:

str

ampal.dssp.tag_dssp_data(assembly, loop_assignments=(' ', 'B', 'S', 'T'))[source]

Adds output data from DSSP to an Assembly.

A dictionary will be added to the tags dictionary of each residue called dssp_data, which contains the secondary structure definition, solvent accessibility phi and psi values from DSSP. A list of regions of continuous secondary assignments will also be added to each Polypeptide.

The tags are added in place, so nothing is returned from this function.

Parameters:
  • assembly (ampal.Assembly) – An Assembly containing some protein.

  • loop_assignments (tuple or list) – A tuple containing the DSSP secondary structure identifiers to that are classed as loop regions.

ampal.geometry module

Basic geometric operations that AMPAL requires.

Functions that have names preceded by an underscore indicate that they are not to be used by users. However, all of these functions will generally be defined using cdef and so will not be imported anyway.

class ampal.geometry.Axis(start, end)

Bases: Curve

Axis: a straight line curve in 3-D space, defined by its start and end points (3-D vectors).

arc_length(t, t_ref=0)

arc length between point t and point t0

Notes

This is a trivial function in this particular case. Implemented for compatability with Curve class. Returns arc length along the curve between the points specified by t and t_ref. If t_ref > t then the returned value is negative.

Parameters:
  • t (float) –

  • t_ref (float) – t and t_ref define points on curve (See self.point()).

property direction_vector
distance_to_point(point)
get_coords(n_points=10, spacing=1.52, t_ref=0)

List of evenly spaced points along the curve of specified total length.

Parameters:
  • n_points (int) – The number of points to be returned.

  • spacing (float) – arc-length spacing of points to be returned.

  • t_ref (The value of t for specifying the first point in the coordinates list.) –

Returns:

coords – 3D coordinates of evenly spaced points along the curve.

Return type:

list of numpy.arrays

property length
property midpoint
point(t)
t_from_arc_length(arc_length, t_ref=0)
Returns parameter t that corresponds to a point at distance=arc_length along the curve from the curve

at t_ref.

Notes

This is a trivial function in this particular case. Implemented for compatability with Curve class.

Parameters:
  • arc_length (float) – distance along curve

  • t_ref (float) –

Returns:

Value of parameter t that corresponds to a point at distance=arc_length along the curve from the point defined by parameter t_ref.

Return type:

float

property unit_binormal

Binormal defined as (Tangent X Normal).

property unit_normal

Normal is perpendicular to tangent. See notes.

Notes

Since the curvature of a straight line is zero, the normal and binormal are not well-defined. https://en.wikipedia.org/wiki/Differential_geometry_of_curves#Special_Frenet_vectors_and_generalized_curvatures Here we have chosen to use the following convention:

If the Axis is parallel to the z-axis (i.e. unit_tangent = [0, 0, 1]),

then the normal is [1, 0, 0] (x-axis) and the binormal is [0, 1, 0] (y-axis).

If the Axis is not parallel to the z-axis, then we can rotate the z-axis so it aligns with our Axis. The normal is the same rotation appled to the x-axis. It’s helpful here to think of the right-hand rule. https://upload.wikimedia.org/wikipedia/commons/thumb/d/d2/Right_hand_rule_cross_product.svg/2000px-Right_hand_rule_cross_product.svg.png https://en.wikipedia.org/wiki/Cross_product#Definition Using the right-hand rule, the Axis direction vector is your index finger, the normal vector is the middle finger and the binormal is the thumb.

property unit_tangent

Tangent points along the direction of the line.

class ampal.geometry.Curve

Bases: object

Not intended to be instantiated. This is a base class for other Curve classes to inherit from.

Notes

Any curve must have a point method, which returns the coordinates of a point, given a parameter t.

point(t)
class ampal.geometry.HelicalCurve(a, b, handedness='r', t0=0)

Bases: Curve

HelicalCurve, parameterised by a and b.

property alpha

Angle between tangent to helix and the helix axis, in degrees.

classmethod alpha_and_pitch(alpha, pitch, handedness='r', t0=0, radians=False)
classmethod alpha_and_radius(alpha, radius, handedness='r', t0=0, radians=False)
arc_length(t, t_ref=None)

arc length between point t and point t0

Notes

Returns arc length along the curve between the points specified by t and t_ref. If t_ref > t then the returned value is negative.

Parameters:
  • t (float) –

  • t_ref (float) – t and t_ref define points on curve (See self.point()).

property axis
curvature()
get_coords(n_points=10, spacing=1.52, t_ref=None)

List of evenly spaced points along the curve of specified total length.

Parameters:
  • n_points (int) – The number of points to be returned.

  • spacing (float) – arc-length spacing of points to be returned.

  • axis_start (list or tuple or numpy.array of length 3) –

  • axis_end (list or tuple or numpy.array of length 3) – Coordinates returned will be for a helical path winding around an axis specified by axis_end - axis_start.

  • t_ref (The value of t for specifying the first point in the coordinates list.) –

Returns:

coords – 3D coordinates of evenly spaced points along the curve.

Return type:

list of numpy.arrays

property pitch
classmethod pitch_and_radius(pitch, radius, handedness='r', t0=0)
point(t)
property radius
t_from_arc_length(arc_length, t_ref=None)
Returns parameter t that corresponds to a point at distance=arc_length along the curve from the curve

at t_ref.

Parameters:
  • arc_length (float) – distance along curve

  • t_ref (float) –

Returns:

Value of parameter t that corresponds to a point at distance=arc_length along the curve from the point defined by parameter t_ref.

Return type:

float

torsion()
unit_binormal(t)
unit_normal(t)

Unit vector normal to the curve. Together with unit_tangent, defined the osculating plane at t

Notes

See this article for more information on the tangent, normal and binormal vectors. https://en.wikipedia.org/wiki/Differential_geometry_of_curves#Special_Frenet_vectors_and_generalized_curvatures

unit_tangent(t)

Unit vector pointing along the direction of the curve.

x(t)
y(t)
z(t)
class ampal.geometry.Quaternion(r, i, j, k, is_rotation=False)

Bases: object

classmethod angle_and_axis(angle, axis, radians=False)
apply_conjugation(other)
property as_array
as_rotation(radians=False, force=False)
property conjugate
distance(other)
extract_angle_and_axis(radians=False)

Returns the angle and axis of the rotation that the Quaternion represents.

property norm
classmethod real_and_vector(r, vector, is_rotation=False)
rotate_vector(v, point=None)
property vector
ampal.geometry.angle_between_vectors()

Calculate angle between two vectors using the dot product.

Notes

All angle-based functions in tools_geometry will return values in degrees by default. The radians flag can be used to return values in radians. If a or b is very close to the origin (norm < 1e-7), returns 0.0.

Parameters:
  • a (list or tuple or numpy.array.) –

  • b (list or tuple or numpy.array.) –

  • radians (bool, optional) – If True, returns angle in radians. If False, returns angle in degrees. Default value is False.

Returns:

angle – The angle between vectors v1 and v2.

Return type:

double

ampal.geometry.apply_unit_quaternion()

Applies the rotation defined by the unit quaternion to a 3D vector.

Parameters:
  • vector ([float, float, float]) – 3D coordinate.

  • u_in (double[3]) – Imaginary component of the unit quaternion.

  • s (double) – Real component of the quaternion.

Returns:

r – 3D coordinate.

Return type:

double[3]

ampal.geometry.cartesian_to_cylindrical()

Converts from cartesian to cylindrical coordinates system

Notes

See cylindrical_to_cartesian for more information.

Parameters:
  • x (float) –

  • y (float) –

  • z (float) –

  • radians (bool) – True if returned azimuth is in radians instead of degrees.

Returns:

  • radius (float)

  • azimuth (float)

  • z (float) – Cylindrical coordinate parameters. See cylindrical_to_cartesian for more information.

ampal.geometry.cartesian_to_spherical()

Converts from cartesian to spherical coordinates system

Notes

See spherical_to_cartesian for more information.

Parameters:
  • x (float) –

  • y (float) –

  • z (float) –

  • radians (bool) – True if returned azimuth and zenith are in radians instead of degrees.

Returns:

  • radius (float)

  • azimuth (float)

  • zenith (float) – Spherical coordinate parameters. See spherical_to_cartesian for more information.

ampal.geometry.centre_of_mass()

Find centre of mass of list of points.

Notes

Uniform weighting applied if masses not supplied.

Parameters:
  • points (list(tuple, list or numpy.array)) –

  • masses (list(float)) –

Returns:

com – Centre of mass of points.

Return type:

numpy.array

Raises:

ValueError – if len(masses) != len(points)

ampal.geometry.closest_distance()

Minimun distance between two lists of points in R^3.

Distances are calculated pairwise between points in points1 and points2. The minimum of all these distance values is returned.

Parameters:
  • points1 (list or tuple or numpy.array of triples.) – A number of points in R^3.

  • points2 (list or tuple or numpy.array of triples.) – A number of points in R^3.

Returns:

min_dist – The minimum distance between any point from points1 and any point from point2.

Return type:

float

ampal.geometry.closest_points_on_lines()

Notes

p0, p1, q0, q1 describe two lines, L1 (p0 -> p1) and L2 (q0 -> q1). L1 (2) is the infinite extension of the line that passes through a (q0) and has direction p0->p1 (q0->q1). Two points point_1, and point_2 on L1 and L2 respectively, are returned. The distance between point_1 and point_2 is smaller than any other points on L1 and L2. Code adapted from http://geomalgorithms.com/a07-_distance.html, where a more detailed explanation is available.

Parameters:
  • p0 (list or tuple or numpy.array) –

  • p1 (list or tuple or numpy.array) –

  • q0 (list or tuple or numpy.array) –

  • q1 (list or tuple or numpy.array) –

  • segments (bool) – If True, the points returned must lie on line segments p0->p1 and q0->q1, rather than their infinite extensions.

Returns:

  • point_1 (numpy.array)

  • point_2 (numpy.array)

ampal.geometry.cylindrical_to_cartesian()

Converts from cylindrical to cartesian coordinates system

Notes

Cylindrical coordinates are here defined by a radius, an azimuth angle and a z value. The azimuth is the angle in the x-y plane, with counter-clockwise direction being positive. So: azimuth = 0 => lies on x-axis, azimuth = 90 => lies in y axis, zenith = 180 => lies antiparallel to x-axis. The radius and the z value are the radius and height of the cylinder, respectively. See https://en.wikipedia.org/wiki/Cylindrical_coordinate_system and/or http://mathworld.wolfram.com/CylindricalCoordinates.html for more information.

Parameters:
  • radius (float) –

  • azimuth (float) –

  • z (float) –

  • radians (bool) – True if azimuth is given in radians instead of degrees.

Returns:

cartesian coordinates corresponding to the input cylindrical coordinates.

Return type:

numpy.array([x, y, z])

ampal.geometry.dihedral()

Calculates the dihedral angle defined by four points in R^3.

The points a, b and c form a plane in R^3. The points b, c and d form a plane in R^3. The dihedral angle is the angle between these two planes.

Notes

See http://en.wikipedia.org/wiki/Dihedral_angle and https://en.wikipedia.org/wiki/Atan2 for more information. Dihedral between 4 points. a, b, c, d are all 3-vectors (x,y,z) coordinates. (a,b,c) and (b,c,d) must both be non-colinear.

Parameters:
  • a (list or tuple or numpy.array) –

  • b (list or tuple or numpy.array) –

  • c (list or tuple or numpy.array) –

  • d (list or tuple or numpy.array) –

  • radians (bool, optional) – If True, returns angle in radians. If False, returns angle in degrees. Default value is False.

Returns:

angle – The dihedral angle formed between the points a, b, c and d.

Return type:

double

ampal.geometry.distance()

Calculates the Euclidian distance between two points in R^3.

Parameters:
  • a (list or tuple or numpy.array) –

  • b (list or tuple or numpy.array) –

Returns:

dist – Euclidean distance between a and b.

Return type:

double

ampal.geometry.find_foot()

Finds the perpendicular foot of a point (p) and a line (ab).

The points a and b define a line in R^3. There is a point on this line, q, that has the shortest Euclidian distance from the point p. This point is the perpendicular foot, and is returned by this function.

Notes

See http://mathworld.wolfram.com/Point-LineDistance3-Dimensional.html for a useful work-through of the maths. Parameterise the line ab by t, so that a point on the line is of the form: a + (b-a)*t. Define f(t) = [distance(t, p)]^2. Want t such that f(t) is at a minimum. Solve f’(t) = 0 for t. Substituting this value of t into a + (b-a)*t yields the desired point q.

Parameters:
  • a (list or tuple or numpy.array) –

  • b (list or tuple or numpy.array) – Together, the points a and b define a line in R^3.

  • p (list or tuple or numpy.array) – A point in R^3.

Returns:

q – The point on the line between a and b that is closest to q.

Return type:

list or tuple or numpy.array

ampal.geometry.find_foot_on_plane()

Finds the perpendicular foot of a point (p) and a plane defined by points a, b and c.

The points a, b and c must define a plane in R^3 (i.e. they cannot be collinear). There is a point on this plane, q, that has the shortest Euclidian distance from the point p. This point is the perpendicular foot, and is returned by this function.

Notes

This page:

http://math.stackexchange.com/questions/723937/find-the-point-on-a-plane-3x-4y-z-1-that-is-closest-to-1-0-1 gives a work-through of this kind of problem for a specific example.

The solution uses the Lagrangian multiplier: https://en.wikipedia.org/wiki/Lagrange_multiplier We look for the minimum values of the distance function f (squared to make it easier!)

f(x, y, z) = (x - p[0])^2 + (y - p[1])^2 + (z - p[2])^2

subject to the constraint g(x, y, z) = 0, where g is the equation of the plane. The Langrangian multiplier L = f(x, y, z) - lambda*(g(x, y, z)). We calculate the partial deriviatves of this (wrt x, y, z and lambda). Solve first for lambda gives the equation for the lagrangian_lambda value in this function. Then substitute lambda back in to get the solution for (x, y, z) = q.

Parameters:
  • a (list or tuple or numpy.array) –

  • b (list or tuple or numpy.array) –

  • c (list or tuple or numpy.array) – Together, the points a, b and c define a plane in R^3.

  • p (list or tuple or numpy.array) – A point in R^3.

Returns:

s – The point on the plane that is closest to q.

Return type:

numpy.array

ampal.geometry.find_limits()
ampal.geometry.find_transformations()

Returns geometrical transformations required to align the vector from s1 -> e1 with the vector from s2 -> e2

Notes

Applying the transformations returned to the line from s1 -> e2, will align it with the line from s2 -> e2. After transformation, s1 will have been moved onto s2. Apply the rotation (given by angle, axis and point) first and then the translation (given by translation).

Parameters:
  • s1 ([float, float, float]) – 3D coordinate. Start of line 1

  • e1 ([float, float, float]) – 3D coordinate. End of line 1

  • s2 ([float, float, float]) – 3D coordinate. Start of line 2

  • e2 ([float, float, float]) – 3D coordinate. End of line 2

  • radians (bool) – If True, angle will be returned in radians.

Returns:

  • translation (numpy.array(3)) – 3D vector. Translation to be applied to s1 to move it onto s2.

  • angle (float) – angle for rotation (default: degrees).

  • axis (numpy.array(3)) – axis about which to rotate line 1 to align it with line 2.

  • point (numpy.array(3)) – a point lying on the rotation axis.

ampal.geometry.gen_sectors()

Separates the coordinates into overlapping sectors.

ampal.geometry.intersection_start_and_end()

Find start and end indices of the part of the list of points ‘subject’ that will intersect with ‘reference’

Notes

This is a function designed to find assemblies from within bundles of helices. The idea here is to take two lists of points, a subject and a reference. We want to find the part of the subject (this could be all of it, or none of it) that would intersect the line joining the points in the reference, if the line joining the points of subject was the prinicipal axis of a cylinder of infinite radius. If subject was a list of ten points joining (0, 0 , 0) and (0, 0, 10), and reference was a list of 1000 points joining (10, 0, 7) and (10, 0, 12) this function would return (7, 10). If the subject and reference for the above example were swapped around, however, the function would return (0, 599).

Parameters:
  • subject (list of numpy.array) – A list of points in R^3.

  • reference (list of numpy.array) – A list of points in R^3.

Returns:

  • start_index (int) – The index of the first point in subject that is part of the intersection described above.

  • end_index (int) – The index of the final point in subject that is part of the intersection described above.

ampal.geometry.is_acute()

Return True if angle between two vectors is less than 90 degrees.

Parameters:
  • a (list or tuple or numpy.array) –

  • b (list or tuple or numpy.array) –

Returns:

True if angle between v1 and v2 is less than 90 degrees. False otherwise.

Return type:

bool

ampal.geometry.minimal_distance_between_lines()

Finds the minimal distance between two the two lines joining p0 to p1 and q0 to q1.

The points p0, p1, q0 and d are used to construct two lines in R^3. The first line L1 is from p0 to p1, the second line L2 is from q0 to q1. The minimal distance between these lines is the length of the line that connects them and is perpendicular to both L1 and L2.

Notes

See https://en.wikipedia.org/wiki/Skew_lines#Distance for more information. Also http://geomalgorithms.com/a07-_distance.html

Parameters:
  • p0 (list or tuple or numpy.array) –

  • p1 (list or tuple or numpy.array) –

  • q0 (list or tuple or numpy.array) –

  • q1 (list or tuple or numpy.array) –

  • segments (bool) – If True, the points returned must lie on line segments p0->p1 and q0->q1, rather than their infinite extensions.

Returns:

dist – The distance between the two lines.

Return type:

float

ampal.geometry.points_on_a_circle()

List of n uniformly distributed (x, y) coordinates on the circumference of a circle.

Parameters:
  • n (int) – Number of points to return.

  • radius (float) –

  • centre (tuple or list or numpy.array) –

  • rotation (float) – Angle in degrees by which all points will be rotated. rotation = 0 means that the line from the centre to the first point is parallel to the x-axis. rotation > 0 => anti-clockwise rotation

Returns:

points – (x, y) coordinates of uniformly distributed points on the circumference of the circle

Return type:

list(2-tuples)

ampal.geometry.radius_of_circumcircle()

Returns the radius the circumcircle of a triangle, when supplied with three points (p, q, r).

Notes

See http://www.analyzemath.com/Geometry/Circumcircle/CircumcircleRadius.html for more information. Three (non-colinear) points define a traingle. The circumcircle is the circle that passes though these points.

Parameters:
  • p (list or tuple or numpy.array) –

  • q (list or tuple or numpy.array) –

  • r (list or tuple or numpy.array) –

Returns:

radius

Return type:

float

ampal.geometry.rmsd()

Returns average distance between points in points1 and their counterparts in points2.

Parameters:
  • points1 (list of [float, float, float]) –

  • points2 (list of [float, float, float]) –

Returns:

rmsd

Return type:

float

:raises ValueError : if len(points1) != len(points2):

ampal.geometry.rotation_matrix()

Find rotation matrix for rotating about unit_vect by angle.

Notes

See https://en.wikipedia.org/wiki/Rotation_matrix for more information. In particular, the section headed ‘Rotation matrix from axis and angle’.

Parameters:
  • angle (float) – Angle about which to rotate.

  • unit_vect (list or tuple or numpy.array) – Vector about which to rotate.

  • point (list or tuple or numpy.array, optional.) – Specify a point about which to rotate. If point is None, then the rotation will be about the origin. Default value is None.

Returns:

m – A (4 x 4) rotation matrix.

Return type:

numpy.array

ampal.geometry.spherical_to_cartesian()

Converts from spherical to cartesian coordinates system

Notes

Spherical coordinates are here defined by a radius, an azimuth angle and a zenith angle. The zenith is the angle between the point and the z-axis. So: zenith = 0 => lies on z-axis, zenith = 90 => lies in x-y plane, zenith = 180 => lies antiparallel to z-axis. The azimuth is the angle in the x-y plane, with counter-clockwise direction being positive. So: azimuth = 0 => lies on x-axis, azimuth = 90 => lies in y axis, zenith = 180 => lies antiparallel to x-axis. See https://en.wikipedia.org/wiki/Spherical_coordinate_system and/or http://mathworld.wolfram.com/SphericalCoordinates.html for more information. If zenith = 0, then x and y are both 0 regardless of azimuth value. If azimuth = 0, then y is 0 regardless of zenith value.

Parameters:
  • radius (float) –

  • azimuth (float) –

  • zenith (float) –

  • radians (bool) – True if azimuth and zenith are given in radians instead of degrees.

Returns:

cartesian coordinates corresponding to the input spherical coordinates.

Return type:

numpy.array([x, y, z])

ampal.geometry.unit_vector()

Return vector normalised by length.

Parameters:

a (list or tuple or numpy.array) – The vector to be normalised.

Returns:

v_out – A numpy.array that is a unit vector.

Return type:

numpy.array

Raises:

ZeroDivisionError – If the vector has length 0.

ampal.interactions module

Contains code for analysing chemical interactions in AMPAL objects.

class ampal.interactions.CovalentBond(a, b, dist)[source]

Bases: Interaction

Defines a covalent bond.

property a

One Atom involved in the covalent bond.

property b

One Atom involved in the covalent bond.

class ampal.interactions.HydrogenBond(donor, acceptor, dist, ang_d, ang_a)[source]

Bases: NonCovalentInteraction

Defines a hydrogen bond in terms of a donor and an acceptor.

Parameters:
  • donor (ampal.Atom) – The donor Atom in the interaction.

  • acceptor (ampal.Atom) – The acceptor atom in the interaction.

  • dist (float) – The distance between Atom a and b.

  • ang_a (float) – Angle between the acceptor and the interaction vector.

  • ang_d (float) – Angle between the donor and the interaction vector.

donor

The donor Atom in the interaction.

Type:

ampal.Atom

acceptor

The acceptor atom in the interaction.

Type:

ampal.Atom

dist

The distance between Atom a and b.

Type:

float

ang_a

Angle between the acceptor and the interaction vector.

Type:

float

ang_d

Angle between the donor and the interaction vector.

Type:

float

property acceptor_monomer

The acceptor Monomer in the interaction.

property donor_monomer

The donor Monomer in the interaction.

class ampal.interactions.Interaction(a, b, dist)[source]

Bases: object

A container for all types of interaction with donor and acceptor.

Parameters:
  • a (ampal.Atom) – A member of a pairwise interaction.

  • b (ampal.Atom) – A member of a pairwise interaction.

  • dist (float) – The distance between a and b.

a

A member of a pairwise interaction.

Type:

ampal.Atom

b

A member of a pairwise interaction.

Type:

ampal.Atom

dist

The distance between Atom a and b.

Type:

float

class ampal.interactions.NonCovalentInteraction(donor, acceptor, dist)[source]

Bases: Interaction

A container for all non-covalent interaction.

Parameters:
  • donor (ampal.Atom) – The donor Atom in the interaction.

  • acceptor (ampal.Atom) – The acceptor atom in the interaction.

  • dist (float) – The distance between Atom a and b.

donor

The donor Atom in the interaction.

Type:

ampal.Atom

acceptor

The acceptor atom in the interaction.

Type:

ampal.Atom

dist

The distance between Atom a and b.

Type:

float

property acceptor

The acceptor in the interaction.

property donor

The donor Atom in the interaction.

ampal.interactions.covalent_bonds(atoms, threshold=1.1)[source]

Returns all the covalent bonds in a list of Atom pairs.

Notes

Uses information ELEMENT_DATA, which can be accessed directly through this module i.e. isambard.ampal.interactions.ELEMENT_DATA.

Parameters:
  • atoms ([(Atom, Atom)]) – List of pairs of Atoms.

  • threshold (float, optional) – Allows deviation from ideal covalent bond distance to be included. For example, a value of 1.1 would allow interactions up to 10% further from the ideal distance to be included.

ampal.interactions.find_covalent_bonds(ampal, max_range=2.2, threshold=1.1, tag=True)[source]

Finds all covalent bonds in the AMPAL object.

Parameters:
  • ampal (AMPAL Object) – Any AMPAL object with a get_atoms method.

  • max_range (float, optional) – Used to define the sector size, so interactions at longer ranges will not be found.

  • threshold (float, optional) – Allows deviation from ideal covalent bond distance to be included. For example, a value of 1.1 would allow interactions up to 10% further from the ideal distance to be included.

  • tag (bool, optional) – If True, will add the covalent bond to the tags dictionary of each Atom involved in the interaction under the covalent_bonds key.

ampal.interactions.generate_bond_subgraphs_from_break(bond_graph, atom1, atom2)[source]

Splits the bond graph between two atoms to producing subgraphs.

Notes

This will not work if there are cycles in the bond graph.

Parameters:
  • bond_graph (networkx.Graph) – Graph of covalent bond network

  • atom1 (isambard.ampal.Atom) – First atom in the bond.

  • atom2 (isambard.ampal.Atom) – Second atom in the bond.

Returns:

subgraphs – A list of subgraphs generated when a bond is broken in the covalent bond network.

Return type:

[networkx.Graph]

ampal.interactions.generate_covalent_bond_graph(covalent_bonds)[source]

Generates a graph of the covalent bond network described by the interactions.

Parameters:

covalent_bonds ([CovalentBond]) – List of CovalentBond.

Returns:

bond_graph – A graph of the covalent bond network.

Return type:

networkx.Graph

ampal.ligands module

AMPAL objects that represent ligands.

class ampal.ligands.Ligand(mol_code, atoms=None, monomer_id=' ', insertion_code=' ', is_hetero=False, parent=None)[source]

Bases: Monomer

Monomer that represents a Ligand.

Notes

All Monomers that do not have dedicated classes are represented using the Ligand class.

Parameters:
  • mol_code (str) – PDB molecule code that represents the monomer.

  • atoms (OrderedDict, optional) – OrderedDict containing Atoms for the Monomer. OrderedDict is used to maintain the order items were added to the dictionary.

  • monomer_id (str, optional) – String used to identify the residue.

  • insertion_code (str, optional) – Insertion code of monomer, used if reading from pdb.

  • is_hetero (bool, optional) – True if is a hetero atom in pdb. Helps with PDB formatting.

atoms

OrderedDict containing Atoms for the Monomer. OrderedDict is used to maintain the order items were added to the dictionary.

Type:

OrderedDict

mol_code

PDB molecule code that represents the Ligand.

Type:

str

insertion_code

Insertion code of Ligand, used if reading from pdb.

Type:

str

is_hetero

True if is a hetero atom in pdb. Helps with PDB formatting.

Type:

bool

self.states

Contains an OrderedDicts containing atom information for each state available for the Ligand.

Type:

dict

id

String used to identify the residue.

Type:

str

parent

A reference to the LigandGroup containing this Ligand.

Type:

Polymer or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

class ampal.ligands.LigandGroup(monomers=None, polymer_id=' ', parent=None, sl=2)[source]

Bases: Polymer

A container for Ligand Monomers.

Parameters:
  • monomers (Monomer or [Monomer], optional) – Monomer or list containing Monomer objects to form the Polymer().

  • polymer_id (str, optional) – An ID that the user can use to identify the Polymer. This is used when generating a pdb file using Polymer().pdb.

  • parent (ampal.Assembly, optional) – Reference to Assembly containing the Polymer.

  • sl (int, optional) – The default smoothing level used when calculating the backbone primitive.

ampal.nucleic_acid module

Contains AMPAL objects representing nucleic acids.

class ampal.nucleic_acid.Nucleotide(atoms=None, mol_code='UNK', monomer_id=' ', insertion_code=' ', is_hetero=False, parent=None)[source]

Bases: Monomer

Represents a nucleotide base.

Parameters:
  • atoms (OrderedDict, optional) – OrderedDict containing Atoms for the Nucleotide. OrderedDict is used to maintain the order items were added to the dictionary.

  • mol_code (str, optional) – One or three letter code that represents the Nucleotide.

  • monomer_id (str, optional) – String used to identify the Nucleotide.

  • insertion_code (str, optional) – Insertion code of Nucleotide, used if reading from pdb.

  • is_hetero (bool, optional) – True if is a hetero atom in pdb. Helps with PDB formatting.

  • parent (ampal.Polynucleotide, optional) – Reference to Polynucleotide containing the Nucleotide.

mol_code

PDB molecule code that represents the Nucleotide.

Type:

str

insertion_code

Insertion code of Nucleotide, used if reading from pdb.

Type:

str

is_hetero

True if is a hetero atom in pdb. Helps with PDB formatting.

Type:

bool

states

Contains an OrderedDicts containing atom information for each state available for the Nucleotide.

Type:

dict

id

String used to identify the Nucleotide.

Type:

str

reference_atom

The key that corresponds to the reference Atom. This is used by various functions, for example backbone primitives are calculated using the Atom defined using this key.

Type:

str

parent

A reference to the Polynucleotide containing this Nucleotide.

Type:

Polynucleotide or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

Raises:

ValueError – Raised if mol_code is not length 1 or 3.

class ampal.nucleic_acid.Polynucleotide(monomers=None, polymer_id=' ', parent=None, sl=2)[source]

Bases: Polymer

Polymer type object that represents a Polynucleotide.

Parameters:
  • monomers (Nucleotide or [Nucleotide], optional) – Nucleotide or list containing Nucleotide objects to form the Polynucleotide.

  • polymer_id (str, optional) – An ID that the user can use to identify the Polynucleotide. This is used when generating a pdb file using Polynucleotide().pdb.

  • parent (ampal.Assembly, optional) – Reference to Assembly containing the Polynucleotide.

  • sl (int, optional) – The default smoothing level used when calculating the backbone primitive.

id

Polynucleotide ID

Type:

str

parent

Reference to Assembly containing the Polynucleotide

Type:

ampal.Assembly or None

molecule_type

A description of the type of Polymer i.e. Protein, DNA etc.

Type:

str

ligands

A LigandGroup containing all the Ligands associated with this Polynucleotide chain.

Type:

ampal.LigandGroup

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

sl

The default smoothing level used when calculating the backbone primitive.

Type:

int

Raises:

TypeErrorPolymer type objects can only be initialised empty or using a Monomer.

property sequence

Returns the sequence of the Polynucleotide as a string.

Returns:

sequence – String of the monomer sequence of the Polynucleotide.

Return type:

str

ampal.pdb_parser module

Contains code for parsing PDB files.

class ampal.pdb_parser.PdbParser(pdb, path=True, pdb_id='', ignore_end=False)[source]

Bases: object

Parses a PDB file and produces and AMPAL Assembly.

Parameters:
  • pdb (str) – Either a path to a PDB file or a string containing PDB format structural data.

  • path (bool, optional) – If true, flags pdb as a path and not a PDB string.

  • pdb_id (str, optional) – Identifier for the Assembly.

  • ignore_end (bool, optional) – If false, parsing of the file will stop when an “END” record is encountered.

proc_functions

Keys are PDB labels i.e. “ATOM” or “END”, values are the function that processes that specific line.

Type:

dict

id

Identifier for the Assembly.

Type:

str

pdb_lines

Input PDB split into line.

Type:

[str]

new_labels

Indicates if new atom or residue labels have been found while parsing the PDB file.

Type:

bool

state

Current state being appended to. This is used on multi-state files like NMR structures.

Type:

int

pdb_parse_tree

This dictionary represents the parse tree of the PDB file. Each line of the structure is broken down into a key, the entry label, and a value, the data.

Type:

dict

current_line

The line that is currently being parsed.

Type:

int

ignore_end

If false, parsing of the file will stop when an “END” record is encountered.

Type:

bool, optional

change_state()[source]

Increments current state and adds a new dict to the parse tree.

static check_for_non_canonical(residue)[source]

Checks to see if the residue is non-canonical.

end()[source]

Processes an “END” record.

gen_states(monomer_data, parent)[source]

Generates the states dictionary for a Monomer.

monomer_datalist

A list of atom data parsed from the input PDB.

parentampal.Monomer

Monomer used to assign parent on created Atoms.

make_ampal()[source]

Generates an AMPAL object from the parse tree.

Notes

Will create an Assembly if there is a single state in the parese tree or an AmpalContainer if there is more than one.

parse_pdb_file()[source]

Runs the PDB parser.

proc_atom()[source]

Processes an “ATOM” or “HETATM” record.

proc_chain(chain_info, parent)[source]

Converts a chain into a Polymer type object.

Parameters:
  • chain_info ((set, OrderedDict)) – Contains a set of chain labels and atom records.

  • parent (ampal.Assembly) – Assembly used to assign parent on created Polymer.

Raises:
  • ValueError – Raised if multiple or unknown atom types found within the same chain.

  • AttributeError – Raised if unknown Monomer type encountered.

proc_line_coordinate(line)[source]

Extracts data from columns in ATOM/HETATM record.

proc_monomer(monomer_info, parent, mon_cls=False)[source]

Processes a records into a Monomer.

Parameters:
  • monomer_info ((set, OrderedDict)) – Labels and data for a monomer.

  • parent (ampal.Polymer) – Polymer used to assign parent on created Monomer.

  • mon_cls (Monomer class or subclass, optional) – A Monomer class can be defined explicitly.

proc_state(state_data, state_id)[source]

Processes a state into an Assembly.

Parameters:
  • state_data (dict) – Contains information about the state, including all the per line structural data.

  • state_id (str) – ID given to Assembly that represents the state.

ampal.pdb_parser.load_pdb(pdb, path=True, pdb_id='', ignore_end=False)[source]

Converts a PDB file into an AMPAL object.

Parameters:
  • pdb (str) – Either a path to a PDB file or a string containing PDB format structural data.

  • path (bool, optional) – If true, flags pdb as a path and not a PDB string.

  • pdb_id (str, optional) – Identifier for the Assembly.

  • ignore_end (bool, optional) – If false, parsing of the file will stop when an “END” record is encountered.

Returns:

ampal – AMPAL object that contains the structural information from the PDB file provided. If the PDB file has a single state then an Assembly will be returned, otherwise an AmpalContainer will be returned.

Return type:

ampal.Assembly or ampal.AmpalContainer

ampal.protein module

AMPAL objects that represent protein.

class ampal.protein.Polypeptide(monomers=None, polymer_id=' ', parent=None, sl=2)[source]

Bases: Polymer

Container for Residues, inherits from Polymer.

Parameters:
  • monomers (Residue or [Residue], optional) – Residue or list containing Residue objects to form the Polypeptide.

  • polymer_id (str, optional) – An ID that the user can use to identify the Polypeptide. This is used when generating a pdb file using Polypeptide().pdb.

  • parent (ampal.Assembly, optional) – Reference to Assembly containing the Polymer.

  • sl (int, optional) – The default smoothing level used when calculating the backbone primitive.

id

Polypeptide ID

Type:

str

parent

Reference to Assembly containing the Polypeptide

Type:

ampal.Assembly or None

molecule_type

A description of the type of Polymer i.e. Protein, DNA etc.

Type:

str

ligands

A LigandGroup containing all the Ligands associated with this Polypeptide chain.

Type:

ampal.LigandGroup

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

sl

The default smoothing level used when calculating the backbone primitive.

Type:

int

Raises:

TypeErrorPolymer type objects can only be initialised empty or using a Monomer.

property backbone

Returns a new Polymer containing only the backbone atoms.

Notes

Metadata is not currently preserved from the parent object. Sequence data is retained, but only the main chain atoms are retained.

Returns:

bb_poly – Polymer containing only the backbone atoms of the original Polymer.

Return type:

Polypeptide

property backbone_bond_angles

Dictionary containing backbone bond angles as lists of floats.

Returns:

bond_angles – Keys are n_ca_c, ca_c_o, ca_c_n and c_n_ca, referring to the N-CA-C, CA-C=O, CA-C-N and C-N-CA angles respectively. Values are lists of floats : the bond angles in degrees. The lists of n_ca_c, ca_c_o are of length k for a Polypeptide containing k Residues. The list of ca_c_n and c_n_ca are of length k-1 for a Polypeptide containing k Residues (These angles are across the peptide bond, and are therefore formed between successive Residue pairs).

Return type:

dict

property backbone_bond_lengths

Dictionary containing backbone bond lengths as lists of floats.

Returns:

bond_lengths – Keys are n_ca, ca_c, c_o and c_n, referring to the N-CA, CA-C, C=O and C-N bonds respectively. Values are lists of floats : the bond lengths in Angstroms. The lists of n_ca, ca_c and c_o are of length k for a Polypeptide containing k Residues. The list of c_n bonds is of length k-1 for a Polypeptide containing k Residues (C-N formed between successive Residue pairs).

Return type:

dict

c_join(other, psi=-40.76, omega=-178.25, phi=-65.07, o_c_n_angle=None, c_n_ca_angle=None, c_n_length=None, relabel=True)[source]

Joins other to self at the C-terminus via a peptide bond.

Notes

This function directly modifies self. It does not return a new object.

Parameters:
  • other (Residue or Polypeptide) –

  • psi (float, optional) – Psi torsion angle (degrees) between final Residue of self and first Residue of other.

  • omega (float, optional) – Omega torsion angle (degrees) between final Residue of self and first Residue of other.

  • phi (float, optional) – Phi torsion angle (degrees) between final Residue of self and first Residue of other.

  • o_c_n_angle (float or None, optional) – Desired angle between O, C (final Residue of self) and N (first Residue of other) atoms. If None, default value is taken from ideal_backbone_bond_angles.

  • c_n_ca_angle (float or None, optional) – Desired angle between C (final Residue of self) and N, CA (first Residue of other) atoms. If None, default value is taken from ideal_backbone_bond_angles.

  • c_n_length (float or None, optional) – Desired peptide bond length between final Residue of self and first Residue of other. If None, default value is taken from ideal_backbone_bond_lengths.

  • relabel (bool, optional) – If True, relabel_all is run on self before returning.

Raises:

TypeError: – If other is not a Residue or a Polypeptide.

property fasta

Generates sequence data for the protein in FASTA format.

get_slice_from_res_id(start, end)[source]

Returns a new Polypeptide containing the Residues in start/end range.

Parameters:
  • start (str) – string representing start residue id (PDB numbering)

  • end (str) – string representing end residue id (PDB numbering)

Returns:

slice_polymer – Polymer containing the residue range specified by start-end

Return type:

Polymer

property isoelectric_point

Returns the isoelectric point of the Assembly.

property molar_extinction_280

Returns the extinction co-efficient of the Assembly at 280 nm.

property molecular_weight

Returns the molecular weight of the Assembly in Daltons.

n_join(other, psi=-40.76, omega=-178.25, phi=-65.07, o_c_n_angle=None, c_n_ca_angle=None, c_n_length=None, relabel=True)[source]

Joins other to self at the N-terminus via a peptide bond.

Notes

This function directly modifies self. It does not return a new object.

Parameters:
  • other (Residue or Polypeptide) –

  • psi (float) – Psi torsion angle (degrees) between final Residue of other and first Residue of self.

  • omega (float) – Omega torsion angle (degrees) between final Residue of other and first Residue of self.

  • phi (float) – Phi torsion angle (degrees) between final Residue of other and first Residue of self.

  • o_c_n_angle (float or None) – Desired angle between O, C (final Residue of other) and N (first Residue of self) atoms. If None, default value is taken from ideal_backbone_bond_angles.

  • c_n_ca_angle (float or None) – Desired angle between C (final Residue of other) and N, CA (first Residue of self) atoms. If None, default value is taken from ideal_backbone_bond_angles.

  • c_n_length (float or None) – Desired peptide bond length between final Residue of other and first Residue of self. If None, default value is taken from ideal_backbone_bond_lengths.

  • relabel (bool) – If True, relabel_all is run on self before returning.

Raises:

TypeError: – If other is not a Residue or a Polypeptide

property primitive

Primitive of the backbone.

Notes

This is the average of the positions of all the CAs in frames of sl Residues.

radii_of_curvature()[source]

List of radius of curvature values along the Polypeptide.

rise_per_residue()[source]

List of rise per residue values along the Polypeptide. .. rubric:: Notes

Calculated from Polypeptide.primitive.

property sequence

Returns the sequence of the Polymer as a string.

Returns:

sequence – String of the Residue sequence of the Polypeptide.

Return type:

str

tag_ca_geometry(force=False, reference_axis=None, reference_axis_name='ref_axis')[source]

Tags each Residue with rise_per_residue, radius_of_curvature and residues_per_turn.

Parameters:
  • force (bool, optional) – If True the tag will be run even if Residues are already tagged.

  • reference_axis (list(numpy.array or tuple or list), optional) – Coordinates to feed to geometry functions that depend on having a reference axis.

  • reference_axis_name (str, optional) – Used to name the keys in tags at Polypeptide and Residue level.

tag_sidechain_dihedrals(force=False)[source]

Tags each monomer with side-chain dihedral angles. For residues that do not have any rotamers (Alanine and Glycine) the rotamer is tagged as 0 and chi_angles as None.

force: bool, optional

If True the tag will be run even if Residues are already tagged.

tag_torsion_angles(force=False)[source]

Tags each Monomer of the Polymer with its omega, phi and psi torsion angle.

Parameters:

force (bool, optional) – If True the tag will be run even if Residues are already tagged.

valid_backbone_bond_angles(atol=20)[source]

True if all backbone bond angles are within atol degrees of their expected values.

Notes

Ideal bond angles taken from [1].

References

Parameters:

atol (float, optional) – Tolerance value in degrees for the absolute deviation away from ideal backbone bond angles.

valid_backbone_bond_lengths(atol=0.1)[source]

True if all backbone bonds are within atol Angstroms of the expected distance.

Notes

Ideal bond lengths taken from [1].

References

Parameters:

atol (float, optional) – Tolerance value in Angstoms for the absolute deviation away from ideal backbone bond lengths.

class ampal.protein.Residue(atoms=None, mol_code='UNK', monomer_id=' ', insertion_code=' ', is_hetero=False, parent=None)[source]

Bases: Monomer

Represents a amino acid Residue.

Parameters:
  • atoms (OrderedDict, optional) – OrderedDict containing Atoms for the Monomer. OrderedDict is used to maintain the order items were added to the dictionary.

  • mol_code (str, optional) – One or three letter code that represents the monomer.

  • monomer_id (str, optional) – String used to identify the residue.

  • insertion_code (str, optional) – Insertion code of monomer, used if reading from pdb.

  • is_hetero (bool, optional) – True if is a hetero atom in pdb. Helps with PDB formatting.

  • parent (ampal.Polypeptide, optional) – Reference to Polypeptide containing the Residue.

mol_code

PDB molecule code that represents the Residue.

Type:

str

insertion_code

Insertion code of Residue, used if reading from pdb.

Type:

str

is_hetero

True if is a hetero atom in pdb. Helps with PDB formatting.

Type:

bool

states

Contains an OrderedDicts containing atom information for each state available for the Residue.

Type:

dict

id

String used to identify the residue.

Type:

str

reference_atom

The key that corresponds to the reference atom. This is used by various functions, for example backbone primitives are calculated using the atom defined using this key.

Type:

str

parent

A reference to the Polypeptide containing this Residue.

Type:

Polypeptide or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

Raises:

ValueError – Raised if mol_code is not length 1 or 3.

property backbone

Returns a new Residue containing only the backbone atoms.

Returns:

bb_monomerResidue containing only the backbone atoms of the original Monomer.

Return type:

Residue

Raises:

IndexError – Raise if the atoms dict does not contain the backbone atoms (N, CA, C, O).

property centroid

Calculates the centroid of the residue.

Returns:

centroid – Returns a 3D coordinate for the residue unless a CB atom is not available, in which case None is returned.

Return type:

numpy.array or None

Notes

Uses the definition of the centroid from Huang et al [2].

References

property side_chain

List of the side-chain atoms (R-group).

Notes

Returns empty list for glycine.

Returns:

side_chain_atoms

Return type:

list(Atoms)

property unique_id

Generates a tuple that uniquely identifies a Monomer in an Assembly.

Notes

The unique_id will uniquely identify each monomer within a polymer. If each polymer in an assembly has a distinct id, it will uniquely identify each monomer within the assembly.

The hetero-flag is defined as in Biopython as a string that is either a single whitespace in the case of a non-hetero atom, or ‘H_’ plus the name of the hetero-residue (e.g. ‘H_GLC’ in the case of a glucose molecule), or ‘W’ in the case of a water molecule.

For more information, see the Biopython documentation or this Biopython wiki page: http://biopython.org/wiki/The_Biopython_Structural_Bioinformatics_FAQ

Returns:

unique_id – unique_id[0] is the polymer_id unique_id[1] is a triple of the hetero-flag, the monomer id (residue number) and the insertion code.

Return type:

tuple

ampal.protein.align(target, mobile, target_i=0, mobile_i=0)[source]

Aligns one Polypeptide (mobile) to another (target).

Notes

This function directly modifies atoms of the mobile Polypeptide! It does not return a new object.

Parameters:
  • target (Polypeptide) – Polypeptide to be aligned to.

  • mobile (Polypeptide) – Polypeptide to be moved during alignment.

  • target_i (int, optional) – Index of Residue in target to align to.

  • mobile_i (int, optional) – Index of Residue in mobile to be aligned.

ampal.protein.flat_list_to_dummy_chain(atom_list, atom_group_s=1)[source]

Converts flat list of coordinates into dummy C-alpha carbons

Parameters:
  • atom_list ([Atom]) – Flat list of co-ordinates.

  • atom_group_s (int, optional) – Size of atom groups.

Returns:

polymerPolymer object containing atom coord converted Monomers with ‘DUM’ atom name.

Return type:

Polypeptide

ampal.protein.flat_list_to_polymer(atom_list, atom_group_s=4)[source]

Takes a flat list of atomic coordinates and converts it to a Polymer.

Parameters:
  • atom_list ([Atom]) – Flat list of coordinates.

  • atom_group_s (int, optional) – Size of atom groups.

Returns:

polymerPolymer object containing atom coords converted Monomers.

Return type:

Polypeptide

Raises:

ValueError – Raised if atom_group_s != 4 or 5

ampal.pseudo_atoms module

Contains AMPAL objects representing pseudo atoms.

class ampal.pseudo_atoms.Primitive(monomers=None, polymer_id=' ', parent=None, sl=2)[source]

Bases: PseudoGroup

A backbone path composed of PseudoAtoms.

Parameters:
  • pseudo_atoms (OrderedDict, optional) – OrderedDict containing Atoms for the PsuedoMonomer. OrderedDict is used to maintain the order items were added to the dictionary.

  • mol_code (str, optional) – One or three letter code that represents the PsuedoMonomer.

  • monomer_id (str, optional) – String used to identify the PsuedoMonomer.

  • insertion_code (str, optional) – Insertion code of PsuedoMonomer, used if reading from pdb.

  • is_hetero (bool, optional) – True if is a hetero atom in pdb. Helps with PDB formatting.

  • parent (ampal.PseudoGroup, optional) – Reference to PseudoGroup containing the PsuedoMonomer.

mol_code

PDB molecule code that represents the Nucleotide.

Type:

str

insertion_code

Insertion code of Nucleotide, used if reading from pdb.

Type:

str

is_hetero

True if is a hetero atom in pdb. Helps with PDB formatting.

Type:

bool

states

Contains an OrderedDicts containing atom information for each state available for the Nucleotide.

Type:

dict

id

String used to identify the Nucleotide.

Type:

str

reference_atom

The key that corresponds to the reference Atom. This is used by various functions, for example backbone primitives are calculated using the Atom defined using this key.

Type:

str

parent

A reference to the Polynucleotide containing this Nucleotide.

Type:

Polynucleotide or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

Raises:

ValueError – Raised if mol_code is not length 1 or 3.

property coordinates

Returns the backbone coordinates for the Primitive.

classmethod from_coordinates(coordinates)[source]

Creates a Primitive from a list of coordinates.

radii_of_curvature()[source]

The radius of curvature at each point on the Polymer primitive.

Notes

Each element of the returned list is the radius of curvature, at a point on the Polymer primitive. Element i is the radius of the circumcircle formed from indices [i-1, i, i+1] of the primitve. The first and final values are None.

rise_per_residue()[source]

The rise per residue at each point on the Primitive.

Notes

Each element of the returned list is the rise per residue, at a point on the Primitive. Element i is the distance between primitive[i] and primitive[i + 1]. The final value is None.

class ampal.pseudo_atoms.PseudoAtom(coordinates, name='', occupancy=1.0, bfactor=1.0, charge=' ', parent=None)[source]

Bases: Atom

Object containing 3D coordinates and name.

Notes

Used to represent pseudo atoms (e.g. centre_of_mass) in ISAMBARD.

Parameters:
  • coordinates (3D Vector (tuple, list, numpy.array)) – Position of PseudoAtom in 3D space.

  • element (str) – Element of PseudoAtom.

  • atom_id (str) – Identifier for PseudoAtom, usually a number.

  • res_label (str, optional) – Label used in Monomer to refer to the PseudoAtom type i.e. “CA” or “OD1”.

  • occupancy (float, optional) – The occupancy of the PseudoAtom.

  • bfactor (float, optional) – The bfactor of the PseudoAtom.

  • charge (str, optional) – The point charge of the PseudoAtom.

  • state (str) – The state of this PseudoAtom. Used to identify PseudoAtoms with a number of conformations.

  • parent (ampal.Monomer, optional) – A reference to the Monomer containing this PseudoAtom.

id

Identifier for PseudoAtom, usually a number.

Type:

str

res_label

Label used in PseudoGroup to refer to the Atom type i.e. “CA” or “OD1”.

Type:

str

element

Element of Atom.

Type:

str

parent
A reference to the PseudoGroup containing this PseudoAtom.

number of conformations.

Type:

ampal.PseudoAtom

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

class ampal.pseudo_atoms.PseudoGroup(monomers=None, polymer_id=' ', parent=None, sl=2)[source]

Bases: Polymer

Container for PseudoMonomer, inherits from Polymer.

Parameters:
  • monomers (PseudoAtom or [PseudoGroup], optional) – PseudoMonomer or list containing PseudoMonomer objects to form the PseudoGroup.

  • polymer_id (str, optional) – An ID that the user can use to identify the PseudoGroup. This is used when generating a pdb file using PseudoGroup().pdb.

  • parent (ampal.Assembly, optional) – Reference to Assembly containing the PseudoGroup.

  • sl (int, optional) – The default smoothing level used when calculating the backbone primitive.

id

PseudoGroup ID

Type:

str

parent

Reference to Assembly containing the PseudoGroup

Type:

ampal.Assembly or None

molecule_type

A description of the type of Polymer i.e. Protein, DNA etc.

Type:

str

ligands

A LigandGroup containing all the Ligands associated with this PseudoGroup chain.

Type:

ampal.LigandGroup

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

sl

The default smoothing level used when calculating the backbone primitive.

Type:

int

Raises:

TypeErrorPolymer type objects can only be initialised empty or using a Monomer.

class ampal.pseudo_atoms.PseudoMonomer(pseudo_atoms=None, mol_code='UNK', monomer_id=' ', insertion_code=' ', parent=None)[source]

Bases: Monomer

Represents a collection of PsuedoAtoms.

Parameters:
  • pseudo_atoms (OrderedDict, optional) – OrderedDict containing Atoms for the PsuedoMonomer. OrderedDict is used to maintain the order items were added to the dictionary.

  • mol_code (str, optional) – One or three letter code that represents the PsuedoMonomer.

  • monomer_id (str, optional) – String used to identify the PsuedoMonomer.

  • insertion_code (str, optional) – Insertion code of PsuedoMonomer, used if reading from pdb.

  • is_hetero (bool, optional) – True if is a hetero atom in pdb. Helps with PDB formatting.

  • parent (ampal.PseudoGroup, optional) – Reference to PseudoGroup containing the PsuedoMonomer.

mol_code

PDB molecule code that represents the Nucleotide.

Type:

str

insertion_code

Insertion code of Nucleotide, used if reading from pdb.

Type:

str

is_hetero

True if is a hetero atom in pdb. Helps with PDB formatting.

Type:

bool

states

Contains an OrderedDicts containing atom information for each state available for the Nucleotide.

Type:

dict

id

String used to identify the Nucleotide.

Type:

str

reference_atom

The key that corresponds to the reference Atom. This is used by various functions, for example backbone primitives are calculated using the Atom defined using this key.

Type:

str

parent

A reference to the Polynucleotide containing this Nucleotide.

Type:

Polynucleotide or None

tags

A dictionary containing information about this AMPAL object. The tags dictionary is used by AMPAL to cache information about this object, but is also intended to be used by users to store any relevant information they have.

Type:

dict

Raises:

ValueError – Raised if mol_code is not length 1 or 3.

property pdb

Generates a PDB string for the PseudoMonomer.

Module contents