The Molecular Entity
This document describes the EntityHandle and the related classes.
The Handle Classes
-
CreateEntity ()
Creates a new entity. The created entity is empty, that is, it does not
contain any atoms, residues, chains, bonds or torsions. To populate the
entity, use an entity editor.
-
class
EntityHandle
The entity class represents a molecular structure. Such a structure is in
general made up of one or more chains of residues, which in turn are
formed by one or more atoms.
The interface of entities is tailored to biological macromolecules, but this
does not prevent it to be used for molecules in general: An entity also
represent a ligand or a collection of water molecules - hence the rather
generic name.
-
chains
List of all chains of this entity. The chains are in the same
order they have been added, for example the order they have been listed
in a PDB file.
Also available as GetChainList() . To access a single chain, use
FindChain() .
This property is read-only.
-
chain_count
Number of chains. Read-only. See GetChainCount() .
-
residues
List of all residues of this entity. The returned residues are first
sorted by chain and then from N- to C-terminus.
This property is read-only.
Example usage:
for res in ent.residues:
print(res.name, res.atom_count)
Also available as GetResidueList() . To access a single residue, use
FindResidue() .
-
residue_count
Number of residues. Read-only. See GetResidueCount() .
-
atoms
Get list of all atoms of this entity. To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
atom_count
Number of atoms. Read-only. See GetAtomCount() .
-
bounds
Axis-aligned bounding box of the entity. Read-only
-
mass
The total mass of the entity in Dalton. Also available as GetMass() .
-
center_of_atoms
Center of atoms, that is the average atom position of the entity.
Use center_of_mass for the the mass-weighted center of the entity.
Also available as GetCenterOfAtoms() .
-
center_of_mass
Center of mass of the entity.
Also available as GetCenterOfMass()
-
geometric_center
Mid-point of the axis aligned bounding box of the entity.
Also available as GetGeometricCenter()
-
positions
Equivalent to calling GetPositions() with sort_by_index = True. This
property is read-only and only available if OpenStructure was compiled with
an enabled USE_NUMPY flag (see here for details).
-
valid
Validity of handle.
-
GetName ()
Returns: | Name associated to this entity. |
Return type: | str |
-
SetName (name)
Parameters: | name (str ) – Sets this as new name to be associated to this entity. |
-
FindChain (chain_name)
Get chain by name. See also chains
Parameters: | chain_name (str) – Chain identifier, e.g. “A” |
Returns: | A valid ChainHandle , if the entity contains a
chain with the given name, an invalid
ChainHandle otherwise. |
-
GetChainList ()
See chains
-
GetChainCount ()
See chain_count
-
FindResidue (chain_name, res_num)
Get residue by chain name and residue number. See also
GetResidueList()
Parameters: |
- chain_name (str) – Chain identifier, e.g. “A”
- res_num (
ResNum ) – residue number
|
Returns: | A valid ResidueHandle if the chain exists and
the chain contains a residue of the given residue
number, an invalid ResidueHandle otherwise.
|
-
GetResidueList ()
See residues
-
GetResidueCount ()
See residue_count
-
FindAtom (chain_name, res_num, atom_name)
Get atom by chain name, residue number and atom name. See also
atoms
Parameters: |
- chain_name (str) – Chain identifier, e.g. “A”
- res_num (
ResNum ) – residue number
- atom_name (str) – atom name, e.g. CA
|
Returns: | A valid AtomHandle if an atom matching the
parameters could be found, an invalid
AtomHandle otherwise
|
-
GetAtomList ()
See atoms
-
GetAtomCount ()
See atom_count
-
EditXCS ([edit_mode=mol.EditMode.UNBUFFERED_EDIT])
Request XCSEditor for editing the external coordinate system. This
call will fail when there are pending changes of the internal coordinate
system.
Parameters: | edit_mode (mol.EditMode) – Must be EditMode.BUFFERED_EDIT or
EditMode.UNBUFFERED_EDIT. For more details, see the
editor documentation. |
Returns: | XCSEditor |
-
EditICS ([edit_mode=mol.EditMode.UNBUFFERED_EDIT])
Request ICSEditor for editing the internal coordinate system, such
as torsions, bond lengths and angle between two bonds. This call will fail
when there are pending changes of the external coordinate system.
Parameters: | edit_mode (mol.EditMode) – Must be EditMode.BUFFERED_EDIT or
EditMode.UNBUFFERED_EDIT. For more details, see the
editor documentation. |
Returns: | ICSEditor |
-
Select (query, flags=0)
Perform a selection on the entity. The result of the selection is an
EntityView which (usually) contains only a subset of chains,
residues, atoms and bonds of the original entity.
Example Usage:
# select calpha atoms of peptides
calphas = ent.Select('aname=CA and peptide=true')
# select atoms in a box of size 10, centred at the origin
in_box = ent.Select('x=-5:5 and y=-5:5 and z=-5:5')
Parameters: |
- query (
Query / str ) – The query to be executed.
- flags (
int / QueryFlag ) – An ORed together combination of QueryFlag
|
Returns: | An entity view .
|
Raises: | QueryError when the query could not be executed due
to syntactic errors.
|
-
CreateFullView ()
Creates an entity view containing all chains, residues, atoms and bonds of
this entity.
# these two lines are identical
full=ent.Select('')
full=ent.CreateFullView()
-
CreateEmptyView ()
Creates an entity view pointing to this entity, but otherwise empty. This
method is usually the starting point for manual entity view creation, e.g.:
view=ent.CreateEmtpyView()
for atom in ent.atoms:
if ComplicatedPredicate(atom):
view.AddAtom(atom)
-
Copy ()
Creates a deep copy of the entity.
Returns: | A new EntityHandle that is an exact copy of this entity. |
Note
alternative atom positions are not handled yet.
-
GetMass ()
See mass
-
GetCenterOfAtoms ()
See center_of_atoms
-
GetCenterOfMass ()
See center_of_mass
-
GetGeometricCenter ()
See geometric_center
-
GetPositions (sort_by_index=True)
Returns: | Array of atom positions for this entity. |
Return type: | numpy.array (shape [atom_count , 3]) |
Parameters: | sort_by_index – If True, the atoms are sorted by their
index . Otherwise, they are sorted
as they appear in the atoms list. |
This method is only available if OpenStructure was compiled with an enabled
USE_NUMPY flag (see here for details).
-
FindWithin (pos, radius)
Find all atoms in sphere of given radius centred at pos. To turn the
returned list of atoms into an EntityView , use
CreateViewFromAtomList() .
Parameters: |
- pos (
Vec3 ) – Center of sphere
- radius (float) – The radius of the sphere
|
Returns: | AtomHandleList (list of AtomHandle )
|
-
IsValid ()
See valid
-
class
ChainHandle
A chain of one or more residues . Chains are always
part of an entity.
-
atoms
Get list of all atoms of this chain. To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
bounds
Axis-aligned bounding box of the chain. Read-only
-
center_of_atoms
Center of atoms (not mass weighted). Also available as
GetCenterOfAtoms() .
-
center_of_mass
Center of mass. Also available as GetCenterOfMass() .
-
geometric_center
Mid-point of the axis aligned bounding box of the chain.
Also available as GetGeometricCenter()
-
description
Details about the chain. Not categorised, just text.
-
in_sequence
Whether the residue numbers are in ascending order. For example:
chain=ent.FindChain("A")
print(chain.residues) # [A.GLY1, A.GLY2, A.GLY4A, A.GLY4B]
print(chain.in_sequence) # prints true
chain=ent.FindChain("B")
print(chain.residues) # [B.GLY1, B.GLY4, B.GLY3]
print(chain.in_sequence) # prints false
-
is_oligosaccharide
Indicates if the chain is an oligosaccharide, a branched, non-linear entity
of multiple sugars. Also available as IsOligosaccharide() .
-
is_polymer
Indicates if a chain is a polymer. True for polypeptides, polynucleotides,
polysaccharides, oligosaccharides and branched chains. Also available as
IsPolymer() .
-
is_polynucleotide
Indicates if a chain is a nucleic acid. Also available as
IsPolynucleotide() .
-
is_polypeptide
Indicates if a chain is a protein. Also available as IsPolypeptide() .
-
is_polysaccharide
Indicates if a chain is a polysaccharide. Also available as
IsPolysaccharide() .
-
mass
The total mass of this chain in Dalton. Also available as GetMass()
-
name
The chain name. The name uniquely identifies the chain in the entity. In
most cases, the chain name is one character. This is restriction of the PDB
file format. However, you are free to use longer names as long as you don’t
want to save them as PDB files
This property is read-only. To change the name, use an XCSEditor .
Also available as GetName()
-
residue_count
Number of residues. Read-only. See GetResidueCount() .
-
residues
List of all residues of this chain. The residues are sorted from N- to
C-terminus. Usually the residue numbers are in ascending order
(see in_sequence ).
This property is read-only.
Example usage:
chain=ent.FindChain("A")
for res in chain.residues:
print(res.name, res.atom_count)
Also available as GetResidueList() . To access a single residue, use
FindResidue() .
-
type
Describes the type of the chain.
-
hash_code
A unique identifier for this chain. Note that a deep copy of an entity (see
EntityHandle.Copy() ) will have chains with differing identifiers.
Shallow copies of the entity preserve the identifier.
Also available as GetHashCode() .
-
valid
Validity of handle.
-
FindResidue (res_num)
Get residue by residue number. See also residues .
Parameters: | res_num (ResNum ) – residue number |
Returns: | A valid ResidueHandle if the chain contains
a residue with matching residue number, an invalid
ResidueHandle otherwise. |
-
GetResidueList ()
See residues
-
GetResidueCount ()
See residue_count
-
FindAtom (res_num, atom_name)
Get atom by residue number and atom name. See also atoms
Parameters: |
- res_num (
ResNum ) – residue number
- atom_name (str) – atom name, e.g. CA
|
Returns: | A valid AtomHandle if an atom matching the
parameters could be found, an invalid
AtomHandle otherwise
|
-
GetAtomList ()
See atoms
-
GetName ()
See name
-
GetType ()
See type
-
IsOligosaccharide ()
See is_oligosaccharide
-
IsPolymer ()
See is_polymer
-
IsPolynucleotide ()
See is_polynucleotide
-
IsPolypeptide ()
See is_polypeptide
-
IsPolysaccharide ()
See is_polysaccharide
-
GetDescription ()
See description
-
GetHashCode ()
See hash_code
-
IsValid ()
See valid
-
class
ResidueHandle
The residue is either used to represent complete molecules or building blocks
in a polymer, e.g. in a protein, DNA or RNA. A residue consists of one or
more atoms . Residues are always part of a
ChainHandle , even if they are ligands or water molecules where the
concept of a chain does not apply.
-
name
The residue name is usually a str of 3 characters, e.g. GLY for
glycine or ALA for alanine, but may be shorter, e.g. G for guanosine,
or longer for structures loaded from formats other than PDB.
Also available as GetName() .
This property is read-only. To change the name of the residue, use
RenameResidue() .
-
qualified_name
The qualified name consists of a residue identifier and chain name.
For a glycine with residue number 2 of chain A, the qualified name is
“A.GLY2”. Also available as GetQualifiedName() .
-
number
The number of this residue. The residue number has a numeric part and an
insertion-code. This property is read-only. Also available as
GetNumber() .
-
one_letter_code
For amino acids, and nucleotides the one_letter_code is an alpha-numeric
character. For unknown or more exotic residues, the one letter code is set
to ‘?’.
Example
This code-snippet shows how to get the sequence string from a list of
residues.
print(''.join([r.one_letter_code for r in chain.residues]))
-
atoms
Get list of all atoms of this residue. To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
atom_count
Number of atoms. Read-only. See GetAtomCount() .
-
bond_count
Number of bonds. Read-only. See GetBondCount() .
-
bounds
Axis-aligned bounding box of the residue. Read-only.
-
mass
The total mass of this residue in Dalton. Also available as GetMass() .
-
center_of_atoms
Center of atoms (not mass weighted). Also available as
GetCenterOfAtoms() .
-
center_of_mass
Center of mass. Also available as GetCenterOfMass()
-
geometric_center
Mid-point of the axis aligned bounding box of the residue.
-
chain
The chain this residue belongs to. Read-only. Also available as
GetChain()
-
phi_torsion
The PHI dihedral angle between this residue and the next. For residues
that are not amino acids, residues that do not have all atoms required or
residues that do not have bonds between the four atoms involved in the
torsion, the PHI torsion is an invalid handle.
Read-only. Also available as GetPhiTorsion()
-
psi_torsion
The PSI dihedral angle between this residue and the previous. For residues
that are not amino acids, residues that do not have all atoms required or
residues that do not have bonds between the four atoms involved in the
torsion, the PSI torsion is an invalid handle.
Read-only. Also available as GetPsiTorsion()
-
omega_torsion
The OMEGA dihedral angle between this residue and the previous. For residues
that are not amino acids, residues that do not have all atoms required or
residues that do not have bonds between the four atoms involved in the
torsion, the OMEGA torsion is an invalid handle.
Read-only. Also available as GetOmegaTorsion()
-
chem_class
The chemical class of the residue.
-
chem_type
The chemical type of the residue. The type is only properly set if a
compound library is used.
-
sec_structure
The secondary structure of the residue.
-
is_ligand
Whether the residue is a ligand. When loading PDB structures, this property
is set based on the HET records. This also means, that this property will
most likely not be set properly for all except PDB files coming from
pdb.org. When loading MMCIF structures, this property is set based on
seqres information and differs from PDB structures.
Also available as IsLigand() , SetIsLigand() .
-
is_protein
Whether the residue is considered to be part of a protein. This is set when
loading a structure if the residue forms a feasible peptide bond to the
previous or next residue (see IsBondFeasible() ). Also
available as IsProtein() , SetIsProtein() . In contrast to
IsPeptideLinking() this excludes residues which are not connected to
neighbouring residues such as CA-only residues or badly positioned ones.
-
peptide_linking
Whether residue can form peptide bonds. This is determined based on
chem_class which is set when loading the structure.
-
nucleotide_linking
Whether residue can form nucleotide bonds. This is determined based on
chem_class which is set when loading the structure.
-
index
Residue index (starting at 0) within chain.
-
central_atom
Central atom used for rendering traces. For peptides, this is usually
the CA atom. For nucleotides, this is usually the P atom.
-
central_normal
Normal computed for central_atom . Only defined for peptides and
nucleotides if all required atoms available. Otherwise, the (1,0,0) vector
is returned.
-
hash_code
A unique identifier for this residue. Note that a deep copy of an entity (see
EntityHandle.Copy() ) will have residues with differing identifiers.
Shallow copies of the entity preserve the identifier. Also available as
GetHashCode() .
-
valid
Validity of handle.
-
next
Residue after this one in the same chain. Invalid handle returned if there
is no next residue. Residues are ordered as in ChainHandle.residues
independently on whether they are connected or not (see
InSequence() to check for connected residues).
-
prev
Residue before this one in the same chain. Otherwise same behaviour as
next .
-
FindAtom (atom_name)
Get atom by atom name. See also atoms
Parameters: | atom_name (str) – atom name, e.g. CA |
Returns: | A valid AtomHandle if an atom with the given
name could be found, an invalid AtomHandle
otherwise |
-
HasAltAtoms ()
Returns whether the residue has any alternative atom groups
-
HasAltAtomGroup (group_name)
Returns whether the residue has a certain alternative atom group
Parameters: | group_name (str ) – Group of interest |
Return type: | bool |
-
GetAltAtomGroupNames ()
Returns names of all alternative atom groups in residue
-
GetCurrentAltGroupName ()
Returns the currently active alternative atom group, empty
str if residue has no alternative atom groups
-
SwitchAtomPos (group_name)
Switches the atoms in residue to the specified alternative atom group
Parameters: | group_name (str ) – Group of interest |
Return type: | bool |
Returns: | Whether the switch was successful (e.g. False if no such
group exists) |
-
GetName ()
See name
-
GetQualifiedName ()
See qualified_name
-
GetNumber ()
See number
-
GetOneLetterCode ()
See one_letter_code
-
GetAtomList ()
See atoms
-
GetAtomCount ()
See atom_count
-
GetBondCount ()
See bond_count
-
GetBounds ()
See bounds
-
GetMass ()
See mass
-
GetCenterOfAtoms ()
See center_of_atoms
-
GetCenterOfMass ()
See center_of_mass
-
GetGeometricCenter ()
See geometric_center
-
GetChain ()
See chain
-
GetPhiTorsion ()
See phi_torsion
-
GetPsiTorsion ()
See psi_torsion
-
GetOmegaTorsion ()
See omega_torsion
-
GetChemClass ()
See chem_class
-
GetChemType ()
See chem_type
-
GetSecStructure ()
See sec_structure
-
IsLigand ()
See is_ligand
-
IsProtein ()
See is_protein
-
IsPeptideLinking ()
See peptide_linking
-
IsNucleotideLinking ()
See nucleotide_linking
-
GetIndex ()
See index
-
GetCentralAtom ()
-
SetCentralAtom ()
See central_atom
-
GetCentralNormal ()
See central_normal
-
GetHashCode ()
See hash_code
-
IsValid ()
See valid
-
SetIsLigand ()
Set the IsLigand() flag explicitly.
Parameters: | ligand (bool) – Whether this residue is a ligand or not |
-
SetIsProtein ()
Set the IsProtein() flag explicitly.
Parameters: | protein (bool) – Whether this residue is a protein or not |
-
GetNext ()
See next
-
GetPrev ()
See next
-
class
AtomHandle
Represents an atom in a molecular structure. Atoms are always part of a
residue.
-
name
The atom name, usually a capitalized string consisting of the element and
an position indicator, e.g. G1. In general, the atom name uniquely
identifies an atom within a residue. However, this is not enforced and there
may be more than one atoms with the same name.
-
qualified_name
The qualified name consists of the atom name as well as a residue
identifier and chain name. For CA of a glycine with residue number 2 of
chain A, the qualified name is “A.GLY2.CA”.
-
element
The atom’s element. By convention in Openstructure, this is the chemical
symbol in uppercase, but this is not strictly enforced and may be a non-
existing element or an empty string. Also available as GetElement() .
Read-only.
-
mass
The atom’s mass in Dalton. Also available as GetMass() . Read-only.
-
pos
The atom’s position in global coordinates, transformed by the entity
transformation. Also available as GetPos() .
Read-only.
-
original_pos
The atom’s untransformed position in global coordinates. Also available as
GetOriginalPos() .
Read-only.
-
radius
The van-der-Waals radius of the atom. Also available as GetRadius() .
Read/write.
-
occupancy
The atom’s occupancy in the range 0 to 1. Read/write. Also available as
GetOccupancy() , SetOccupancy() .
-
b_factor
The atom’s temperature factor. Read/write. Also available as
GetBFactor() , SetBFactor() .
-
charge
The atom’s charge. Also available as GetCharge() , SetCharge() .
-
residue
The residue this atom belongs to. Read-only.
-
is_hetatom
Indicates whether this atom is a HETATM. When loading a structure from a
PDB file, all non-standard aminoacid and nucleotide atoms as well as ligands
are marked as HETATM.
-
bonds
List of bonds this atom is involved in. The bonds are in no particular
order. Also available as GetBondList() .
-
index
Atom index (starting at 0) within entity.
-
hash_code
A unique identifier for this atom. Note that a deep copy of an entity (see
EntityHandle.Copy() ) will have atoms with differing identifiers.
Shallow copies of the entity preserve the identifier. Atom views on a handle
have different identifiers, but the atom view handles (see
AtomView.handle ) have the same identifier. Also available
as GetHashCode() .
-
valid
Validity of handle.
-
FindBondToAtom (other_atom)
Finds and returns the bond formed between this atom and other_atom. If no
bond exists between the two atoms, an invalid BondHandle is
returned.
-
GetBondCount ()
-
-
GetBondList ()
See bonds
-
GetBondPartners ()
Get list of atoms this atom is bonded with. This method exists as a shortcut
to determine all the bonding partners of an atom and avoids code like this:
bond_parters=[]
for bond in atom.bonds:
if bond.first==atom:
bond_partners.append(bond.second)
else:
bond_partners.append(bond.first)
Return type: | AtomHandleList |
-
GetCharge ()
See charge
-
GetElement ()
See element
-
GetEntity ()
The entity this atom belongs to
-
GetHashCode ()
See hash_code
-
GetIndex ()
See index
-
GetMass ()
See mass
-
GetName ()
See name
-
GetOriginalPos ()
-
-
GetPos ()
See pos
-
GetQualifiedName ()
See qualified_name
-
GetRadius ()
See radius
-
GetResidue ()
See residue
-
IsHetAtom ()
See is_hetatom
-
IsValid ()
See valid
-
class
BondHandle
Represents a chemical bond between two atoms (first and second).
-
first
-
second
Atoms involved in the bond. No assumptions about the order should be made.
With the internal coordinate system enabled, first and second may even be
swapped when rebuilding the internal connectivity tree. Also available as
GetFirst() and GetSecond() .
-
pos
Midpoint between the two atoms (transformed coordinates). Also available as
GetPos() .
-
length
Length of the bond. Also available as GetLength() .
-
bond_order
The bond order. Possible values:
1 - single bond
2 - double bond
3 - triple bond
4 - aromatic bond
Also available as GetBondOrder() .
-
hash_code
A unique identifier for this bond handle. Also available as GetHashCode() .
-
valid
Validity of handle. Also available as IsValid() .
-
GetFirst ()
See first
-
GetSecond ()
See second
-
GetPos ()
See pos
-
GetLength ()
See length
-
GetBondOrder ()
See bond_order
-
GetOther (other_atom)
Get the other atom. Returns the one of the two atoms that does not match
the given one.
-
SetBondOrder (order)
Set the bond order. See GetBondOrder() .
Parameters: | order (int ) – The bond order |
See bond_order
-
IsValid ()
See valid
-
GetHashCode ()
See hash_code
The View Classes
-
class
EntityView
An entity view represents a structural subset of an EntityHandle . For
an introduction, see Introduction to the mol Module.
-
handle
The underlying handle of the entity view. Also
available as GetHandle() .
-
chains
List of all chains of this entity. The chains are in the same
order they have been added.
Also available as GetChainList() . To access a single chain, use
FindChain() .
This property is read-only.
-
chain_count
Number of chains. Read-only. See GetChainCount() .
-
residues
List of all residues of this entity. The returned residues are first
sorted by chain and then from N- to C-terminus.
This property is read-only.
Example usage:
for res in ent.residues:
print(res.name, res.atom_count)
Also available as GetResidueList() . To access a single residue, use
FindResidue() .
-
residue_count
Number of residues. Read-only. See GetResidueCount() .
-
atoms
Get list of all atoms of this entity. To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
atom_count
Number of atoms. Read-only. See GetAtomCount() .
-
bounds
See EntityHandle.bounds
-
mass
See EntityHandle.mass
-
center_of_atoms
See EntityHandle.center_of_atoms
-
center_of_mass
See EntityHandle.center_of_mass
-
geometric_center
See EntityHandle.geometric_center
-
valid
Validity of view.
-
GetName ()
-
-
SetName (name)
-
-
CreateEmptyView ()
See EntityHandle.CreateEmptyView()
-
CreateFullView ()
Returns a copy of this entity. Provided for duck-typing <http://en.wikipedia.org/wiki/Duck_typing> purposes.
-
AddChain (chain_handle[, view_add_flags])
Add chain to view. By default, only the chain is added to the view, but not
its residues and atoms. This behaviour can be changed by passing in an
appropriate set of view_add_flags.
The following example just adds a chain without residues and atoms:
pdb = ost.io.LoadPDB(<PDB file name>)
v = pdb.CreateEmptyView()
v.AddChain(pdb.chains[0])
To get the chain with residues and atoms, go like:
pdb = ost.io.LoadPDB(<PDB file name>)
v = pdb.CreateEmptyView()
v.AddChain(pdb.chains[0], ost.mol.INCLUDE_ALL)
Note that the view above still lacks bonds which can be added with the
AddAllInclusiveBonds() method.
-
AddResidue (residue_handle[, view_add_flags])
Add residue to view. If the residue’s chain is not already part of the
view, it will be added. By default, only the residue is added, but not its
atoms. This behaviour can be modified by passing in an appropriate
combination of ViewAddFlag .
-
AddAtom (atom_handle[, view_add_flags])
Add the atom to view. The chain and residue of the atom are added to the
view if they haven’t already been added.
-
AddBond (bond_handle)
Add bond to the view.
Parameters: | bond_handle (BondHandle) – Bond to be added |
Returns: | True if the bond has been added, false if the bond was already
present in the view. |
Return type: | bool |
-
AddAllInclusiveBonds ()
Convenience method to add all bonds to the view for which both of the bond
partners are in the view. This method is useful when manually assembling
views.
-
RemoveChain (chain)
Remove chain and all its residues from the entity view.
Parameters: | chain (ChainView) – The chain to be removed. May be invalid |
-
RemoveResidue (residue)
Remove residue from view.
Parameters: | residue (ResidueView) – The residue to be removed. May be invalid |
-
RemoveAtom (atom)
Remove atom from view
Parameters: | atom (AtomView) – The atom to be removed. May be invalid |
-
GetAngle (atom1, atom2, atom3)
Parameters: |
|
Return type: | float
|
-
FindWithin (pos, radius)
Find all atoms that are within radius of the given position. See
EntityHandle.FindWithin() .
Parameters: |
- pos (
Vec3 ) – Center of sphere
- radius (float) – The radius of the sphere
|
Returns: | AtomHandleList (list of AtomHandle )
|
-
FindChain (chain_name)
Find chain by name.
Parameters: | chain_name (str) – Chain identifier, e.g. “A” |
Returns: | The chain if present in the view, an invalid ChainView
otherwise |
Return type: | ChainView |
-
FindResidue (residue)
Find residue view of pointing to the given handle.
Parameters: | residue (ResidueHandle) – Residue handle |
Returns: | The residue view pointing the the handle, or an invalid handle if the residue is not part of the view |
Return type: | ResidueView |
-
FindAtom (chain_name, res_num, atom_name)
Parameters: |
- chain_name (str) – The chain name
- res_num (
ResNum or int ) – The residue number
- atom_name (str) – The name of the atom
|
Return type: | AtomView
|
-
Select (query, flags=0)
Perform selection on entity view. See EntityHandle.Select() .
-
Copy ()
Returns a deep copy of the entity view. The effect is identical to calling
-
GetMass ()
Get total mass of view.
-
GetCenterOfMass ()
Get the center of mass. For a non-mass-weighted center, see
GetCenterOfAtoms() .
-
GetGeometricCenter ()
Get the geometric center (the center of the axis-aligned
bounding box).
-
GetGeometricStart ()
-
-
GetCenterOfAtoms ()
Calculates the center of atoms. For a mass-weighted center, use
GetCenterOfMass() .
-
GetBondCount ()
Get number of bonds
-
GetBondList ()
See bonds
-
GetHandle ()
See handle
-
GetGeometricEnd ()
-
-
GetChainList ()
See chains
-
GetChainCount ()
See chain_count
-
GetResidueList ()
See residues
-
GetResidueCount ()
See residue_count
-
GetAtomList ()
See atoms
-
GetAtomCount ()
See atom_count
-
IsValid ()
See valid
-
class
ChainView
A view representation of a ChainHandle . Mostly, the same
functionality is provided as for the handle.
-
handle
The chain handle this view points to. Also available as GetHandle() .
-
hash_code
A unique identifier for this chain view. Note, that this is not the same as
for the chain handle (see ChainHandle.hash_code ).
Also available as GetHashCode() .
-
name
The chain name. The name uniquely identifies the chain in the entity. In
most cases, the chain name is one character. This is a restriction of the
PDB file format. However, you are free to use longer names as long as you
don’t want to save them as PDB files.
This property is read-only. To change the name, use an XCSEditor .
Also available as GetName()
-
residues
List of all residues of this chain. The residues are sorted from N- to
C-terminus. Usually the residue numbers are in ascending order
(see in_sequence ).
This property is read-only.
Example usage:
chain=ent.FindChain("A")
for res in chain.residues:
print(res.name, res.atom_count)
Also available as GetResidueList() . To access a single residue, use
FindResidue() .
-
in_sequence
Whether the residue numbers are in ascending order. For example:
chain=ent.FindChain("A")
print(chain.residues) # [A.GLY1, A.GLY2, A.GLY4A, A.GLY4B]
print(chain.in_sequence) # prints true
chain=ent.FindChain("B")
print(chain.residues) # [B.GLY1, B.GLY4, B.GLY3]
print(chain.in_sequence) # prints false
Note that the value of in_sequence is independent from the value of
ChainHandle.in_sequence .
-
atoms
Get list of all atoms of this chain. To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
bounds
Axis-aligned bounding box of the chain. Read-only
-
mass
The total mass of this chain in Dalton. Also available as GetMass()
-
center_of_atoms
See ChainHandle.center_of_atoms
-
center_of_mass
See ChainHandle.center_of_mass
-
geometric_center
See ChainHandle.geometric_center
-
valid
Validity of view.
-
AddAtom (atom_handle[, view_add_flags])
Add atom to the view. If the residue of the atom is not already part of the
view, it will be added. If the atom does not belong to chain, the result is
undefined.
-
AddResidue (residue_handle[, view_add_flags])
Add residue to the view. If the atom does not belong to chain, the result is
undefined. By default, only the residue, but no atoms are added to the
view. To change the behavior, pass in a suitable combination of
ViewAddFlag .
-
FindAtom (res_num, atom_name)
Find atom with the given residue number and atom name.
Parameters: |
- res_num (
ResNum ) – The residue number
- atom_name (str) – The atom name
|
Returns: | The atom view, or an invalid atom view if no atom with the given
parameters is part of the view.
|
Return type: | AtomView
|
-
FindResidue (res_num)
Find residue by residue number.
Parameters: | res_num (ResNum ) – |
Return type: | ResidueView |
Returns: | The residue view, or an invalid residue view if no residue with
the given residue number is in the view. |
-
GetCenterOfAtoms ()
See center_of_atoms
-
GetCenterOfMass ()
See center_of_mass
-
GetEntity ()
The parent entity.
-
GetGeometricCenter ()
See geometric_center
-
GetHandle ()
See handle
-
GetHashCode ()
See hash_code
-
GetMass ()
See mass
-
GetName ()
See name
-
GetResidueByIndex (index)
Returns the residue view at index. If index is negative or bigger than
the number of residues minus one, an invalid ResidueView is returned.
Parameters: | index (int) – The index |
Return type: | ResidueView |
-
GetResidueList ()
See residues
-
InSequence ()
See in_sequence
-
IsValid ()
See valid
-
RemoveResidue (residue)
Remove residue from chain. If the residue is not part of the view, the
chain is left unchanged.
Parameters: | residue (ResidueView ) – The residue view. May be invalid |
Return type: | None |
-
RemoveResidues ()
Remove all residues from this chain view
-
Select (query, flags=0)
Perform query on chain view. See EntityHandle.Select() .
-
class
ResidueView
A view representation of a ResidueHandle . Mostly, the same
functionality is provided as for the handle.
-
handle
The residue handle this view points to. Also available as
GetHandle() .
-
hash_code
A unique identifier for this residue view. Note, that this is not the same as
for the residue handle (see ResidueHandle.hash_code ).
Also available as GetHashCode() .
-
name
-
qualified_name
-
number
-
one_letter_code
-
bounds
-
mass
-
center_of_atoms
-
center_of_mass
-
geometric_center
-
phi_torsion
-
psi_torsion
-
omega_torsion
-
chem_class
-
chem_type
-
sec_structure
-
is_ligand
-
is_protein
-
peptide_linking
-
nucleotide_linking
-
central_atom
-
central_normal
-
valid
-
next
-
prev
See the respective attributes in ResidueHandle .
-
atoms
Get list of all atoms of this residue included in the view.
To access a single atom, use
FindAtom() .
This property is read-only. Also available as GetAtomList()
-
atom_count
Number of atoms included in the view. Read-only. See GetAtomCount() .
-
chain
The chain this residue belongs to. Read-only. Also available as
GetChain()
-
index
Residue index (starting at 0) within chain view.
-
RemoveAtom (atom_view)
Remove atom from residue and all associated bonds. If the atom is not part
of the view, the residue view is left unchanged.
Parameters: | atom_view (AtomView ) – The atom to be removed. May be an invalid handle |
Return type: | None |
-
GetHandle ()
See handle
-
GetHashCode ()
See hash_code
-
GetName ()
-
GetQualifiedName ()
-
GetNumber ()
-
GetOneLetterCode ()
-
GetBounds ()
-
GetMass ()
-
GetCenterOfAtoms ()
-
GetCenterOfMass ()
-
GetGeometricCenter ()
-
GetPhiTorsion ()
-
GetPsiTorsion ()
-
GetOmegaTorsion ()
-
GetChemClass ()
-
GetChemType ()
-
GetSecStructure ()
-
IsLigand ()
-
SetIsLigand ()
-
IsProtein ()
-
IsPeptideLinking ()
-
IsNucleotideLinking ()
-
GetCentralAtom ()
-
SetCentralAtom ()
-
GetCentralNormal ()
-
IsValid ()
-
GetNext ()
-
GetPrev ()
See the respective methods in ResidueHandle .
-
GetChain ()
See chain
-
FindAtom (atom_name)
Find atom by name. Returns the atom, or an invalid atom if no atoms with
the given name is part of the view.
Parameters: | atom_name (str) – |
Return type: | AtomView |
-
IsAtomIncluded (atom_handle)
Returns true if the given atom is part of the view, false if not.
Parameters: | atom_handle (AtomHandle ) – |
Return type: | bool |
-
AddAtom (atom_handle[, flags])
Add atom to the view.
-
GetAtomList ()
See atoms
-
GetAtomCount ()
See atom_count
-
GetIndex ()
See index
-
Select (query, flags=0)
Perform selection on residue view. See EntityHandle.Select() .
-
class
AtomView
A view representation of an AtomHandle . Mostly, the same
functionality is provided as for the handle.
-
handle
The underlying AtomHandle of the atom view. Also available as
GetHandle() .
-
hash_code
A unique identifier for this atom view. Note, that this is not the same as
for the atom handle (see AtomHandle.hash_code ).
Also available as GetHashCode() .
-
GetHandle ()
See handle
-
GetHashCode ()
See hash_code
Functions
Boolean Operators
-
Intersection (view_a, view_b)
Calculates and returns the intersection of view_a and view_b. view_a
and view_b must be views of the same entity.
-
Difference (view_a, view_b)
Calculates and returns the difference between view_a and view_b. view_a
and view_b must be views of the same entity.The returned view will
contain atoms, residues, chains and bonds that are in view_a and not in
view_b.
-
Union (view_a, view_b)
Calculates and returns the union of view_a and view_b. view_a
and view_b must be views of the same entity.The returned view will contain
all atoms, residues, chains and bonds that are either in view_a, view_b
or part of both views.
Residue Numbering
-
class
ResNum (num, ins_code='\0')
Number for a residue. The residue number has a numeric part and an (optional)
insertion-code. You can work with this object as if it was an integer and
comparison will look first at the numeric part and then the insertion-code.
All access to existing objects is read-only. Openstructure supports a range
of (-8388608 to 8388607) for the numeric part. However, the PDB format only
supports a range of (-999, 9999). This becomes relevant when a structure is
saved in PDB format where an IOException is raised if the PDB range is not
respected.
Parameters: |
- num (
int ) – Numeric part of residue number.
- ins_code (
str ) – Alpha-numeric part of residue number (optional insertion
code). Only first character kept.
|
-
num
Numeric part of residue number.
-
ins_code
Alpha-numeric part of residue number (insertion code). Single character.
-
GetNum ()
-
-
GetInsCode ()
-
ChainType
A ChainType fills the ChainHandle.type attribute. Different types are
described in the ChainType enum. The type is set with an editor in
EditorBase.SetChainType() . Further convenience functions are described
here.
-
class
ChainType
The ChainType enum enumerates all types defined by the PDB for the MMCif file
format. Following values are supported:
CHAINTYPE_POLY , CHAINTYPE_NON_POLY , CHAINTYPE_WATER ,
CHAINTYPE_POLY_PEPTIDE_D , CHAINTYPE_POLY_PEPTIDE_L ,
CHAINTYPE_POLY_DN , CHAINTYPE_POLY_RN , CHAINTYPE_POLY_SAC_D ,
CHAINTYPE_POLY_SAC_L , CHAINTYPE_POLY_DN_RN ,
CHAINTYPE_UNKNOWN , CHAINTYPE_MACROLIDE ,
CHAINTYPE_CYCLIC_PSEUDO_PEPTIDE , CHAINTYPE_POLY_PEPTIDE_DN_RN ,
CHAINTYPE_BRANCHED , CHAINTYPE_OLIGOSACCHARIDE ,
CHAINTYPE_N_CHAINTYPES
Where CHAINTYPE_N_CHAINTYPES holds the number of different types available.
-
ChainTypeFromString (identifier)
Create a ChainType item for a given string.
Parameters: | identifier (str ) – String to request a type for. |
Raises: | runtime_error if identifier is unrecognised. |
Returns: | ChainType |
-
StringFromChainType (type)
Return the String identifier for a given type.
Parameters: | type (ChainType ) – To be translated |
Raises: | runtime_error if type is unrecognised. |
Returns: | str |
ViewAddFlag
-
class
ViewAddFlag
Defines flags controlling behaviour of routines adding handles to views:
INCLUDE_ATOMS - Include all atoms when adding a residue handle to a view
INCLUDE_RESIDUES - Include all residues when adding a chain to a view
INCLUDE_CHAINS - Include all chains when creating a new entity view
INCLUDE_ALL = INCLUDE_ATOMS | INCLUDE_RESIDUES |
INCLUDE_CHAINS - Convenience flags to include all substructures
CHECK_DUPLICATES - If set, it will be checked that no duplicates are
created when adding a new handle
Flags can be ORed to combine them.
SecStructure
-
class
SecStructure (type)
Defines a secondary structure type following the types defined by DSSP.
-
IsHelical ()
Returns: | True, if the set type is any type of helix (i.e. ALPHA_HELIX,
PI_HELIX or THREE_TEN_HELIX) |
-
IsExtended ()
Returns: | True, if the set type is any type of beta sheet (i.e. EXTENDED,
BETA_BRIDGE) |
-
IsCoil ()
Returns: | True, if the set type is any type of coil (i.e. TURN, BEND or COIL) |
-
__str__ ()
-
-
class
SecStructure. Type
Enumerates all popssible secondary structure types distinguished by DSSP.
Their values with the corresponding character code are listed here:
ALPHA_HELIX = ‘H’
PI_HELIX = ‘I’
THREE_TEN_HELIX = ‘G’
EXTENDED = ‘E’
|
BETA_BRIDGE = ‘B’
TURN = ‘T’
BEND = ‘S’
COIL = ‘C’
|
ChemClass
-
class
ChemClass (chem_class)
The chemical class is used to broadly categorize residues based on their
chemical properties. For example, peptides belong to some PEPTIDE_LINKING
class. Possible values as constant variable names and as characters:
PEPTIDE_LINKING = ‘P’
D_PEPTIDE_LINKING = ‘D’
L_PEPTIDE_LINKING = ‘L’
RNA_LINKING = ‘R’
DNA_LINKING = ‘S’
NON_POLYMER = ‘N’
|
L_SACCHARIDE = ‘X’
D_SACCHARIDE = ‘Y’
SACCHARIDE = ‘Z’
WATER = ‘W’
UNKNOWN = ‘U’
|
The constants are defined directly within the mol module.
Python can implicitly convert characters to objects of this type. Note however
that only the first character of a str object is considered!
Parameters: | chem_class (str ) – Chemical class to set. |
-
IsPeptideLinking ()
Returns: | True, if set class is PEPTIDE_LINKING, D_PEPTIDE_LINKING or
L_PEPTIDE_LINKING |
-
IsNucleotideLinking ()
Returns: | True, if set class is RNA_LINKING or DNA_LINKING |
ChemType
-
class
ChemType
The chemical type of a residue is a classification of all compounds obtained
from the PDB component dictionary. For example, ions belong to the class IONS,
amino acids to AMINOACIDS. Possible values as constant variable names and as
characters:
IONS = ‘I’
NONCANONICALMOLS = ‘M’
SACCHARIDES = ‘S’
NUCLEOTIDES = ‘N’
AMINOACIDS = ‘A’
|
COENZYMES = ‘E’
WATERCOORDIONS = ‘C’
DRUGS = ‘D’
WATERS = ‘W’
UNKNOWN = ‘U’
|
Python can implicitly convert characters to objects of this type. Note however
that only the first character of a str object is considered!
Parameters: | chem_type (str ) – Chemical type to set. |
-
IsIon ()
Returns: | True, if set type is IONS or WATERCOORDIONS |
-
IsNucleotide ()
Returns: | True, if set type is NUCLEOTIDES |
-
IsSaccharide ()
Returns: | True, if set type is SACCHARIDES |
-
IsAminoAcid ()
Returns: | True, if set type is AMINOACIDS |
-
IsCoenzyme ()
Returns: | True, if set type is COENZYMES |
-
IsDrug ()
Returns: | True, if set type is DRUGS |
-
IsNonCanonical ()
Returns: | True, if set type is NONCANONICALMOLS |
-
IsWater ()
Returns: | True, if set type is WATERS |
-
IsKnown ()
Returns: | True, if set type is not UNKNOWN |
|
Contents
Search
Enter search terms or a module, class or function name.
Previous topic
mol – Molecular structures and surfaces
Next topic
Editors
You are here
|