Bio.PDB.Atom module
Atom class, used in Structure objects.
- class Bio.PDB.Atom.Atom(name: str, coord: ndarray, bfactor: float | None, occupancy: float | None, altloc: str, fullname: str, serial_number, element: str | None = None, pqr_charge: float | None = None, radius: float | None = None)
Bases:
object
Define Atom class.
The Atom object stores atom name (both with and without spaces), coordinates, B factor, occupancy, alternative location specifier and (optionally) anisotropic B factor and standard deviations of B factor and positions.
In the case of PQR files, B factor and occupancy are replaced by atomic charge and radius.
- __init__(name: str, coord: ndarray, bfactor: float | None, occupancy: float | None, altloc: str, fullname: str, serial_number, element: str | None = None, pqr_charge: float | None = None, radius: float | None = None)
Initialize Atom object.
- Parameters:
name (string) – atom name (eg. “CA”). Note that spaces are normally stripped.
coord (NumPy array (Float0, length 3)) – atomic coordinates (x,y,z)
bfactor (number) – isotropic B factor
occupancy (number) – occupancy (0.0-1.0)
altloc (string) – alternative location specifier for disordered atoms
fullname (string) – full atom name, including spaces, e.g. “ CA “. Normally these spaces are stripped from the atom name.
element (uppercase string (or None if unknown)) – atom element, e.g. “C” for Carbon, “HG” for mercury,
pqr_charge (number) – atom charge
radius (number) – atom radius
- __eq__(other)
Test equality.
- __ne__(other)
Test inequality.
- __gt__(other)
Test greater than.
- __ge__(other)
Test greater or equal.
- __lt__(other)
Test less than.
- __le__(other)
Test less or equal.
- __hash__()
Return atom full identifier.
- __repr__()
Print Atom object as <Atom atom_name>.
- __sub__(other)
Calculate distance between two atoms.
- Parameters:
other (L{Atom}) – the other atom
Examples
This is an incomplete but illustrative example:
distance = atom1 - atom2
- strictly_equals(other: _AtomT, compare_coordinates: bool = False) bool
Compare this atom to the other atom using a strict definition of equality.
Indicates whether the atoms have the same name, B factor, occupancy, alternate location indicator (altloc), fullname, element, charge, and radius. If
compare_coordinates
is true, then the coordinates are also compared.- Parameters:
other (Atom) – The atom to compare this atom with
compare_coordinates (bool) – Whether to compare the coordinates of the atoms
- Returns:
Whether the atoms are strictly equal
- Return type:
bool
- set_serial_number(n)
Set serial number.
- set_bfactor(bfactor: float | None)
Set isotroptic B factor.
- set_coord(coord: ndarray)
Set coordinates.
- set_altloc(altloc: str)
Set alternative location specifier.
- set_occupancy(occupancy: float | None)
Set occupancy.
- set_sigatm(sigatm_array: ndarray | None)
Set standard deviation of atomic parameters.
The standard deviation of atomic parameters consists of 3 positional, 1 B factor and 1 occupancy standard deviation.
- Parameters:
sigatm_array (NumPy array (length 5)) – standard deviations of atomic parameters.
- set_siguij(siguij_array: ndarray | None)
Set standard deviations of anisotropic temperature factors.
- Parameters:
siguij_array (NumPy array (length 6)) – standard deviations of anisotropic temperature factors.
- set_anisou(anisou_array: ndarray | None)
Set anisotropic B factor.
- Parameters:
anisou_array (NumPy array (length 6)) – anisotropic B factor.
- set_charge(pqr_charge: float | None)
Set charge.
- set_radius(radius: float | None)
Set radius.
- flag_disorder()
Set the disordered flag to 1.
The disordered flag indicates whether the atom is disordered or not.
- is_disordered() int
Return the disordered flag (1 if disordered, 0 otherwise).
- set_parent(parent)
Set the parent residue.
- Arguments:
parent - Residue object
- detach_parent()
Remove reference to parent.
- get_sigatm() ndarray | None
Return standard deviation of atomic parameters.
- get_siguij() ndarray | None
Return standard deviations of anisotropic temperature factors.
- get_anisou() ndarray | None
Return anisotropic B factor.
- get_serial_number()
Return the serial number.
- get_name() str
Return atom name.
- get_id() str
Return the id of the atom (which is its atom name).
- get_full_id()
Return the full id of the atom.
The full id of an atom is a tuple used to uniquely identify the atom and consists of the following elements: (structure id, model id, chain id, residue id, atom name, altloc)
- get_coord() ndarray
Return atomic coordinates.
- get_bfactor() float | None
Return B factor.
- get_occupancy() float | None
Return occupancy.
- get_fullname() str
Return the atom name, including leading and trailing spaces.
- get_altloc() str
Return alternative location specifier.
- get_level() str
Return level.
- get_charge() float | None
Return charge.
- get_radius() float | None
Return radius.
- transform(rot: ndarray, tran: ndarray)
Apply rotation and translation to the atomic coordinates.
- Parameters:
rot (3x3 NumPy array) – A right multiplying rotation matrix
tran (size 3 NumPy array) – the translation vector
Examples
This is an incomplete but illustrative example:
from numpy import pi, array from Bio.PDB.vectors import Vector, rotmat rotation = rotmat(pi, Vector(1, 0, 0)) translation = array((0, 0, 1), 'f') atom.transform(rotation, translation)
- get_vector() Vector
Return coordinates as Vector.
- Returns:
coordinates as 3D vector
- Return type:
Bio.PDB.Vector class
- copy()
Create a copy of the Atom.
Parent information is lost.
- class Bio.PDB.Atom.DisorderedAtom(id: str)
Bases:
DisorderedEntityWrapper
Contains all Atom objects that represent the same disordered atom.
One of these atoms is “selected” and all method calls not caught by DisorderedAtom are forwarded to the selected Atom object. In that way, a DisorderedAtom behaves exactly like a normal Atom. By default, the selected Atom object represents the Atom object with the highest occupancy, but a different Atom object can be selected by using the disordered_select(altloc) method.
- __init__(id: str)
Create DisorderedAtom.
- Arguments:
id - string, atom name
- __iter__()
Iterate through disordered atoms.
- __repr__()
Return disordered atom identifier.
- center_of_mass()
Return the center of mass of the DisorderedAtom as a numpy array.
Assumes all child atoms have the same mass (same element).
- disordered_get_list() list[Atom]
Return list of atom instances.
Sorts children by altloc (empty, then alphabetical).
- disordered_add(atom)
Add a disordered atom.
- disordered_remove(altloc: str)
Remove a child atom altloc from the DisorderedAtom.
- Arguments:
altloc - name of the altloc to remove, as a string.
- transform(rot: ndarray, tran: ndarray)
Apply rotation and translation to all children.
See the documentation of Atom.transform for details.