SphericalHarmonics#

The spherical harmonic classes store the respective definitions and basis functions, as well as their derivative and inverse. Please refer to the Spherical Harmonic Definitions page for more specific information on the mathematical definitions of spherical harmonics.

Classes:

SphericalHarmonicDefinition([n_max, ...])

Class storing the (discrete) definition of spherical harmonics.

SphericalHarmonics(n_max, coordinates[, ...])

Compute spherical harmonic basis matrices, their inverses, and gradients.

class spharpy.SphericalHarmonicDefinition(n_max=0, basis_type='real', normalization='N3D', channel_convention='ACN', condon_shortley='auto')[source]#

Class storing the (discrete) definition of spherical harmonics.

This class can serve as a container to create related objects, e.g., spherical harmonic basis matrices for given sampling points, transforms, or other spherical harmonic related data and computations.

Parameters:
  • n_max (int, optional) – Maximum spherical harmonic order. Must be an integer greater or equal to 0. The default is 0.

  • basis_type (str) – Type of spherical harmonic basis, either 'complex' or 'real'. The default is 'real'.

  • normalization (str, optional) – Normalization convention, either 'N3D', 'NM', 'SN3D', 'SNM', or 'maxN'. 'maxN' is only supported up to 3rd order. The default is 'N3D'.

  • channel_convention (str, optional) – Channel ordering convention, either 'ACN' or 'FuMa'. 'FuMa' is only supported up to 3rd order. The default is 'ACN'.

  • condon_shortley (bool, str, optional) – Whether to include the Condon-Shortley phase term. If True, Condon-Shortley is included, if False it is not included. The default is 'auto', which corresponds to True for complex basis_type and False for real basis_type.

Attributes:

basis_type

Get or set the spherical harmonic basis type.

channel_convention

Get or set the channel ordering convention.

condon_shortley

Get or set the Condon-Shortley phase term.

n_max

Get or set the spherical harmonic order.

normalization

Get or set the normalization convention.

property basis_type#

Get or set the spherical harmonic basis type.

property channel_convention#

Get or set the channel ordering convention.

property condon_shortley#

Get or set the Condon-Shortley phase term.

property n_max#

Get or set the spherical harmonic order.

property normalization#

Get or set the normalization convention.

class spharpy.SphericalHarmonics(n_max, coordinates, basis_type='real', normalization='N3D', channel_convention='ACN', inverse_method='auto', condon_shortley='auto')[source]#

Compute spherical harmonic basis matrices, their inverses, and gradients.

The the basis functions that are used to build the matrices can be configured with the parameters described below.

Parameters:
  • n_max (int) – Maximum spherical harmonic order. Must be an integer greater or equal to 0.

  • coordinates (Coordinates, SamplingSphere) – The sampling points for which the matrices are computed.

  • basis_type (str, optional) – Type of spherical harmonic basis, either 'complex' or 'real'. The default is 'real'.

  • normalization (str, optional) – Normalization convention, either 'N3D', 'NM', 'SN3D', 'SNM', or 'maxN'. 'maxN' is only supported up to 3rd order. The default is 'N3D'.

  • channel_convention (str, optional) – Channel ordering convention, either 'ACN' or 'FuMa'. 'FuMa' is only supported up to 3rd order. The default is 'ACN'.

  • inverse_method (str, None, optional) –

    Method for computing the inverse transform:

    • 'quadrature': compute the inverse via numerical quadrature. Note that this requires coordinates to be a SamplingSphere object with weights summing to \(4\pi\).

    • 'pseudo_inverse': compute the inverse via the pseudo-inverse. Note that this requires coordinates to be a SamplingSphere object.

    • None: denotes that the inverse is not defined and cannot be computed.

    • 'auto': If coordinates are a SamplingSphere, use 'quadrature' when applicable and 'pseudo_inverse' otherwise. If coordinates are a Coordinates object None is used.

    The default is 'auto'.

  • condon_shortley (bool or str, optional) – Whether to include the Condon-Shortley phase term. If True, Condon-Shortley is included, if False it is not included. The default is 'auto', which corresponds to True for complex basis_type and False for real basis_type.

Attributes:

basis

Get the spherical harmonic basis matrix.

basis_gradient_phi

Get the gradient of the basis matrix with respect to the azimuth/polar angle.

basis_gradient_theta

Get the gradient of the basis matrix with respect to the colatitude/elevation angle.

basis_inv

Get the inverse basis matrix.

basis_type

Get or set the spherical harmonic basis type.

channel_convention

Get or set the channel ordering convention.

condon_shortley

Get or set the Condon-Shortley phase term.

coordinates

Get or set the coordinates.

inverse_method

Get or set the inverse transform method.

n_max

Get or set the spherical harmonic order.

normalization

Get or set the normalization convention.

Methods:

from_definition(definition, coordinates[, ...])

Create SphericalHarmonics instance from SphericalHarmonicDefinition.

property basis#

Get the spherical harmonic basis matrix.

property basis_gradient_phi#

Get the gradient of the basis matrix with respect to the azimuth/polar angle.

See the coordinates documentation for the definition of the angles.

property basis_gradient_theta#

Get the gradient of the basis matrix with respect to the colatitude/elevation angle.

See the coordinates documentation for the definition of the angles.

property basis_inv#

Get the inverse basis matrix.

property basis_type#

Get or set the spherical harmonic basis type.

property channel_convention#

Get or set the channel ordering convention.

property condon_shortley#

Get or set the Condon-Shortley phase term.

property coordinates#

Get or set the coordinates.

classmethod from_definition(definition, coordinates, inverse_method='auto')[source]#

Create SphericalHarmonics instance from SphericalHarmonicDefinition.

Parameters:
  • definition (SphericalHarmonicDefinition) – The spherical harmonic definition.

  • coordinates (Coordinates, SamplingSphere) – The sampling points for which the matrices are computed.

  • inverse_method (str, None, optional) –

    Method for computing the inverse transform:

    • 'quadrature': compute the inverse via numerical quadrature. Note that this requires coordinates to be a SamplingSphere object with weights summing to \(4\pi\).

    • 'pseudo_inverse': compute the inverse via the pseudo-inverse. Note that this requires coordinates to be a SamplingSphere object.

    • None: denotes that the inverse is not defined and cannot be computed. - 'auto': If coordinates are a SamplingSphere, use 'quadrature' when applicable and 'pseudo_inverse' otherwise. If coordinates are a Coordinates object None is used.

    The default is 'auto'.

Returns:

A new SphericalHarmonics instance initialized with the parameters from the definition object.

Return type:

SphericalHarmonics

property inverse_method#

Get or set the inverse transform method.

property n_max#

Get or set the spherical harmonic order.

property normalization#

Get or set the normalization convention.