The Molecular Entity¶
This document describes the EntityHandle and the related classes.
The Handle Classes¶
- ost.mol.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.
Returns: The newly created EntityHandle
- class ost.mol.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.
- 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().
returns: A list of residue handles
- 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()
Type : A list of atom handles
- bounds¶
Axis-aligned bounding box of the entity. Read-only
Type : ost.geom.AlignedCuboid
- center_of_mass¶
Center of mass. Also available as GetCenterOfMass()
Type : Vec3
- center_of_atoms¶
Center of atoms (not mass-weighted). Also available as GetCenterOfAtoms().
Type : Vec3
- 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.
- 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 (mol.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.
- 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 (mol.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
- EditXCS([edit_mode=mol.EditMode.BUFFERED_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.BUFFERED_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:
Parameters: 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()
Returns: EntityView
- 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)
Returns: EntityView
- 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.
- GetCenterOfAtoms()¶
Get center of atoms, that is the average atom position of the entity. Use GetCenterOfMass() to calculate the mass-weighted center of the entity.
Returns: Vec3
- GetCenterOfMass()¶
Calculates the center of mass of the entity. Use GetCenterOfAtoms() to calculate the non-mass-weighted center of the entity.
Returns: Vec3
- GetGeometricCenter()¶
Calculates the mid-point of the axis aligned bounding box of the entity.
Returns: Vec3
- FindWithin(pos, radius)¶
Find all atoms in sphere of given radius centered 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: A list of atom handles
- class ost.mol.ChainHandle¶
A chain of one or more residues. Chains are always part of an entity.
- 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()
Type : str
- 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().
Returns: A list of residue handles
- 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
- 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()
Type : A list of atom handles
- bounds¶
Axis-aligned bounding box of the chain. Read-only
Type : ost.geom.AlignedCuboid
- mass¶
The total mass of this chain in Dalton. Also available as GetMass()
Type : float
- center_of_atoms¶
Center of atoms (not mass weighted). Also available as GetCenterOfAtoms().
Type : Vec3
- FindResidue(res_num)¶
Get residue by residue number. See also residues
Parameters: res_num (mol.ResNum) – residue number Returns: A valid ResidueHandle if the chain contains a residue with matching residue number, an invalid ResidueHandle otherwise.
- GetResidueList()¶
Get list of all residues of this chain. For peptide chains, the residues are usually ordered from N- to C-terminus.To access a single residue, use FindResidue().
Returns: A list of residue handles
- FindAtom(res_num, atom_name)¶
Get atom by residue number and atom name. See also atoms
Parameters: - res_num (mol.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
- class ost.mol.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.
This property is read-only. To change the name of the residue, use EditorBase.SetResidueName().
- 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().
Type : ResNum
- 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])
Type : str
- bounds¶
Axis-aligned bounding box of the residue. Read-only.
Type : ost.geom.AlignedCuboid
- mass¶
The total mass of this residue in Dalton. Also available as GetMass().
Type : float
- center_of_mass¶
Center of mass. Also available as GetCenterOfMass()
Type : Vec3
- center_of_atoms¶
Center of atoms (not mass weighted). Also available as GetCenterOfAtoms().
Type : Vec3
- chain¶
The chain this residue belongs to. Read-only. Also available as GetChain()
Type : ChainHandle
- 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()
Type : TorsionHandle
- 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()
Type : TorsionHandle
- chem_class¶
The chemical class of a residue is used to broadly categorize residues based on their chemical properties. For example, peptides belong to the LPeptideLinking or DPeptideLinking classes.
- 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.
- 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
- GetAtomList()¶
Get list of all atoms of this residue. To access a single atom, use FindAtom().
- GetCenterOfAtoms()¶
See center_of_atoms
- GetCenterOfMass()¶
See center_of_mass
- GetPhiTorsion()¶
See phi_torsion
- GetPsiTorsion()¶
See psi_torsion
- class ost.mol.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.
Type : str
- qualified_name¶
The qualified name consists of the atom name as well as a unique residue identifier and chain name. For CA of a glycine with residue number 2 of chain A, the qualified name is “A.GLY2.CA”.
Type : str
- element¶
The atom’s element. Note that this may return an empty string. Also available as GetElement(). Read-only.
Type : str
- radius¶
The van-der-Waals radius of the atom. Also available as GetRadius(). Read/write.
Type : float
- occupancy¶
The atom’s occupancy in the range 0 to 1. Read/write. Also available as meth:GetOccupancy, SetOccupancy(). :type: float
- b_factor¶
The atom’s temperature factor. Read/write. Also available as GetBFactor(), SetBFactor().
Type : float
- charge¶
The atom’s charge. Also available as GetCharge(), SetCharge().
Type : float
- residue¶
The residue this atom belongs to. Read-only.
Type : ResidueHandle
- 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().
Type : list of bond handles
- 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.
Parameters: other_atom (AtomHandle) – The other atom Return type: BondHandle
- GetBondCount()¶
Return type: int
- 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
- GetEntity()¶
The entity this atom belongs to
Return type: EntityHandle
- GetHashCode()¶
Returns a unique identifier for this atom. :rtype: int
- GetIndex()¶
Returns the index of the atom.
Return type: int
- GetQualifiedName()¶
See qualified_name
Return type: str
- GetResidue()¶
See residue
Return type: ResidueHandle
- IsHetAtom()¶
See is_hetatom
Return type: bool
- IsValid()¶
See valid
Return type: bool
The View Classes¶
- class ost.mol.EntityView¶
An entity view represents a structural subset of an EntityHandle. For an introduction ,see Introduction to the mol Module.
- 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.
- 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().
Type : A list of ResidueViews
- 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()
Type : A list of AtomViews
- bounds¶
Axis-aligned bounding box of the entity view. Read-only.
Type : ost.geom.AlignedCuboid
- handle¶
The underlying handle of the entity view. Also available as GetHandle().
Type : EntityHandle
- CreateEmptyView()¶
See EntityHandle.CreateEmptyView()
Return type: EntityView
- CreateFullView()¶
Returns a copy of this entity. Provided for duck-typing <http://en.wikipedia.org/wiki/Duck_typing> purposes.
Return type: EntityView
- 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 behavior can be changed by passing in an appropriate set of view_add_flags
Parameters: - chain_handle (ChainHandle) –
- view_add_flags (int) – An ORed together combination of view_add_flags.
Return type: class:ChainView
- 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 view_add_flags.
Parameters: - residue_handle (ResidueHandle) –
- view_add_flags (int) –
Return type: class:ResidueView
- 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.
Parameters: - atom_handle (AtomHandle) – The atom handle
- view_add_flags (int) – An ORed together combination of ViewAddFlags
Return type: class:AtomView
- 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: - atom1 (AtomHandle) –
- atom2 (AtomHandle) –
- atom3 (AtomHandle) –
Return type: float
- FindWithin(pos, radius)¶
Find all atoms that are within radius of the given position. See EntityHandle.FindWithin().
Parameters: - pos (Vec3) –
- radius (float) –
Return type: class:AtomViewList
- FindChain(chain_name)¶
Find chain by name.
Parameters: chain_name (str) – Returns: The chain if present in the view, an invalid ChainView otherwise Return type: class: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: class: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: class:AtomView
- Select(query[, flags])¶
Perform selection on entity view. See EntityHandle.Select().
Parameters: - query ((Query or str)) – The query
- flags (int) – An ORed together combination of ViewAddFlags
Return type: EntityView
- Copy()¶
Returns a deep copy of the entity view. The effect is identical to calling
the_copy=view.Select(')
Return type: EntityView
- GetMass()¶
Get total mass of view.
Return type: float
- GetCenterOfMass()¶
Get the center of mass. For a non-mass-weighted center, see GetCenterOfAtoms().
Return type: Vec3
- GetGeometricCenter()¶
Get the geometric center (the center of the axis-aligned bounding box).
Return type: Vec3
- chain_count¶
Number of chains. Read-only. See GetChainCount().
Type : int
- residue_count¶
Number of residues. Read-only. See GetResidueCount().
Type : int
- atom_count¶
Number of atoms. Read-only. See GetAtomCount().
Type : int
- GetCenterOfAtoms()¶
Calculates the center of atoms. For a mass-weighted center, use GetCenterOfMass().
Return type: Vec3
- GetBondCount()¶
Get number of bonds :rtype: int
- GetChainCount()¶
Get number chains. See chain_count
Return type: int
- GetResidueCount()¶
See residue_count
Return type: int
- GetBondList()¶
See bonds
Return type: BondHandleList
- GetResidueList()¶
Return type: class:ResidueViewList
- GetAtomCount()¶
Get number of atoms. See :attr`atom_count`. :rtype: int
- class ost.mol.ChainView¶
- 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()
Type : str
- 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 : A list of residue views
- 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
- 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()
Type : A list of atom handles
- bounds¶
Axis-aligned bounding box of the chain. Read-only
Type : ost.geom.AlignedCuboid
- center_of_mass¶
Center of mass. Also available as GetCenterOfMass()
Type : Vec3
- center_of_atoms¶
Center of atoms (not mass weighted). Also available as GetCenterOfAtoms().
Type : Vec3
- handle¶
The chain handle this view points to. Also available as GetHandle().
Type : ChainHandle
- in_sequence
Whether the residue numbers are in ascending order. Note that the value of in_sequence is independent from the value of ChainHandle.in_sequence.
Type : bool
- 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.
Parameters: - atom_handle (AtomHandle) – The atom to be added
- view_add_flags (int) – An ORed together combination of ViewAddFlags
Return type: AtomView
- 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 view_add_flags.
Parameters: - residue_handle (ResidueHandle) –
- view_add_flags (int) –
Return type:
- 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()¶
See entity
- GetGeometricCenter()¶
See geometric_center
- 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
- 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])¶
Perform query on chain. This will return an entity view containing atoms and bonds of the residue view that match the residue. In case no atom matches the predicate, an empty view is returned.
Parameters: - query (Query or str) – The query
- flags (int) – An ORed together combination of query flags
Return type:
- class ost.mol.ResidueView¶
- handle¶
The residue handle this view points to. Also available as GetHandle().
Type : ResidueHandle
- 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.
This property is read-only. To change the name of the residue, use EditorBase.SetResidueName().
- 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().
Type : ResNum
- 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])
Type : str
- bounds¶
Axis-aligned bounding box of the residue view. Read-only
Type : ost.geom.AlignedCuboid
- center_of_mass¶
Center of mass. Also available as GetCenterOfMass()
Type : Vec3
- center_of_atoms¶
Center of atoms (not mass weighted). Also available as GetCenterOfAtoms().
Type : Vec3
- chain¶
The chain this residue belongs to. Read-only. Also available as GetChain()
Type : ChainView
- handle
The residue handle this view points to
Type : ResidueHandle
- 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
- 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
- GetIndex()¶
See index
- GetCenterOfMass()¶
See center_of_mass
- IsAtomIncluded(atom_handle)¶
Returns true if the given atom is part of the view, false if not.
Parameters: atom_handle (AtomHandle) – Return type: bool
- GetGeometricCenter()¶
See geometric_center
- AddAtom(atom_handle[, flags])¶
Add atom to the view. :param atom_handle: :type atom_handle: AtomHandle :param flags: An ORed together combination of ViewAddFlags. :type flags: int :rtype: AtomView
- GetCenterOfAtoms()¶
See center_of_atoms
- GetAtomList()¶
See atoms
- Select(query[, flags])¶
Perform selection on residue view. This method will return an entity view containing all atoms and bonds of the residue matching the predicate. In case no atom matches the predicate, this will return an empty view.
Parameters: - query (Query or str) – The query
- flags (int) – An ORed together combination of query flags
Return type:
Functions¶
Boolean Operators¶
- ost.mol.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.
Parameters: - view_a (EntityView) – first view
- view_b (EntityView) – second view
- ost.mol.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.
Parameters: - view_a (EntityView) – first view
- view_b (EntityView) – second view
- ost.mol.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.
Parameters: - view_a (EntityView) – first view
- view_b (EntityView) – second view