Bio.PDB.vectors module

Vector class, including rotation-related functions.

Bio.PDB.vectors.m2rotaxis(m)

Return angles, axis pair that corresponds to rotation matrix m.

The case where m is the identity matrix corresponds to a singularity where any rotation axis is valid. In that case, Vector([1, 0, 0]), is returned.

Bio.PDB.vectors.vector_to_axis(line, point)

Vector to axis method.

Return the vector between a point and the closest point on a line (ie. the perpendicular projection of the point on the line).

Parameters
  • line (L{Vector}) – vector defining a line

  • point (L{Vector}) – vector defining the point

Bio.PDB.vectors.rotaxis2m(theta, vector)

Calculate left multiplying rotation matrix.

Calculate a left multiplying rotation matrix that rotates theta rad around vector.

Parameters
  • theta (float) – the rotation angle

  • vector (L{Vector}) – the rotation axis

Returns

The rotation matrix, a 3x3 Numeric array.

Examples

>>> from numpy import pi
>>> from Bio.PDB.vectors import rotaxis2m
>>> from Bio.PDB.vectors import Vector
>>> m = rotaxis2m(pi, Vector(1, 0, 0))
>>> Vector(1, 2, 3).left_multiply(m)
<Vector 1.00, -2.00, -3.00>
Bio.PDB.vectors.rotaxis(theta, vector)

Calculate left multiplying rotation matrix.

Calculate a left multiplying rotation matrix that rotates theta rad around vector.

Parameters
  • theta (float) – the rotation angle

  • vector (L{Vector}) – the rotation axis

Returns

The rotation matrix, a 3x3 Numeric array.

Examples

>>> from numpy import pi
>>> from Bio.PDB.vectors import rotaxis2m
>>> from Bio.PDB.vectors import Vector
>>> m = rotaxis2m(pi, Vector(1, 0, 0))
>>> Vector(1, 2, 3).left_multiply(m)
<Vector 1.00, -2.00, -3.00>
Bio.PDB.vectors.refmat(p, q)

Return a (left multiplying) matrix that mirrors p onto q.

Returns

The mirror operation, a 3x3 Numeric array.

Examples

>>> from Bio.PDB.vectors import refmat
>>> p, q = Vector(1, 2, 3), Vector(2, 3, 5)
>>> mirror = refmat(p, q)
>>> qq = p.left_multiply(mirror)
>>> print(q)
<Vector 2.00, 3.00, 5.00>
>>> print(qq)
<Vector 1.21, 1.82, 3.03>
Bio.PDB.vectors.rotmat(p, q)

Return a (left multiplying) matrix that rotates p onto q.

Parameters
  • p (L{Vector}) – moving vector

  • q (L{Vector}) – fixed vector

Returns

rotation matrix that rotates p onto q

Return type

3x3 Numeric array

Examples

>>> from Bio.PDB.vectors import rotmat
>>> p, q = Vector(1, 2, 3), Vector(2, 3, 5)
>>> r = rotmat(p, q)
>>> print(q)
<Vector 2.00, 3.00, 5.00>
>>> print(p)
<Vector 1.00, 2.00, 3.00>
>>> p.left_multiply(r)
<Vector 1.21, 1.82, 3.03>
Bio.PDB.vectors.calc_angle(v1, v2, v3)

Calculate angle method.

Calculate the angle between 3 vectors representing 3 connected points.

Parameters

v2, v3 (v1,) – the tree points that define the angle

Returns

angle

Return type

float

Bio.PDB.vectors.calc_dihedral(v1, v2, v3, v4)

Calculate dihedral angle method.

Calculate the dihedral angle between 4 vectors representing 4 connected points. The angle is in ]-pi, pi].

Parameters

v2, v3, v4 (v1,) – the four points that define the dihedral angle

class Bio.PDB.vectors.Vector(x, y=None, z=None)

Bases: object

3D vector.

__init__(self, x, y=None, z=None)

Initialize the class.

__repr__(self)

Return vector 3D coordinates.

__neg__(self)

Return Vector(-x, -y, -z).

__add__(self, other)

Return Vector+other Vector or scalar.

__sub__(self, other)

Return Vector-other Vector or scalar.

__mul__(self, other)

Return Vector.Vector (dot product).

__truediv__(self, x)

Return Vector(coords/a).

__pow__(self, other)

Return VectorxVector (cross product) or Vectorxscalar.

__getitem__(self, i)

Return value of array index i.

__setitem__(self, i, value)

Assign values to array index i.

__contains__(self, i)

Validate if i is in array.

norm(self)

Return vector norm.

normsq(self)

Return square of vector norm.

normalize(self)

Normalize the Vector object.

Changes the state of self and doesn’t return a value. If you need to chain function calls or create a new object use the normalized method.

normalized(self)

Return a normalized copy of the Vector.

To avoid allocating new objects use the normalize method.

angle(self, other)

Return angle between two vectors.

get_array(self)

Return (a copy of) the array of coordinates.

left_multiply(self, matrix)

Return Vector=Matrix x Vector.

right_multiply(self, matrix)

Return Vector=Vector x Matrix.

copy(self)

Return a deep copy of the Vector.

Bio.PDB.vectors.homog_rot_mtx(angle_rads: float, axis: str)<built-in function array>

Generate a 4x4 single-axis numpy rotation matrix.

Parameters
  • angle_rads (float) – the desired rotation angle in radians

  • axis (char) – character specifying the rotation axis

Bio.PDB.vectors.set_Z_homog_rot_mtx(angle_rads: float, mtx: numpy.ndarray)

Update existing Z rotation matrix to new angle.

Bio.PDB.vectors.set_Y_homog_rot_mtx(angle_rads: float, mtx: numpy.ndarray)

Update existing Y rotation matrix to new angle.

Bio.PDB.vectors.set_X_homog_rot_mtx(angle_rads: float, mtx: numpy.ndarray)

Update existing X rotation matrix to new angle.

Bio.PDB.vectors.homog_trans_mtx(x: float, y: float, z: float)<built-in function array>

Generate a 4x4 numpy translation matrix.

Parameters

y, z (x,) – translation in each axis

Bio.PDB.vectors.set_homog_trans_mtx(x: float, y: float, z: float, mtx: numpy.ndarray)

Update existing translation matrix to new values.

Bio.PDB.vectors.homog_scale_mtx(scale: float)<built-in function array>

Generate a 4x4 numpy scaling matrix.

Parameters

scale (float) – scale multiplier

Bio.PDB.vectors.get_spherical_coordinates(xyz: <built-in function array>)Tuple[float, float, float]

Compute spherical coordinates (r, azimuth, polar_angle) for X,Y,Z point.

Parameters

xyz (array) – column vector (3 row x 1 column numpy array)

Returns

tuple of r, azimuth, polar_angle for input coordinate

Bio.PDB.vectors.coord_space(a0: numpy.ndarray, a1: numpy.ndarray, a2: numpy.ndarray, rev: bool = False)Tuple[numpy.ndarray, Union[numpy.ndarray, NoneType]]

Generate transformation matrix to coordinate space defined by 3 points.

New coordinate space will have:

acs[0] on XZ plane acs[1] origin acs[2] on +Z axis

Parameters
  • column array x3 acs (numpy) – X,Y,Z column input coordinates x3

  • rev (bool) – if True, also return reverse transformation matrix (to return from coord_space)

Returns

4x4 numpy array, x2 if rev=True

Bio.PDB.vectors.multi_rot_Z(angle_rads: numpy.ndarray)numpy.ndarray

Create [entries] numpy Z rotation matrices for [entries] angles.

Parameters
  • entries – int number of matrices generated.

  • angle_rads – numpy array of angles

Returns

entries x 4 x 4 homogeneous rotation matrices

Bio.PDB.vectors.multi_rot_Y(angle_rads: numpy.ndarray)numpy.ndarray

Create [entries] numpy Y rotation matrices for [entries] angles.

Parameters
  • entries – int number of matrices generated.

  • angle_rads – numpy array of angles

Returns

entries x 4 x 4 homogeneous rotation matrices