spharpy.transforms#

Rotation for data in the spherical harmonic domain. Please refer to the Spherical Harmonic Definitions page for background information on spherical harmonics.

Classes:

SphericalHarmonicRotation(quat, *args, **kwargs)

Class for rotations of coordinates and spherical harmonic expansions.

Functions:

rotation_z_axis(n_max, angle)

Rotation matrix for complex spherical harmonics around the z-axis by a given angle.

rotation_z_axis_real(n_max, angle)

Rotation matrix for real-valued spherical harmonics around the z-axis by a given angle.

wigner_d_function(n, m_dash, m, beta)

Wigner-d function for rotations around the y-axis.

wigner_d_rotation(n_max, alpha, beta, gamma)

Wigner-D rotation matrix for Euler rotations by angles (alpha, beta, gamma) around the (z,y,z)-axes.

wigner_d_rotation_real(n_max, alpha, beta, gamma)

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.SphericalHarmonicRotation(quat, *args, **kwargs)[source]#

Bases: Rotation

Class for rotations of coordinates and spherical harmonic expansions.

The class extends the scipy.spatial.transform.Rotation class and provides methods to express the rotation as spherical harmonic rotation matrices as well as to apply the rotation to spherical harmonic data objects. The spherical harmonic rotation matrices are computed based on the Wigner-D functions as described in [1] and [2] for complex and real-valued spherical harmonics, respectively.

To create a rotation object, use one of the available from_* methods, e.g., scipy.spatial.transform.Rotation.from_euler or scipy.spatial.transform.Rotation.from_quat for creation using Euler angles or quaternions, respectively. Also see the examples below.

Examples

The following examples demonstrate how to create a rotation and apply it to spherical harmonic coefficients and corresponding spherical harmonic data containers using the implemented methods and operators.

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import spharpy

Define the convention and order of the spherical harmonics and a set of coefficients. Note that the generated coefficients correspond to a dipole oriented along the y-axis.

>>> definition = spharpy.SphericalHarmonicDefinition(n_max=1)
>>> coefficients = np.array([0, 1, 0, 0])

Define a rotation by 45 degrees around the z-axis based on the Euler angles.

>>> angle = np.pi / 4
>>> R = spharpy.transforms.SphericalHarmonicRotation.from_euler(
>>>     'z', angle)

The spherical harmonic rotation matrix can be obtained and applied to the spherical harmonic coefficients using the following:

>>> D = R.as_spherical_harmonic_matrix(definition)
>>> rotated_coefficients = D @ coefficients
>>> print(f"Rotated coefficients: {np.round(rotated_coefficients, 2)}")
... # Rotated coefficients: [ 0.   0.71 0.   -0.71]

To visualize the effect of the rotation, we can expand the series expansion on the unit sphere and plot the original and rotated data:

>>> sampling = spharpy.samplings.equal_area(0, n_points=250)
>>> Y = spharpy.SphericalHarmonics.from_definition(
>>>     definition, coordinates=sampling)
...
>>> _, axs = plt.subplots(
>>>     1, 2, subplot_kw={'projection': '3d'}, figsize=(5, 2.5),
>>>     layout='constrained')
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, Y.basis @ coefficients, ax=axs[0], colorbar=False)
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, Y.basis @ rotated_coefficients, ax=axs[1],
>>>     colorbar=False)

(Source code, png, hires.png, pdf)

../_images/spharpy-transforms-1.png

The rotation can also be applied to spherical harmonic data objects using the overloaded multiplication operator. This includes the class object itself as well as SphericalHarmonicFrequencyData, SphericalHarmonicTimeData, and SphericalHarmonicSignal.

The following example demonstrates the application to an arbitrary SphericalHarmonicFrequencyData object containing the same series expansion in frequency_data.freq as above:

>>> frequency_data = spharpy.SphericalHarmonicFrequencyData(
>>>     np.atleast_2d(coefficients).T, frequencies=1e3,
>>>     basis_type=definition.basis_type,
>>>     normalization=definition.normalization,
>>>     channel_convention=definition.channel_convention,
>>>     condon_shortley=definition.condon_shortley)

The rotation can now be applied using the apply method

>>> rotated_frequency_data = R.apply(frequency_data)

or the * operator

>>> rotated_frequency_data = R * frequency_data

The effect of the rotation can again be visualized by evaluating the series expansion on the unit sphere:

>>> _, axs = plt.subplots(
>>>     1, 2, subplot_kw={'projection': '3d'}, figsize=(5, 2.5),
>>>     layout='constrained')
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, np.squeeze(Y.basis @ frequency_data.freq),
>>>     ax=axs[0], colorbar=False)
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, np.squeeze(Y.basis @ rotated_frequency_data.freq),
>>>     ax=axs[1], colorbar=False)

(Source code, png, hires.png, pdf)

../_images/spharpy-transforms-2.png

Multiple rotations can be combined by multiplying the respective SphericalHarmonicRotation objects.

>>> rotated_frequency_data = R * R * frequency_data

This corresponds to a rotation of 90 degrees around the z-axis.

>>> _, axs = plt.subplots(
>>>     1, 2, subplot_kw={'projection': '3d'}, figsize=(5, 2.5),
>>>     layout='constrained')
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, np.squeeze(Y.basis @ frequency_data.freq),
>>>     ax=axs[0], colorbar=False)
>>> spharpy.plot.balloon_wireframe(
>>>     sampling, np.squeeze(Y.basis @ rotated_frequency_data.freq),
>>>     ax=axs[1], colorbar=False)

(Source code, png, hires.png, pdf)

../_images/spharpy-transforms-3.png

Note

Currently, only spherical harmonic rotations matrices following the ACN indexing and N3D normalization are supported.

References

Methods:

apply(target)

Apply the rotation to spherical harmonic data object.

as_spherical_harmonic_matrix(...)

Express the rotation as spherical harmonic rotation matrices.

apply(target: SphericalHarmonicTimeData | SphericalHarmonicFrequencyData | SphericalHarmonicSignal) SphericalHarmonicTimeData | SphericalHarmonicFrequencyData | SphericalHarmonicSignal[source]#

Apply the rotation to spherical harmonic data object.

Note that the spherical harmonic definition of the target object is used to generate a fitting spherical harmonic rotation matrix.

Parameters:

target (SphericalHarmonicTimeData, SphericalHarmonicFrequencyData, SphericalHarmonicSignal) – Spherical harmonic series expansion to be rotated.

Returns:

The rotated signal.

Return type:

SphericalHarmonicTimeData, SphericalHarmonicFrequencyData, SphericalHarmonicSignal

as_spherical_harmonic_matrix(spherical_harmonic_definition: SphericalHarmonicDefinition) ndarray[source]#

Express the rotation as spherical harmonic rotation matrices.

Supports complex and real-valued spherical harmonics.

Parameters:

spherical_harmonic_definition (SphericalHarmonicDefinition) – Definition of the spherical harmonics.

Returns:

Stack of block-diagonal rotation matrices with individual shapes \([(n_{max}+1)^2, (n_{max}+1)^2]\).

Return type:

np.ndarray, complex or float

Note

At the moment, only spherical harmonic rotations matrices following the ACN indexing and N3D normalization are supported.

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 [3].

\[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

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 [4].

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

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 [5].

Parameters:
  • n (int) – order

  • m_dash (int) – degree

  • m (int) – degree

  • beta (float) – Rotation angle

Returns:

Wigner-d symbol

Return type:

float

References

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 [6]. 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

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 [7] 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