spharpy.transforms¶
Classes:
|
Class for rotations of coordinates and data. |
Functions:
|
Rotation matrix for complex spherical harmonics around the z-axis by a given angle. |
|
Rotation matrix for real-valued spherical harmonics around the z-axis by a given angle. |
|
Wigner-d function for rotations around the y-axis. |
|
Wigner-D rotation matrix for Euler rotations by angles (alpha, beta, gamma) around the (z,y,z)-axes. |
|
Wigner-D rotation matrix for Euler rotations for real-valued spherical harmonics by angles (alpha, beta, gamma) around the (z,y,z)-axes. |
- class spharpy.transforms.RotationSH(quat, n_max=0, *args, **kwargs)[source]¶
Bases:
scipy.spatial.transform.rotation.RotationClass for rotations of coordinates and data.
Methods:
__init__(quat[, n_max])Initialize
apply(coefficients[, type])Apply the rotation to L sets of spherical harmonic coefficients
as_spherical_harmonic([type])Export the rotation operations as a spherical harmonic rotation matrices.
from_euler(n_max, seq, angles[, degrees])Initialize from Euler angles.
from_matrix(n_max, matrix, **kwargs)Initialize from rotation matrix.
from_quat(n_max, quat, **kwargs)Initialize from quaternions.
from_rotvec(n_max, rotvec[, degrees])Initialize from rotation vectors.
Attributes:
The spherical harmonic order used for spherical harmonic rotation matrices.
- __init__(quat, n_max=0, *args, **kwargs)[source]¶
Initialize
- Parameters
quat (array_like, shape (N, 4) or (4,)) – Each row is a (possibly non-unit norm) quaternion in scalar-last (x, y, z, w) format. Each quaternion will be normalized to unit norm.
n_max (int) – The spherical harmonic order
- Returns
The rotation object with spherical harmonic order n_max.
- Return type
Note
Initializing using the constructor is not advised. Always use the respective
from_quatmethod.
- apply(coefficients, type='real')[source]¶
Apply the rotation to L sets of spherical harmonic coefficients
- Parameters
coefficients (array, complex, shape \(((n_max+1)^2, L)\)) – L sets of spherical harmonic coefficients with a respective order \(((n_max+1)^2\)
- Returns
The rotated data
- Return type
array, complex
- as_spherical_harmonic(type='real')[source]¶
Export the rotation operations as a spherical harmonic rotation matrices. Supports complex and real-valued spherical harmonics.
- Parameters
real (string, optional) – Spherical harmonic definition. Can either be ‘complex’ or ‘real’, by default ‘real’ is used.
- Returns
Stack of block-diagonal rotation matrices.
- Return type
array, complex or float
- classmethod from_euler(n_max, seq, angles, degrees=False, **kwargs)[source]¶
Initialize from Euler angles.
Rotations in 3-D can be represented by a sequence of 3 rotations around a sequence of axes. In theory, any three axes spanning the 3-D Euclidean space are enough. In practice, the axes of rotation are chosen to be the basis vectors. The three rotations can either be in a global frame of reference (extrinsic) or in a body centred frame of reference (intrinsic), which is attached to, and moves with, the object under rotation 3.
- Parameters
n_max (int) – Spherical harmonic order.
seq (str) – Specifies sequence of axes for rotations. Up to 3 characters belonging to the set {‘X’, ‘Y’, ‘Z’} for intrinsic rotations, or {‘x’, ‘y’, ‘z’} for extrinsic rotations. Extrinsic and intrinsic rotations cannot be mixed in one function call.
angles ((float or array_like, shape (N,) or (N, [1 or 2 or 3]))) –
Euler angles specified in radians (degrees is False) or degrees (degrees is True). For a single character seq, angles can be:
a single value
array_like with shape (N,), where each
angle[i]corresponds to a single rotationarray_like with shape (N, 1), where each
angle[i, 0]corresponds to a single rotation
For 2- and 3-character wide seq, angles can be:
array_like with shape (W,) where W is the width of
seq, which corresponds to a single rotation with W axesarray_like with shape (N, W) where each
angle[i]corresponds to a sequence of Euler angles describing a single rotation
degrees (bool, optional) – If True, then the given angles are assumed to be in degrees. Default is False.
- Returns
Object containing the rotation represented by the sequence of rotations around given axes with given angles.
- Return type
References
Examples
>>> from spharpy.transforms import Rotation as R >>> rot = R.from_euler('z', 90, degrees=True) >>> rot.as_spherical_harmonic()
- classmethod from_matrix(n_max, matrix, **kwargs)[source]¶
Initialize from rotation matrix. Rotations in 3 dimensions can be represented with 3 x 3 proper orthogonal matrices 1. If the input is not proper orthogonal, an approximation is created using the method described in 2.
- Parameters
n_max (int) – Spherical harmonic order
matrix ((array_like, shape (N, 3, 3) or (3, 3))) – A single matrix or a stack of matrices, where
matrix[i]is the i-th matrix.
- Returns
Object containing the rotations represented by the rotation matrices.
- Return type
References
- 1
- 2
F. Landis Markley, “Unit Quaternion from Rotation Matrix”, Journal of guidance, control, and dynamics vol. 31.2, pp. 440-442, 2008.
Examples
>>> from spharpy.transforms import Rotation as R >>> rot = R.from_matrix([ ... [0, -1, 0], ... [1, 0, 0], ... [0, 0, 1]]) >>> rot.as_spherical_harmonic()
- classmethod from_quat(n_max, quat, **kwargs)[source]¶
Initialize from quaternions. 3D rotations can be represented using unit-norm quaternions 4.
- Parameters
n_max (int) – Spherical harmonic order
quat ((array_like, shape (N, 4) or (4,))) – Each row is a (possibly non-unit norm) quaternion in scalar-last (x, y, z, w) format. Each quaternion will be normalized to unit norm.
- Returns
Object containing the rotations represented by input quaternions.
- Return type
References
Examples
>>> from spharpy.transforms import Rotation as R >>> rot = R.from_quat([1, 0, 0, 0]) >>> rot.as_spherical_harmonic()
- classmethod from_rotvec(n_max, rotvec, degrees=False, *args, **kwargs)[source]¶
Initialize from rotation vectors. A rotation vector is a 3 dimensional vector which is co-directional to the axis of rotation and whose norm gives the angle of rotation 5.
- Parameters
n_max (int) – Spherical harmonic order
rotvec (array_like, float, shape (L, 3) or (3,)) – A single vector or a stack of vectors, where rot_vec[i] gives the ith rotation vector.
degrees (bool, optional) – Specify if rotation angles are defined in degrees instead of radians, by default False.
- Returns
Object containing the rotations represented by input rotation vectors.
- Return type
References
Examples
>>> from spharpy.transforms import Rotation as R >>> rot = R.from_rotvec(np.pi/2 * np.array([0, 0, 1])) >>> rot.as_spherical_harmonic()
- property n_max¶
The spherical harmonic order used for spherical harmonic rotation matrices.
- spharpy.transforms.rotation_z_axis(n_max, angle)[source]¶
Rotation matrix for complex spherical harmonics around the z-axis by a given angle. The rotation is performed such that positive angles result in a counter clockwise rotation of the data 6.
\[c_{nm}(\theta, \phi + \xi) = e^{-im\xi} c_{nm}(\theta, \phi)\]- Parameters
n_max (integer) – Spherical harmonic order
angle (number) – Rotation angle in radians [0, 2 pi]
- Returns
Diagonal rotation matrix evaluated for the specified angle
- Return type
array_like, complex, shape \(((n_{max}+1)^2, (n_{max}+1)^2)\)
References
- 6
N. A. Gumerov and R. Duraiswami, “Recursions for the computation of multipole translation and rotation coefficients for the 3-d helmholtz equation,” vol. 25, no. 4, pp. 1344–1381, 2003.
Examples
>>> import numpy as np >>> import spharpy >>> n_max = 1 >>> sh_vec = np.array([0, 1, 0, 0]) >>> rotMat = spharpy.transforms.rotation_z_axis(n_max, np.pi/2) >>> sh_vec_rotated = rotMat @ sh_vec
- spharpy.transforms.rotation_z_axis_real(n_max, angle)[source]¶
Rotation matrix for real-valued spherical harmonics around the z-axis by a given angle. The rotation is performed such that positive angles result in a counter clockwise rotation of the data 7.
- Parameters
n_max (integer) – Spherical harmonic order
angle (number) – Rotation angle in radians [0, 2 pi]
- Returns
Block-diagonal Rotation matrix evaluated for the specified angle.
- Return type
array_like, float, shape \(((n_max+1)^2, (n_max+1)^2)\)
References
- 7
M. Kronlachner, “Spatial Transformations for the Alteration of Ambisonic Recordings,” Master Thesis, University of Music and Performing Arts, Graz, 2014.
Examples
>>> import numpy as np >>> import spharpy >>> n_max = 1 >>> sh_vec = np.array([0, 1, 0, 0]) >>> rotMat = spharpy.transforms.rotation_z_axis_real(n_max, np.pi/2) >>> sh_vec_rotated = rotMat @ sh_vec
- spharpy.transforms.wigner_d_function(n, m_dash, m, beta)[source]¶
Wigner-d function for rotations around the y-axis. Convention as defined in 8.
- Parameters
n (int) – order
m_dash (int) – degree
m (int) – degree
beta (float) – Rotation angle
- Returns
Wigner-d symbol
- Return type
float
References
- 8
B. Rafaely, Fundamentals of Spherical Array Processing, 1st ed., vol. 8. Springer-Verlag GmbH Berlin Heidelberg, 2015.
- spharpy.transforms.wigner_d_rotation(n_max, alpha, beta, gamma)[source]¶
Wigner-D rotation matrix for Euler rotations by angles (alpha, beta, gamma) around the (z,y,z)-axes. The implementation follows 9. and rotation is performed such that positive angles result in a counter clockwise rotation of the data.
\[D_{m^\prime,m}^n(\alpha, \beta, \gamma) = e^{-im^\prime\alpha} d_{m^\prime,m}^n(\beta) e^{-im\gamma}\]- Parameters
n_max (int) – Spherical harmonic order
alpha (float) – First z-axis rotation angle
beta (float) – Y-axis rotation angle
gamma (float) – Second z-axis rotation angle
- Returns
Block diagonal rotation matrix
- Return type
array_like, complex,:math:((n_max+1)^2, (n_max+1)^2)
Examples
>>> import numpy as np >>> import spharpy >>> n_max = 1 >>> sh_vec = np.array([0, 0, 1, 0]) >>> rotMat = spharpy.transforms.wigner_d_rotation(n_max, 0, np.pi/4, 0) >>> sh_vec_rotated = rotMat @ sh_vec
References
- 9
B. Rafaely, Fundamentals of Spherical Array Processing, 1st ed., vol. 8. Springer-Verlag GmbH Berlin Heidelberg, 2015.
- spharpy.transforms.wigner_d_rotation_real(n_max, alpha, beta, gamma)[source]¶
Wigner-D rotation matrix for Euler rotations for real-valued spherical harmonics by angles (alpha, beta, gamma) around the (z,y,z)-axes. The implementation follows 10 and the rotation is performed such that positive angles result in a counter clockwise rotation of the data.
- Parameters
n_max (int) – Spherical harmonic order
alpha (float) – First z-axis rotation angle
beta (float) – Y-axis rotation angle
gamma (float) – Second z-axis rotation angle
- Returns
Block diagonal rotation matrix
- Return type
array_like, float, \(((n_max+1)^2, (n_max+1)^2)\)
Examples
>>> import numpy as np >>> import spharpy >>> n_max = 1 >>> sh_vec = np.array([0, 0, 1, 0]) >>> rotMat = spharpy.transforms.wigner_d_rotation_real( >>> n_max, 0, np.pi/4, 0) >>> sh_vec_rotated = rotMat @ sh_vec
References
- 10
M. A. Blanco, M. Flórez, and M. Bermejo, “Evaluation of the rotation matrices in the basis of real spherical harmonics,” Journal of Molecular Structure: THEOCHEM, vol. 419, no. 1–3, pp. 19–27, Dec. 1997, doi: 10.1016/S0166-1280(97)00185-1.