spyrmsd.qcp module

spyrmsd.qcp.K_mtx(M)[source]

Compute symmetric key matrix.

Parameters

M (numpy.ndarray) – Inner product between coordinate matrices

Returns

Symmetric key matrix

Return type

numpy.ndarray

Notes

The symmetric key matrix corresponds to the matrix \(\mathbf{K}\). 2

If \(S_{xy}\) is defined as

\[S_{xy} = \sum_i^N x_{B,i} y_{A,i}\]

then \(\mathbf{K}\) is the \(4\times 4\) symmetric matrix given by

\[\begin{split}\begin{pmatrix} S_{xx} + S_{yy} + S_{zz} & S_{yz} - S_{zy} & S_{zx} - S_{xz} & S_{xy} - S_{yx} \\ & S_{xx} - S_{yy} - S_{zz} & S_{xy} + S_{yx} & S_{zx} + S_{xz}\\ & & -S_{xx} + S_{yy} - S_{zz} & S_{yz} - S_{zy} \\ & & & -S_{xx} - S_{yy} + S_{zz} \\ \end{pmatrix}\end{split}\]
2

D. L. Theobald, Rapid calculation of RMSDs using a quaternion-based characteristic polynomial, Acta Crys. A 61, 478-480 (2005).

spyrmsd.qcp.M_mtx(A: ndarray, B: ndarray) ndarray[source]

Compute inner product between coordinate matrices.

Parameters
  • A (numpy.ndarray) – Coordinates A

  • B (numpy.ndarray) – Coordinates B

Returns

Inner product of the coordinate matrices A and B

Return type

numpy.ndarray

Notes

The inner product of the coordinate matrices A and B corresponds to the matrix \(\mathbf{M}\). 1

If \(S_{xy}\) is defined as

\[S_{xy} = \sum_i^N x_{B,i} y_{A,i}\]

then \(\mathbf{M}\) is the \(3\times 3\) matrix given by

\[\begin{split}\begin{pmatrix} S_{xx} & S_{xy} & S_{xz} \\ S_{yx} & S_{yy} & S_{yz} \\ S_{zx} & S_{zy} & S_{zz} \\ \end{pmatrix}\end{split}\]
1

D. L. Theobald, Rapid calculation of RMSDs using a quaternion-based characteristic polynomial, Acta Crys. A 61, 478-480 (2005).

spyrmsd.qcp.coefficients(M: ndarray, K: ndarray) Tuple[float, float, float][source]

Compute quaternion polynomial coefficients.

Parameters
  • M (numpy.ndarray) – Inner product between coordinate matrices

  • K (numpy.ndarray) – Symmetric key matrix

Returns

Quaternion polynomial coefficients

Return type

Tuple[float, float, float]

Notes

Returns only \(\mathbf{M}\)- and \(\mathbf{K}\)-dependent coefficients are returned. \(c_4=1\) and \(c_3=0\) are not returned.

The \(\mathbf{M}\)- and \(\mathbf{K}\)-dependent quaternion polynomial coefficients are given by

\[c_2 = -2 \text{ tr}\left(\mathbf{M}^T\mathbf{M}\right)\]
\[c_1 = -8 \text{ det}(\mathbf{M})\]
\[c_0 = \text{ det}(\mathbf{K})\]
spyrmsd.qcp.lambda_max(Ga: float, Gb: float, c2: float, c1: float, c0: float) float[source]

Find largest root of the quaternion polynomial.

Parameters
  • Ga (float) – Inner product of structure A

  • Gb – Inner product of structure B

  • c2 – Coefficient \(c_2\) of the quaternion polynomial

  • c1 – Coefficient \(c_1\) of the quaternion polynomial

  • c0 – Coefficient \(c_0\) of the quaternion polynomial

Returns

Largest root of the quaternion polynomial (\(\lambda_\text{max}\))

Return type

float

spyrmsd.qcp.qcp_rmsd(A: ndarray, B: ndarray, atol: float = 1e-09) float[source]

Compute RMSD using the quaternion polynomial method.

Parameters
  • A (numpy.ndarray) – Coordinates of structure A

  • B (numpy.ndarray) – Coordinates of structure B

  • atol (float) – Absolute tolerance parameter (see notes)

Returns

RMSD between structures A and B

Return type

float

Raises

AssertionError – If the shape of structures A and B is different

Notes

If the structures A and B can be superimposed exactly (i.e. they differ only by center-of-mass translations and rotations), we have

\[G_a + G_b = 2 \lambda_\text{max}\]

This means that \(s = G_a + G_bb - 2 * \lambda_\text{max}\) can become negative because of numerical errors and therefore \(\sqrt{s}\) fails. In order to avoid this problem, the final RMSD is set to \(0\) if \(|s| < atol\).