SAF
Loading...
Searching...
No Matches
saf_sh

Files

file  saf_sh.c
 Public source for the Spherical Harmonic Transform and Spherical Array Processing module (SAF_SH_MODULE)
 
file  saf_sh.h
 Main header for the Spherical Harmonic Transform and Spherical Array Processing module (SAF_SH_MODULE)
 
file  saf_sh_internal.c
 Internal source for the Spherical Harmonic Transform and Spherical Array Processing module (SAF_SH_MODULE)
 
file  saf_sh_internal.h
 Internal header for the Spherical Harmonic Transform and Spherical Array Processing module (SAF_SH_MODULE)
 

Macros

#define ORDER2NSH(order)   ((order+1)*(order+1))
 Converts spherical harmonic order to number of spherical harmonic components i.e: (order+1)^2.
 
#define NSH2ORDER(nSH)   ( (int)(sqrt((double)nSH)-0.999) )
 Converts number of spherical harmonic components to spherical harmonic order i.e: sqrt(nSH)-1.
 

Enumerations

enum  ARRAY_CONSTRUCTION_TYPES { ARRAY_CONSTRUCTION_OPEN , ARRAY_CONSTRUCTION_OPEN_DIRECTIONAL , ARRAY_CONSTRUCTION_RIGID , ARRAY_CONSTRUCTION_RIGID_DIRECTIONAL }
 Microphone/Hydrophone array construction types. More...
 
enum  SECTOR_PATTERNS { SECTOR_PATTERN_PWD , SECTOR_PATTERN_MAXRE , SECTOR_PATTERN_CARDIOID }
 Sector pattern designs for directionally-constraining sound-fields [1]. More...
 
enum  ARRAY_SHT_OPTIONS { ARRAY_SHT_DEFAULT , ARRAY_SHT_REG_LS , ARRAY_SHT_REG_LSHD }
 Microphone array to spherical harmonic domain conversion options. More...
 

Functions

void unnorm_legendreP (int n, double *x, int lenX, double *y)
 Calculates unnormalised Legendre polynomials up to order N, for all values in vector x [1].
 
void unnorm_legendreP_recur (int n, float *x, int lenX, float *Pnm_minus1, float *Pnm_minus2, float *Pnm)
 Calculates unnormalised Legendre polynomials up to order N, for all values in vector x.
 
void getSHreal (int order, float *dirs_rad, int nDirs, float *Y)
 Computes real-valued spherical harmonics [1] for each given direction on the unit sphere.
 
void getSHreal_recur (int order, float *dirs_rad, int nDirs, float *Y)
 Computes real-valued spherical harmonics [1] for each given direction on the unit sphere.
 
void getSHreal_part (int order_start, int order_end, float *dirs_rad, int nDirs, float *Y)
 
void getSHcomplex (int order, float *dirs_rad, int nDirs, float_complex *Y)
 Computes complex-valued spherical harmonics [1] for each given direction on the unit sphere.
 
void complex2realSHMtx (int order, float_complex *T_c2r)
 Computes a complex to real spherical harmonic transform matrix.
 
void real2complexSHMtx (int order, float_complex *T_r2c)
 Computes a real to complex spherical harmonic transform matrix.
 
void complex2realCoeffs (int order, float_complex *C_N, int K, float *R_N)
 Converts SH coefficients from the complex to real basis.
 
void getSHrotMtxReal (float R[3][3], float *RotMtx, int L)
 Generates a real-valued spherical harmonic rotation matrix [1] based on a 3x3 rotation matrix (see quaternion2rotationMatrix(), euler2rotationMatrix())
 
void computeVelCoeffsMtx (int sectorOrder, float_complex *A_xyz)
 Computes the matrices which generate the coefficients of a beampattern of order (sectorOrder+1) that is essentially the product of a pattern of order=sectorOrder and a dipole.
 
float computeSectorCoeffsEP (int orderSec, float_complex *A_xyz, SECTOR_PATTERNS pattern, float *sec_dirs_deg, int nSecDirs, float *sectorCoeffs)
 Computes beamforming matrices (sector coefficients) which, when applied to input SH signals, yield energy-preserving (EP) sectors.
 
float computeSectorCoeffsAP (int orderSec, float_complex *A_xyz, SECTOR_PATTERNS pattern, float *sec_dirs_deg, int nSecDirs, float *sectorCoeffs)
 Computes beamforming matrices (sector coefficients) which, when applied to input SH signals, yield amplitude-preserving (EP) sectors.
 
void beamWeightsCardioid2Spherical (int N, float *b_n)
 Generates spherical coefficients for generating cardioid beampatterns.
 
void beamWeightsDolphChebyshev2Spherical (int N, int paramType, float arrayParam, float *b_n)
 Generates beamweights in the SHD for Dolph-Chebyshev beampatterns, with mainlobe and sidelobe control [1].
 
void beamWeightsHypercardioid2Spherical (int N, float *b_n)
 Generates beamweights in the SHD for hypercardioid beampatterns.
 
void beamWeightsMaxEV (int N, float *b_n)
 Generates beamweights in the SHD for maximum energy-vector beampatterns.
 
void beamWeightsVelocityPatternsReal (int order, float *b_n, float azi_rad, float elev_rad, float_complex *A_xyz, float *velCoeffs)
 Generates beamforming coefficients for velocity patterns (REAL)
 
void beamWeightsVelocityPatternsComplex (int order, float *b_n, float azi_rad, float elev_rad, float_complex *A_xyz, float_complex *velCoeffs)
 Generates beamforming coefficients for velocity patterns (COMPLEX)
 
void rotateAxisCoeffsReal (int order, float *c_n, float theta_0, float phi_0, float *c_nm)
 Generates spherical coefficients for a rotated axisymmetric pattern (REAL)
 
void rotateAxisCoeffsComplex (int order, float *c_n, float theta_0, float phi_0, float_complex *c_nm)
 Generates spherical coefficients for a rotated axisymmetric pattern (COMPLEX)
 
void checkCondNumberSHTReal (int order, float *dirs_rad, int nDirs, float *w, float *cond_N)
 Computes the condition numbers for a least-squares SHT.
 
int calculateGridWeights (float *dirs_rad, int nDirs, int order, float *w)
 Computes approximation of quadrature/integration weights for a given grid.
 
void sphPWD_create (void **const phPWD, int order, float *grid_dirs_deg, int nDirs)
 Creates an instance of a spherical harmonic domain implementation of the steer-response power (SRP) approach for computing power-maps, which can then be used for sound-field visualisation/DoA estimation purposes.
 
void sphPWD_destroy (void **const phPWD)
 Destroys an instance of the spherical harmonic domain PWD implementation.
 
void sphPWD_compute (void *const hPWD, float_complex *Cx, int nSrcs, float *P_map, int *peak_inds)
 Computes a power-map based on determining the energy of hyper-cardioid beamformers; optionally, also returning the grid indices corresponding to the N highest peaks (N=nSrcs)
 
void sphMUSIC_create (void **const phMUSIC, int order, float *grid_dirs_deg, int nDirs)
 Creates an instance of the spherical harmonic domain MUSIC implementation, which may be used for computing pseudo-spectrums for visualisation/DoA estimation purposes.
 
void sphMUSIC_destroy (void **const phMUSIC)
 Destroys an instance of the spherical harmonic domain MUSIC implementation.
 
void sphMUSIC_compute (void *const hMUSIC, float_complex *Vn, int nSrcs, float *P_music, int *peak_inds)
 Computes a pseudo-spectrum based on the MUSIC algorithm in the spherical harmonic domain; optionally returning the grid indices corresponding to the N highest peaks (N=nSrcs)
 
void sphESPRIT_create (void **const phESPRIT, int order)
 Creates an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.
 
void sphESPRIT_destroy (void **const phESPRIT)
 Destroys an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.
 
void sphESPRIT_estimateDirs (void *const hESPRIT, float_complex *Us, int K, float *src_dirs_rad)
 Estimates the direction-of-arrival (DoA) based on the ESPRIT-based estimator, in the spherical harmonic domain.
 
void generatePWDmap (int order, float_complex *Cx, float_complex *Y_grid, int nGrid_dirs, float *pmap)
 Generates a powermap based on the energy of a plane-wave decomposition (PWD) (i.e.
 
void generateMVDRmap (int order, float_complex *Cx, float_complex *Y_grid, int nGrid_dirs, float regPar, float *pmap, float_complex *w_MVDR)
 Generates a powermap based on the energy of adaptive Minimum-Variance Distortion-less Response (MVDR) beamformers.
 
void generateCroPaCLCMVmap (int order, float_complex *Cx, float_complex *Y_grid, int nGrid_dirs, float regPar, float lambda, float *pmap)
 (EXPERIMENTAL) Generates a powermap utilising the CroPaC LCMV post-filter described in [1]
 
void generateMUSICmap (int order, float_complex *Cx, float_complex *Y_grid, int nSources, int nGrid_dirs, int logScaleFlag, float *pmap)
 Generates an activity-map based on the sub-space multiple-signal classification (MUSIC) method.
 
void generateMinNormMap (int order, float_complex *Cx, float_complex *Y_grid, int nSources, int nGrid_dirs, int logScaleFlag, float *pmap)
 Generates an activity-map based on the sub-space minimum-norm (MinNorm) method.
 
void arraySHTmatrices (ARRAY_SHT_OPTIONS method, int order, float amp_thresh_dB, float_complex *H_array, float *grid_dirs_deg, int nBins, int nMics, int nGrid, float *w_grid, float_complex *H_sht)
 Computes matrices required to transform array signals into spherical harmonic signals (frequency-domain)
 
void arraySHTfilters (ARRAY_SHT_OPTIONS method, int order, float amp_thresh_dB, float_complex *H_array, float *grid_dirs_deg, int nFFT, int nMics, int nGrid, float *w_grid, float *h_sht)
 Computes filters required to transform array signals into spherical harmonic signals (time-domain)
 
void arraySHTmatricesDiffEQ (float_complex *H_sht, float_complex *DCM, float *freqVector, float alias_freq_hz, int nBins, int order, int nMics, float_complex *H_sht_eq)
 Diffuse-field equalisation of SHT matrices above the spatial aliasing frequency.
 
void cylModalCoeffs (int order, double *kr, int nBands, ARRAY_CONSTRUCTION_TYPES arrayType, double_complex *b_N)
 Calculates the modal coefficients for open/rigid cylindrical arrays.
 
float sphArrayAliasLim (float r, float c, int maxN)
 Returns a simple estimate of the spatial aliasing limit (the kR = maxN rule)
 
void sphArrayNoiseThreshold (int maxN, int Nsensors, float r, float c, ARRAY_CONSTRUCTION_TYPES arrayType, double dirCoeff, float maxG_db, float *f_lim)
 Computes the frequencies (per order), at which the noise of a SHT of a SMA exceeds a specified maximum level.
 
void sphModalCoeffs (int order, double *kr, int nBands, ARRAY_CONSTRUCTION_TYPES arrayType, double dirCoeff, double_complex *b_N)
 Calculates the modal coefficients for open/rigid spherical arrays.
 
void sphScattererModalCoeffs (int order, double *kr, double *kR, int nBands, double_complex *b_N)
 Calculates the modal coefficients for a rigid spherical scatterer with omni-directional sensors.
 
void sphScattererDirModalCoeffs (int order, double *kr, double *kR, int nBands, double dirCoeff, double_complex *b_N)
 Calculates the modal coefficients for a rigid spherical scatterer with directional sensors.
 
void sphDiffCohMtxTheory (int order, float *sensor_dirs_rad, int N_sensors, ARRAY_CONSTRUCTION_TYPES arrayType, double dirCoeff, double *kr, int nBands, double *M_diffcoh)
 Calculates the theoretical diffuse coherence matrix for a spherical array.
 
void diffCohMtxMeas (float_complex *H_array, int nBins, int N_sensors, int nGrid, float *w_grid, float_complex *M_diffcoh)
 Calculates the diffuse coherence matrices for an arbitrary array.
 
void diffCohMtxMeasReal (float *H_array, int N_sensors, int nGrid, float *w_grid, float *M_diffcoh)
 Calculates the diffuse coherence matrices for an array that uses a broad-band real-valued basis.
 
void simulateCylArray (int order, double *kr, int nBands, float *sensor_dirs_rad, int N_sensors, float *src_dirs_deg, int N_srcs, ARRAY_CONSTRUCTION_TYPES arrayType, float_complex *H_array)
 Simulates a cylindrical microphone array, returning the transfer functions for each (plane wave) source direction on the surface of the cylinder.
 
void simulateSphArray (int order, double *kr, double *kR, int nBands, float *sensor_dirs_rad, int N_sensors, float *src_dirs_deg, int N_srcs, ARRAY_CONSTRUCTION_TYPES arrayType, double dirCoeff, float_complex *H_array)
 Simulates a spherical microphone array, returning the transfer functions for each (plane wave) source direction on the surface of the sphere.
 
void evaluateSHTfilters (int order, float_complex *M_array2SH, int nSensors, int nBands, float_complex *H_array, int nDirs, float_complex *Y_grid, float *cSH, float *lSH)
 Generates some objective measures, which evaluate the performance of spatial encoding filters.
 

Detailed Description

Spherical harmonic domain processing module

Macro Definition Documentation

◆ NSH2ORDER

#define NSH2ORDER ( nSH)    ( (int)(sqrt((double)nSH)-0.999) )

Converts number of spherical harmonic components to spherical harmonic order i.e: sqrt(nSH)-1.

Definition at line 55 of file saf_sh.h.

◆ ORDER2NSH

#define ORDER2NSH ( order)    ((order+1)*(order+1))

Converts spherical harmonic order to number of spherical harmonic components i.e: (order+1)^2.

Definition at line 51 of file saf_sh.h.

Enumeration Type Documentation

◆ ARRAY_CONSTRUCTION_TYPES

Microphone/Hydrophone array construction types.

Enumerator
ARRAY_CONSTRUCTION_OPEN 

Open array, omni-directional sensors.

ARRAY_CONSTRUCTION_OPEN_DIRECTIONAL 

Open array, directional sensors.

ARRAY_CONSTRUCTION_RIGID 

Rigid baffle, omni-directional sensors.

ARRAY_CONSTRUCTION_RIGID_DIRECTIONAL 

Rigid baffle, directional sensors.

Definition at line 64 of file saf_sh.h.

◆ ARRAY_SHT_OPTIONS

Microphone array to spherical harmonic domain conversion options.

See also
[1] Jin, C.T., Epain, N. and Parthy, A., 2014. Design, optimization and evaluation of a dual-radius spherical microphone array. IEEE/ACM Transactions on Audio, Speech, and Language Processing, 22(1), pp.193-204.
Enumerator
ARRAY_SHT_DEFAULT 

The default SHT filters are ARRAY_SHT_REG_LS.

ARRAY_SHT_REG_LS 

Regularised least-squares (LS)

ARRAY_SHT_REG_LSHD 

Regularised least-squares (LS) in the SH domain (similar as in [1])

Definition at line 98 of file saf_sh.h.

◆ SECTOR_PATTERNS

Sector pattern designs for directionally-constraining sound-fields [1].

See also
[1] Politis, A., & Pulkki, V. (2016). Acoustic intensity, energy-density and diffuseness estimation in a directionally-constrained region. arXiv preprint arXiv:1609.03409
Enumerator
SECTOR_PATTERN_PWD 

Plane-wave decomposition/Hyper-cardioid.

SECTOR_PATTERN_MAXRE 

Spatially tapered hyper-cardioid, such that it has maximum energy concentrated in the look- direction.

SECTOR_PATTERN_CARDIOID 

Cardioid pattern.

Definition at line 81 of file saf_sh.h.

Function Documentation

◆ arraySHTfilters()

void arraySHTfilters ( ARRAY_SHT_OPTIONS method,
int order,
float amp_thresh_dB,
float_complex * H_array,
float * grid_dirs_deg,
int nFFT,
int nMics,
int nGrid,
float * w_grid,
float * h_sht )

Computes filters required to transform array signals into spherical harmonic signals (time-domain)

Parameters
[in]methodSee ARRAY_SHT_OPTIONS
[in]orderTransform order
[in]amp_thresh_dBMaximum gain amplification (in dB)
[in]H_arrayArray TFs; FLAT: (nFFT/2+1) x nMics x nGrid
[in]grid_dirs_degGrid directions [azi ELEV] degrees; FLAT: nGrid x 2
[in]nFFTFFT size / filter length
[in]nMicsNumber of microphones
[in]nGridNumber of directions in the grid
[in]w_gridIntegration weights (set to NULL if not available)
[out]h_shtThe SHT filters; FLAT: nSH x nMics x nFFT

Definition at line 2154 of file saf_sh.c.

◆ arraySHTmatrices()

void arraySHTmatrices ( ARRAY_SHT_OPTIONS method,
int order,
float amp_thresh_dB,
float_complex * H_array,
float * grid_dirs_deg,
int nBins,
int nMics,
int nGrid,
float * w_grid,
float_complex * H_sht )

Computes matrices required to transform array signals into spherical harmonic signals (frequency-domain)

Parameters
[in]methodSee ARRAY_SHT_OPTIONS
[in]orderTransform order
[in]amp_thresh_dBMaximum gain amplification (in dB)
[in]H_arrayArray TFs; FLAT: nBins x nMics x nGrid
[in]grid_dirs_degGrid directions [azi ELEV] degrees; FLAT: nGrid x 2
[in]nBinsNumber of frequencies
[in]nMicsNumber of microphones
[in]nGridNumber of directions in the grid
[in]w_gridIntegration weights (set to NULL if not available)
[out]H_shtThe SHT matrices; FLAT: nBins x nSH x nMics

Definition at line 1977 of file saf_sh.c.

◆ arraySHTmatricesDiffEQ()

void arraySHTmatricesDiffEQ ( float_complex * H_sht,
float_complex * DCM,
float * freqVector,
float alias_freq_hz,
int nBins,
int order,
int nMics,
float_complex * H_sht_eq )

Diffuse-field equalisation of SHT matrices above the spatial aliasing frequency.

Parameters
[out]H_shtInput SHT matrices; FLAT: nBins x nSH x nMics
[in]DCMDiffuse coh matrices; FLAT: nBins x nMics x nMics
[in]freqVectorFrequency vector; nBins x 1
[in]alias_freq_hzSpatial aliasing frequency, in Hz
[in]nBinsNumber of frequencies
[in]orderTransform order
[in]nMicsNumber of microphones in array
[out]H_sht_eqEqualised SHT matrices; FLAT: nBins x nSH x nMics

Definition at line 2195 of file saf_sh.c.

◆ beamWeightsCardioid2Spherical()

void beamWeightsCardioid2Spherical ( int N,
float * b_n )

Generates spherical coefficients for generating cardioid beampatterns.

For a specific order N of a higher order cardioid of the form D(theta)=(1/2)^N * (1+cos(theta))^N, this function generates the beamweights for the same pattern, but in the SHD. Because the pattern is axisymmetric only the N+1 coefficients of m=0 are returned.

Parameters
[in]NOrder of spherical harmonic expansion
[out]b_nBeamformer weights; (N+1) x 1

Definition at line 822 of file saf_sh.c.

◆ beamWeightsDolphChebyshev2Spherical()

void beamWeightsDolphChebyshev2Spherical ( int N,
int paramType,
float arrayParam,
float * b_n )

Generates beamweights in the SHD for Dolph-Chebyshev beampatterns, with mainlobe and sidelobe control [1].

Because the pattern is axisymmetric only the N+1 coefficients of m=0 are returned.

Warning
NOT IMPLEMENTED YET!
Parameters
[in]NOrder of spherical harmonic expansion
[in]paramType'0' side-lobe level control, '1' mainlobe width control
[in]arrayParamSidelobe level 1/R or mainlobe with 2*a0
[out]b_nBeamformer weights; (N+1) x 1
See also
[1] Koretz, A. and Rafaely, B., 2009. Dolph-Chebyshev beampattern design for spherical arrays. IEEE Transactions on Signal Processing, 57(6), pp.2417-2420.

◆ beamWeightsHypercardioid2Spherical()

void beamWeightsHypercardioid2Spherical ( int N,
float * b_n )

Generates beamweights in the SHD for hypercardioid beampatterns.

The hypercardioid is the pattern that maximises the directivity-factor for a certain SH order (N). The hypercardioid is also the plane-wave decomposition beamformer in the SHD, also called 'regular' because the beamweights are just the SH values on the beam-direction. Since the pattern is axisymmetric only the N+1 coefficients of m=0 are returned.

Parameters
[in]NOrder of spherical harmonic expansion
[out]b_nBeamformer weights; (N+1) x 1

Definition at line 839 of file saf_sh.c.

◆ beamWeightsMaxEV()

void beamWeightsMaxEV ( int N,
float * b_n )

Generates beamweights in the SHD for maximum energy-vector beampatterns.

Generate the beamweights for the a maximum energy-vector beampattern in the SHD. This pattern originates from Ambisonic-related research as it maximises the ambisonic energy-vector, which is essentially the directional centroid of the squared pattern. It can also be seen as the pattern which maximizes the acoustic intensity vector of a diffuse field weighted with this pattern. In practice it is almost the same as a supercardioid which maximizes front- back power ratio for a certain order, and it can be used as such. Because the pattern is axisymmetric only the N+1 coefficients of m=0 are returned. Details for their theory can be found, for example, in [1].

Parameters
[in]NOrder of spherical harmonic expansion
[out]b_nBeamformer weights; (N+1) x 1
See also
[1] Zotter, F., Pomberger, H. and Noisternig, M., 2012. Energy- preserving ambisonic decoding. Acta Acustica united with Acustica, 98(1), pp.37-47.

Definition at line 857 of file saf_sh.c.

◆ beamWeightsVelocityPatternsComplex()

void beamWeightsVelocityPatternsComplex ( int order,
float * b_n,
float azi_rad,
float elev_rad,
float_complex * A_xyz,
float_complex * velCoeffs )

Generates beamforming coefficients for velocity patterns (COMPLEX)

If the sound-field is weighted with an axisymmetric spatial distribution described by the N+1 SH coefficients b_n, then the beamweights capturing the velocity signals for the weighted sound-field are of an order one higher than the weighting pattern, and can be derived from it. This type of beamforming has some applications for spatial sound reproduction and acoustic analysis, see [1].

Parameters
[in]orderOrder of spherical harmonic expansion
[in]b_nAxisymmetric beamformer weights; (order+1) x 1
[in]azi_radOrientation, azimuth in RADIANS
[in]elev_radOrientation, ELEVATION in RADIANS
[in]A_xyzVelocity coefficients; see computeVelCoeffsMtx(); FLAT: (sectorOrder+2)^2 x (sectorOrder+1)^2 x 3
[out]velCoeffsBeamforming coefficients for velocity patterns; FLAT: (order+2)^2 x 3
See also
[1] Politis, A. and Pulkki, V., 2016. Acoustic intensity, energy-density and diffuseness estimation in a directionally-constrained region. arXiv preprint arXiv:1609.03409.

Definition at line 905 of file saf_sh.c.

◆ beamWeightsVelocityPatternsReal()

void beamWeightsVelocityPatternsReal ( int order,
float * b_n,
float azi_rad,
float elev_rad,
float_complex * A_xyz,
float * velCoeffs )

Generates beamforming coefficients for velocity patterns (REAL)

If the sound-field is weighted with an axisymmetric spatial distribution described by the N+1 SH coefficients b_n, then the beamweights capturing the velocity signals for the weighted sound-field are of an order one higher than the weighting pattern, and can be derived from it. This type of beamforming has some applications for spatial sound reproduction and acoustic analysis, see [1].

Parameters
[in]orderOrder of spherical harmonic expansion
[in]b_nAxisymmetric beamformer weights; (order+1) x 1
[in]azi_radOrientation, azimuth in RADIANS
[in]elev_radOrientation, ELEVATION in RADIANS
[in]A_xyzVelocity coefficients; see computeVelCoeffsMtx(); FLAT: (sectorOrder+2)^2 x (sectorOrder+1)^2 x 3
[out]velCoeffsBeamforming coefficients for velocity patterns; FLAT: (order+2)^2 x 3
See also
[1] Politis, A. and Pulkki, V., 2016. Acoustic intensity, energy-density and diffuseness estimation in a directionally-constrained region. arXiv preprint arXiv:1609.03409.

Definition at line 884 of file saf_sh.c.

◆ calculateGridWeights()

int calculateGridWeights ( float * dirs_rad,
int nDirs,
int order,
float * w )

Computes approximation of quadrature/integration weights for a given grid.

Test
test__calculateGridWeights
Parameters
[in]dirs_radGrid directions [azi, INCLINATION] convention, in RADIANS; FLAT: nDirs x 2
[in]nDirsNumber of directions in the grid
[in]orderSupported spherical harmonic order, or -1 to calculate
[out]wIntegration weights; nDirs x 1
Returns
Supported spherical harmonic order, or 0 if not found

Definition at line 1065 of file saf_sh.c.

◆ checkCondNumberSHTReal()

void checkCondNumberSHTReal ( int order,
float * dirs_rad,
int nDirs,
float * w,
float * cond_N )

Computes the condition numbers for a least-squares SHT.

Test
test__checkCondNumberSHTReal()
Parameters
[in]orderOrder of spherical harmonic expansion
[in]dirs_radDirections on the sphere [azi, INCLINATION] convention, in RADIANS; FLAT: nDirs x 2
[in]nDirsNumber of directions
[in]wIntegration weights; nDirs x 1
[out]cond_NCondition numbers; (order+1) x 1

Definition at line 990 of file saf_sh.c.

◆ complex2realCoeffs()

void complex2realCoeffs ( int order,
float_complex * C_N,
int K,
float * R_N )

Converts SH coefficients from the complex to real basis.

Parameters
[in]orderOrder of spherical harmonic expansion
[in]C_NComplex coeffients; FLAT: (order+1)^2 x K
[in]KNumber of columns
[out]R_NReal coefficients; FLAT: (order+1)^2 x K

Definition at line 554 of file saf_sh.c.

◆ complex2realSHMtx()

void complex2realSHMtx ( int order,
float_complex * T_c2r )

Computes a complex to real spherical harmonic transform matrix.

Computes the unitary transformation matrix T_c2r. It expresses the real spherical harmonics with respect to the complex harmonics, so that r_N = T_c2r * y_N, where r_N and y_N is are the real and complex SH vectors, respectively.

Test
test__complex2realSHMtx()
Parameters
[in]orderOrder of spherical harmonic expansion
[out]T_c2rTransformation matrix for complex->real; FLAT: (order+1)^2 x (order+1)^2

Definition at line 490 of file saf_sh.c.

◆ computeSectorCoeffsAP()

float computeSectorCoeffsAP ( int orderSec,
float_complex * A_xyz,
SECTOR_PATTERNS pattern,
float * sec_dirs_deg,
int nSecDirs,
float * sectorCoeffs )

Computes beamforming matrices (sector coefficients) which, when applied to input SH signals, yield amplitude-preserving (EP) sectors.

This partitioning of the sound-field into spatially-localised sectors has been used e.g. for parametric sound-field reproduction in [1] and visualisation in [2,3].

Note
Each sector comprises 1x sector pattern of order "orderSec", and 3x weighted dipoles which are essentially the product of the sector pattern with (unweighted) dipoles, and have the directivity of one higher order (orderSec+1).
Parameters
[in]orderSecOrder of sector patterns
[in]A_xyzVelocity coefficients, see computeVelCoeffsMtx(); FLAT: (sectorOrder+2)^2 x (sectorOrder+1)^2 x 3
[in]patternSee SECTOR_PATTERNS enum for the options
[in]sec_dirs_degSector directions [azi elev], in DEGREES; FLAT: nSecDirs x 2
[in]nSecDirsNumber of sectors
[out]sectorCoeffsThe sector coefficients; FLAT: (nSecDirs*4) x (orderSec+2)^2
Returns
Normalisation coefficient
See also
[1] Politis, A., Vilkamo, J., & Pulkki, V. (2015). Sector-based parametric sound field reproduction in the spherical harmonic domain. IEEE Journal of Selected Topics in Signal Processing, 9(5), 852-866.
[2] McCormack, L., Politis, A., and Pulkki, V. (2019). "Sharpening of angular spectra based on a directional re-assignment approach for ambisonic sound-field visualisation". IEEE International Conference on Acoustics, Speech and Signal Processing (ICASSP).
[3] McCormack, L., Delikaris-Manias, S., Politis, A., Pavlidi, D., Farina, A., Pinardi, D. and Pulkki, V., 2019. Applications of Spatially Localized Active-Intensity Vectors for Sound-Field Visualization. Journal of the Audio Engineering Society, 67(11), pp.840-854.

Definition at line 769 of file saf_sh.c.

◆ computeSectorCoeffsEP()

float computeSectorCoeffsEP ( int orderSec,
float_complex * A_xyz,
SECTOR_PATTERNS pattern,
float * sec_dirs_deg,
int nSecDirs,
float * sectorCoeffs )

Computes beamforming matrices (sector coefficients) which, when applied to input SH signals, yield energy-preserving (EP) sectors.

This partitioning of the sound-field into spatially-localised sectors has been used e.g. for parametric sound-field reproduction in [1] and visualisation in [2,3].

Note
Each sector comprises 1x sector pattern of order "orderSec", and 3x weighted dipoles which are essentially the product of the sector pattern with (unweighted) dipoles, and have the directivity of one higher order (orderSec+1).
Test
test__computeSectorCoeffsEP()
Parameters
[in]orderSecOrder of sector patterns
[in]A_xyzVelocity coefficients, see computeVelCoeffsMtx(); FLAT: (sectorOrder+2)^2 x (sectorOrder+1)^2 x 3
[in]patternSee SECTOR_PATTERNS enum for the options
[in]sec_dirs_degSector directions [azi elev], in DEGREES; FLAT: nSecDirs x 2
[in]nSecDirsNumber of sectors
[out]sectorCoeffsThe sector coefficients; FLAT: (nSecDirs*4) x (orderSec+2)^2
Returns
Normalisation coefficient
See also
[1] Politis, A., Vilkamo, J., & Pulkki, V. (2015). Sector-based parametric sound field reproduction in the spherical harmonic domain. IEEE Journal of Selected Topics in Signal Processing, 9(5), 852-866.
[2] McCormack, L., Politis, A., and Pulkki, V. (2019). "Sharpening of angular spectra based on a directional re-assignment approach for ambisonic sound-field visualisation". IEEE International Conference ' on Acoustics, Speech and Signal Processing (ICASSP).
[3] McCormack, L., Delikaris-Manias, S., Politis, A., Pavlidi, D., Farina, A., Pinardi, D. and Pulkki, V., 2019. Applications of Spatially Localized Active-Intensity Vectors for Sound-Field Visualization. Journal of the Audio Engineering Society, 67(11), pp.840-854.

Definition at line 703 of file saf_sh.c.

◆ computeVelCoeffsMtx()

void computeVelCoeffsMtx ( int sectorOrder,
float_complex * A_xyz )

Computes the matrices which generate the coefficients of a beampattern of order (sectorOrder+1) that is essentially the product of a pattern of order=sectorOrder and a dipole.

It is used in beamWeightsVelocityPatterns(). For a derivation of the matrices refer to [1].

Test
test__computeSectorCoeffsEP()
Parameters
[in]sectorOrderOrder of patterns
[out]A_xyzVelocity coefficients; FLAT: (sectorOrder+2)^2 x (sectorOrder+1)^2 x 3
See also
[1] Politis, A. and Pulkki, V., 2016. Acoustic intensity, energy-density and diffuseness estimation in a directionally-constrained region. arXiv preprint arXiv:1609.03409

Definition at line 671 of file saf_sh.c.

◆ cylModalCoeffs()

void cylModalCoeffs ( int order,
double * kr,
int nBands,
ARRAY_CONSTRUCTION_TYPES arrayType,
double_complex * b_N )

Calculates the modal coefficients for open/rigid cylindrical arrays.

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[out]b_NModal coefficients per kr and 0:order; FLAT: nBands x (order+1)

Definition at line 2266 of file saf_sh.c.

◆ diffCohMtxMeas()

void diffCohMtxMeas ( float_complex * H_array,
int nBins,
int N_sensors,
int nGrid,
float * w_grid,
float_complex * M_diffcoh )

Calculates the diffuse coherence matrices for an arbitrary array.

Parameters
[in]H_arrayArray TFs; FLAT: nBins x N_sensors x nGrid
[in]nBinsNumber of frequencies
[in]N_sensorsNumber of sensors
[in]nGridNumber of directions
[in]w_gridIntegration weights (set to NULL if not available)
[out]M_diffcohDiffuse coherence matrix per frequency; FLAT: nBands x N_sensors x N_sensors

Definition at line 2646 of file saf_sh.c.

◆ diffCohMtxMeasReal()

void diffCohMtxMeasReal ( float * H_array,
int N_sensors,
int nGrid,
float * w_grid,
float * M_diffcoh )

Calculates the diffuse coherence matrices for an array that uses a broad-band real-valued basis.

Parameters
[in]H_arrayArray TFs; FLAT: N_sensors x nGrid
[in]N_sensorsNumber of sensors
[in]nGridNumber of directions
[in]w_gridIntegration weights (set to NULL if not available)
[out]M_diffcohDiffuse coherence matrix; N_sensors x N_sensors

Definition at line 2683 of file saf_sh.c.

◆ evaluateSHTfilters()

void evaluateSHTfilters ( int order,
float_complex * M_array2SH,
int nSensors,
int nBands,
float_complex * H_array,
int nDirs,
float_complex * Y_grid,
float * cSH,
float * lSH )

Generates some objective measures, which evaluate the performance of spatial encoding filters.

This analysis is performed by comparing the spatial resolution of the spherical harmonic components generated by the encoding filters, with the ideal SH components. For more information, the reader is directed to [1,2].

Parameters
[in]orderTransform/encoding order
[in]M_array2SHEncoding matrix per frequency; FLAT: nBands x (order+1)^2 x nSensors
[in]nSensorsNumber of sensors
[in]nBandsNumber of frequency bands/bins
[in]H_arrayMeasured/modelled array responses for many directions; FLAT: nBands x nSensors x nDirs
[in]nDirsNumber of directions the array was measured/modelled
[in]Y_gridSpherical harmonics weights for each grid direction; FLAT: nDirs x (order+1)^2
[out]cSHAbsolute values of the spatial correlation per band and order; FLAT: nBands x (order+1)
[out]lSHLevel difference per band and order; FLAT: nBands x (order+1)
See also
[1] Moreau, S., Daniel, J., Bertet, S., 2006, 3D sound field recording with higher order ambisonics–objective measurements and validation of spherical microphone. In Audio Engineering Society Convention 120.
[2] Politis, A., Gamper, H. (2017). "Comparing Modelled And Measurement- Based Spherical Harmonic Encoding Filters For Spherical Microphone Arrays. In IEEE Workshop on Applications of Signal Processing to Audio and Acoustics (WASPAA).

Definition at line 2846 of file saf_sh.c.

◆ generateCroPaCLCMVmap()

void generateCroPaCLCMVmap ( int order,
float_complex * Cx,
float_complex * Y_grid,
int nGrid_dirs,
float regPar,
float lambda,
float * pmap )

(EXPERIMENTAL) Generates a powermap utilising the CroPaC LCMV post-filter described in [1]

The spatial post-filter is estimated for all directions on the grid, and is used to supress reverb/noise interference that may be present in an MVDR map. Unlike in the paper, the second column for the contraints 'A', is Y.*diag(Cx), rather than utilising a maximum energy beamformer. The post- filters are then applied to the MVDR powermap map derived in the sherical harmonic domain, rather than an MVDR beamformer generated directly in the microphone array signal domain, like in the paper. Otherwise, the algorithm is the same.

Parameters
[in]orderAnalysis order
[in]CxCorrelation/covariance matrix; FLAT: (order+1)^2 x (order+1)^2
[in]Y_gridSteering vectors for each grid direcionts; FLAT: (order+1)^2 x nGrid_dirs
[in]nGrid_dirsNumber of grid directions
[in]regParRegularisation parameter, for diagonal loading of Cx
[in]lambdaParameter controlling how harsh CroPaC is applied, 0..1; 0: fully CroPaC, 1: fully MVDR
[out]pmapResulting CroPaC LCMV powermap; nGrid_dirs x 1
See also
[1] Delikaris-Manias, S., Vilkamo, J., & Pulkki, V. (2016). Signal- dependent spatial filtering based on weighted-orthogonal beamformers in the spherical harmonic domain. IEEE/ACM Transactions on Audio, Speech and Language Processing (TASLP), 24(9), 1507-1519.

Definition at line 1764 of file saf_sh.c.

◆ generateMinNormMap()

void generateMinNormMap ( int order,
float_complex * Cx,
float_complex * Y_grid,
int nSources,
int nGrid_dirs,
int logScaleFlag,
float * pmap )

Generates an activity-map based on the sub-space minimum-norm (MinNorm) method.

Parameters
[in]orderAnalysis order
[in]CxCorrelation/covariance matrix; FLAT: (order+1)^2 x (order+1)^2
[in]Y_gridSteering vectors for each grid direcionts; FLAT: (order+1)^2 x nGrid_dirs
[in]nSourcesNumber of sources present in sound scene
[in]nGrid_dirsNumber of grid directions
[in]logScaleFlag'1' log(pmap), '0' pmap.
[out]pmapResulting MinNorm pseudo-spectrum; nGrid_dirs x 1

Definition at line 1916 of file saf_sh.c.

◆ generateMUSICmap()

void generateMUSICmap ( int order,
float_complex * Cx,
float_complex * Y_grid,
int nSources,
int nGrid_dirs,
int logScaleFlag,
float * pmap )

Generates an activity-map based on the sub-space multiple-signal classification (MUSIC) method.

Parameters
[in]orderAnalysis order
[in]CxCorrelation/covariance matrix; FLAT: (order+1)^2 x (order+1)^2
[in]Y_gridSteering vectors for each grid direcionts; FLAT: (order+1)^2 x nGrid_dirs
[in]nSourcesNumber of sources present in sound scene
[in]nGrid_dirsNumber of grid directions
[in]logScaleFlag'1' log(pmap), '0' pmap.
[out]pmapResulting MUSIC pseudo-spectrum; nGrid_dirs x 1

Definition at line 1868 of file saf_sh.c.

◆ generateMVDRmap()

void generateMVDRmap ( int order,
float_complex * Cx,
float_complex * Y_grid,
int nGrid_dirs,
float regPar,
float * pmap,
float_complex * w_MVDR )

Generates a powermap based on the energy of adaptive Minimum-Variance Distortion-less Response (MVDR) beamformers.

Parameters
[in]orderAnalysis order
[in]CxCorrelation/covariance matrix; FLAT: (order+1)^2 x (order+1)^2
[in]Y_gridSteering vectors for each grid direcionts; FLAT: (order+1)^2 x nGrid_dirs
[in]nGrid_dirsNumber of grid directions
[in]regParRegularisation parameter, for diagonal loading of Cx
[out]pmapResulting MVDR powermap; nGrid_dirs x 1
[out]w_MVDR(Optional) weights will be copied to this, unless it's NULL; FLAT: nSH x nGrid_dirs || NULL

Definition at line 1698 of file saf_sh.c.

◆ generatePWDmap()

void generatePWDmap ( int order,
float_complex * Cx,
float_complex * Y_grid,
int nGrid_dirs,
float * pmap )

Generates a powermap based on the energy of a plane-wave decomposition (PWD) (i.e.

hyper-cardioid) beamformers

Parameters
[in]orderAnalysis order
[in]CxCorrelation/covariance matrix; FLAT: (order+1)^2 x (order+1)^2
[in]Y_gridSteering vectors for each grid direcionts; FLAT: (order+1)^2 x nGrid_dirs
[in]nGrid_dirsNumber of grid directions
[out]pmapResulting PWD powermap; nGrid_dirs x 1

Definition at line 1656 of file saf_sh.c.

◆ getSHcomplex()

void getSHcomplex ( int order,
float * dirs_rad,
int nDirs,
float_complex * Y )

Computes complex-valued spherical harmonics [1] for each given direction on the unit sphere.

The real spherical harmonics are computed WITH the 1/sqrt(4*pi) term. This function employs unnorm_legendreP() and double precision.

Warning
This function assumes [azi, inclination] convention! Note that one may convert from elevation, with: [azi, pi/2-elev].
Test
test__getSHcomplex()
Parameters
[in]orderOrder of spherical harmonic expansion
[in]dirs_radDirections on the sphere [azi, INCLINATION] convention, in RADIANS; FLAT: nDirs x 2
[in]nDirsNumber of directions
[out]YThe SH weights [WITH the 1/sqrt(4*pi)]; FLAT: (order+1)^2 x nDirs
See also
[1] Rafaely, B. (2015). Fundamentals of spherical array processing (Vol. 8). Berlin: Springer.

Definition at line 439 of file saf_sh.c.

◆ getSHreal()

void getSHreal ( int order,
float * dirs_rad,
int nDirs,
float * Y )

Computes real-valued spherical harmonics [1] for each given direction on the unit sphere.

The spherical harmonic values are computed WITH the 1/sqrt(4*pi) term. Compared to getSHreal_recur(), this function uses unnorm_legendreP() and double precision, so is more suitable for being computed in an initialisation stage. This version is indeed slower, but more precise (especially for high orders).

Warning
This function assumes [azi, inclination] convention! Note that one may convert from elevation, with: [azi, pi/2-elev].
Test
test__getSHreal()
Parameters
[in]orderOrder of spherical harmonic expansion
[in]dirs_radDirections on the sphere [azi, INCLINATION] convention, in RADIANS; FLAT: nDirs x 2
[in]nDirsNumber of directions
[out]YThe SH weights [WITH the 1/sqrt(4*pi)]; FLAT: (order+1)^2 x nDirs
See also
[1] Rafaely, B. (2015). Fundamentals of spherical array processing (Vol. 8). Berlin: Springer.

Definition at line 190 of file saf_sh.c.

◆ getSHreal_part()

void getSHreal_part ( int order_start,
int order_end,
float * dirs_rad,
int nDirs,
float * Y )

Definition at line 331 of file saf_sh.c.

◆ getSHreal_recur()

void getSHreal_recur ( int order,
float * dirs_rad,
int nDirs,
float * Y )

Computes real-valued spherical harmonics [1] for each given direction on the unit sphere.

The real spherical harmonics are computed WITH the 1/sqrt(4*pi) term. Compared to getSHreal(), this function uses unnorm_legendreP_recur() and single precision, so is more suitable for being computed in a real-time loop. It sacrifices some precision, and numerical error propogates through the recursion, but it is much faster.

The function also uses static memory buffers for single direction and up to 7th order, which speeds things up considerably for such use cases.

Warning
This function assumes [azi, inclination] convention! Note that one may convert from elevation, with: [azi, pi/2-elev].
Test
test__getSHreal_recur()
Parameters
[in]orderOrder of spherical harmonic expansion
[in]dirs_radDirections on the sphere [azi, INCLINATION] convention, in RADIANS; FLAT: nDirs x 2
[in]nDirsNumber of directions
[out]YThe SH weights [WITH the 1/sqrt(4*pi)]; FLAT: (order+1)^2 x nDirs
See also
[1] Rafaely, B. (2015). Fundamentals of spherical array processing (Vol. 8). Berlin: Springer.

Definition at line 253 of file saf_sh.c.

◆ getSHrotMtxReal()

void getSHrotMtxReal ( float R[3][3],
float * RotMtx,
int L )

Generates a real-valued spherical harmonic rotation matrix [1] based on a 3x3 rotation matrix (see quaternion2rotationMatrix(), euler2rotationMatrix())

The rotation should then be applied as:

outSig = RotMtx * inSig; % where inSig/outSig are: (L+1)^2 x signalLength
Note
The normalisation convention does not matter, since only dipoles are used to rotate dipoles, quadrapoles to rotate quadrapoles etc. So any order-dependent scaling is irrelevant.
Warning
The resulting rotation matrix should be applied to signals which follow the ACN channel ordering convention!
Test
test__getSHrotMtxReal()
Parameters
[in]RThe 3x3 rotation matrix
[in]LOrder of spherical harmonic expansion
[out]RotMtxSH domain rotation matrix; FLAT: (L+1)^2 x (L+1)^2
See also
[1] Ivanic, J., Ruedenberg, K. (1998). Rotation Matrices for Real Spherical Harmonics. Direct Determination by Recursion Page: Additions and Corrections. Journal of Physical Chemistry A, 102(45), 9099–9100.

Definition at line 585 of file saf_sh.c.

◆ real2complexSHMtx()

void real2complexSHMtx ( int order,
float_complex * T_r2c )

Computes a real to complex spherical harmonic transform matrix.

Computes the unitary transformation matrix T_r2c the expresses the complex spherical harmonics with respect to the real harmonics, so that y_N = T_r2c * r_N, where r_N and y_N are the real and complex SH vectors, respectively.

Test
test__real2complexSHMtx()
Parameters
[in]orderOrder of spherical harmonic expansion
[out]T_r2cTransformation matrix for real->complex; FLAT: (order+1)^2 x (order+1)^2

Definition at line 522 of file saf_sh.c.

◆ rotateAxisCoeffsComplex()

void rotateAxisCoeffsComplex ( int order,
float * c_n,
float theta_0,
float phi_0,
float_complex * c_nm )

Generates spherical coefficients for a rotated axisymmetric pattern (COMPLEX)

Parameters
[in]orderOrder of spherical harmonic expansion
[in]c_nCoefficients describing a rotationally symmetric pattern order N, expressed as a sum of spherical harmonics of degree m=0; (N+1) x 1
[in]theta_0POLAR rotation for the pattern, in RADIANS
[in]phi_0Azimuthal rotation for the pattern, in RADIANS
[out]c_nmCoefficients of rotated pattern expressed as a sum of SHs; (N+1)^2 x 1

Definition at line 965 of file saf_sh.c.

◆ rotateAxisCoeffsReal()

void rotateAxisCoeffsReal ( int order,
float * c_n,
float theta_0,
float phi_0,
float * c_nm )

Generates spherical coefficients for a rotated axisymmetric pattern (REAL)

Parameters
[in]orderOrder of spherical harmonic expansion
[in]c_nCoefficients describing a rotationally symmetric pattern order N, expressed as a sum of spherical harmonics of degree m=0; (N+1) x 1
[in]theta_0POLAR rotation for the pattern, in RADIANS
[in]phi_0Azimuthal rotation for the pattern, in RADIANS
[out]c_nmCoefficients of rotated pattern expressed as a sum of SHs; (N+1)^2 x 1

Definition at line 945 of file saf_sh.c.

◆ simulateCylArray()

void simulateCylArray ( int order,
double * kr,
int nBands,
float * sensor_dirs_rad,
int N_sensors,
float * src_dirs_deg,
int N_srcs,
ARRAY_CONSTRUCTION_TYPES arrayType,
float_complex * H_array )

Simulates a cylindrical microphone array, returning the transfer functions for each (plane wave) source direction on the surface of the cylinder.

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[in]sensor_dirs_radSpherical coords of the sensors in RADIANS, [azi ELEV]; FLAT: N_sensors x 2
[in]N_sensorsNumber of sensors
[in]src_dirs_degSpherical coords of the plane waves in DEGREES, [azi ELEV]; FLAT: N_srcs x 2
[in]N_srcsNumber sources (DoAs of plane waves)
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[out]H_arraySimulated array response for each plane wave; FLAT: nBands x N_sensors x N_srcs

Definition at line 2716 of file saf_sh.c.

◆ simulateSphArray()

void simulateSphArray ( int order,
double * kr,
double * kR,
int nBands,
float * sensor_dirs_rad,
int N_sensors,
float * src_dirs_deg,
int N_srcs,
ARRAY_CONSTRUCTION_TYPES arrayType,
double dirCoeff,
float_complex * H_array )

Simulates a spherical microphone array, returning the transfer functions for each (plane wave) source direction on the surface of the sphere.

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*array_radius; nBands x 1
[in]kRwavenumber*scatterer_radius, set to NULL if not needed
[in]nBandsNumber of frequency bands/bins
[in]sensor_dirs_radSpherical coords of the sensors in RADIANS, [azi ELEV]; FLAT: N_sensors x 2
[in]N_sensorsNumber of sensors
[in]src_dirs_degSpherical coords of the plane waves in DEGREES, [azi ELEV]; FLAT: N_srcs x 2
[in]N_srcsNumber sources (DoAs of plane waves)
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[in]dirCoeffOnly for directional (open) arrays, 1: omni, 0.5: card, 0:dipole
[out]H_arraySimulated array response for each plane wave; FLAT: nBands x N_sensors x N_srcs

Definition at line 2768 of file saf_sh.c.

◆ sphArrayAliasLim()

float sphArrayAliasLim ( float r,
float c,
int maxN )

Returns a simple estimate of the spatial aliasing limit (the kR = maxN rule)

Parameters
[in]rArray radius, meters
[in]cSpeed of sound, m/s
[in]maxNOrder
Returns
Spatial aliasing limit estimate, in Hz

Definition at line 2331 of file saf_sh.c.

◆ sphArrayNoiseThreshold()

void sphArrayNoiseThreshold ( int maxN,
int Nsensors,
float r,
float c,
ARRAY_CONSTRUCTION_TYPES arrayType,
double dirCoeff,
float maxG_db,
float * f_lim )

Computes the frequencies (per order), at which the noise of a SHT of a SMA exceeds a specified maximum level.

Computes the frequencies that the noise in the output channels of a spherical microphone array (SMA), after performing the spherical harmonic transform (SHT) and equalisation of the output signals, reaches a certain user-defined threshold maxG_db. The frequencies are computed only at the lower range of each order, where its response decays rapidly, ignoring for example the nulls of an open array at the higher frequencies. The estimation of the limits are based on a linear approximation of the log-log response found e.g. in [1]

Parameters
[in]maxNMaximum order of the array
[in]NsensorsNumber of sensors
[in]rMic radius, meters
[in]cSpeed of sound, m/s
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[in]dirCoeffOnly for directional (open) arrays, 1: omni, 0.5: card, 0:dipole
[in]maxG_dbMax allowed amplification for the noise level, maxG_db = 20*log10(maxG)
[out]f_limNoise limit estimate; (maxN+1) x 1
See also
[1] Politis, A., Vilkamo, J., & Pulkki, V. (2015). Sector-based parametric sound field reproduction in the spherical harmonic domain. IEEE Journal of Selected Topics in Signal Processing, 9(5), 852-866.

Definition at line 2341 of file saf_sh.c.

◆ sphDiffCohMtxTheory()

void sphDiffCohMtxTheory ( int order,
float * sensor_dirs_rad,
int N_sensors,
ARRAY_CONSTRUCTION_TYPES arrayType,
double dirCoeff,
double * kr,
int nBands,
double * M_diffcoh )

Calculates the theoretical diffuse coherence matrix for a spherical array.

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]sensor_dirs_radSpherical coords of the sensors in RADIANS, [azi ELEV]; FLAT: N_sensors x 2
[in]N_sensorsNumber of sensors
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[in]dirCoeffOnly for directional (open) arrays, 1: omni, 0.5: card, 0:dipole
[in]krwavenumber*sensor_radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[out]M_diffcohTheoretical diffuse coherence matrix per frequency; FLAT: N_sensors x N_sensors x nBands

Definition at line 2569 of file saf_sh.c.

◆ sphESPRIT_create()

void sphESPRIT_create ( void **const phESPRIT,
int order )

Creates an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.

The ESPRIT method (in this case, using spherical harmonic input signals) returns the analysed DoAs directly; i.e. without any grid searching/ scanning (like e.g. MUSIC requires...). The DoA estimates are therefore continuous, i.e. not bound to any grid.

This particular implementation is is based on the "3-recurrence relationship" design, detailed in [1].

Test
test__sphESPRIT()
Parameters
[in]phESPRIT(&) address of the ESPRIT DoA estimator handle
[in]orderSpherical harmonic input order
See also
[1] B. Jo and J.-W. Choi, "Parametric direction-of-arrival estimation with three recurrence relations of spherical harmonics," J. Acoust. Soc. Amer.,vol. 145, no. 1, pp. 480–488, Jan. 2019.

Definition at line 1420 of file saf_sh.c.

◆ sphESPRIT_destroy()

void sphESPRIT_destroy ( void **const phESPRIT)

Destroys an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.

Parameters
[in]phESPRIT(&) address of the ESPRIT DoA estimator handle

Definition at line 1495 of file saf_sh.c.

◆ sphESPRIT_estimateDirs()

void sphESPRIT_estimateDirs ( void *const hESPRIT,
float_complex * Us,
int K,
float * src_dirs_rad )

Estimates the direction-of-arrival (DoA) based on the ESPRIT-based estimator, in the spherical harmonic domain.

Note
The "signal subspace" refers to the first K eigenvectors of the spatial correlation matrix, after sorting them such that the eigenvalue are in descending order.
Warning
The number of sources (K) cannot exceed: order^2!
Parameters
[in]hESPRITThe ESPRIT DoA estimator handle
[in]UsSignal subspace; FLAT: (order+1)^2 x K
[in]KNumber of sources
[out]src_dirs_radSource directions, in radians; K x 2

Definition at line 1545 of file saf_sh.c.

◆ sphModalCoeffs()

void sphModalCoeffs ( int order,
double * kr,
int nBands,
ARRAY_CONSTRUCTION_TYPES arrayType,
double dirCoeff,
double_complex * b_N )

Calculates the modal coefficients for open/rigid spherical arrays.

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[in]arrayTypeSee ARRAY_CONSTRUCTION_TYPES enum
[in]dirCoeffOnly for directional (open) arrays, 1: omni, 0.5: card, 0:dipole
[out]b_NModal coefficients per kr and 0:order; FLAT: nBands x (order+1)

Definition at line 2369 of file saf_sh.c.

◆ sphMUSIC_compute()

void sphMUSIC_compute ( void *const hMUSIC,
float_complex * Vn,
int nSrcs,
float * P_music,
int * peak_inds )

Computes a pseudo-spectrum based on the MUSIC algorithm in the spherical harmonic domain; optionally returning the grid indices corresponding to the N highest peaks (N=nSrcs)

Warning
The number of sources should not exceed: floor(nSH/2)!
Parameters
[in]hMUSICsphMUSIC handle
[in]VnNoise subspace; FLAT: nSH x (nSH - nSrcs)
[in]nSrcsNumber of sources (or an estimate of the number of sources)
[in]P_musicPseudo-spectrum (set to NULL if not wanted); nDirs x 1
[in]peak_indsIndices corresponding to the "nSrcs" highest peaks in the pseudo-spectrum (set to NULL if not wanted); nSrcs x 1

Definition at line 1355 of file saf_sh.c.

◆ sphMUSIC_create()

void sphMUSIC_create ( void **const phMUSIC,
int order,
float * grid_dirs_deg,
int nDirs )

Creates an instance of the spherical harmonic domain MUSIC implementation, which may be used for computing pseudo-spectrums for visualisation/DoA estimation purposes.

Note
Subspace approaches such as MUSIC can offer higher spatial resolution than beamforming approaches, such as the steered-response power (PWD), as long as the source signals are not correlated between them and are presented in a reverberant/diffuse sound.
Test
test__sphMUSIC()
Parameters
[in]phMUSIC(&) address of the sphMUSIC handle
[in]orderSpherical harmonic input order
[in]grid_dirs_degScanning grid directions; FLAT: nDirs x 2
[in]nDirsNumber of scanning directions

Definition at line 1284 of file saf_sh.c.

◆ sphMUSIC_destroy()

void sphMUSIC_destroy ( void **const phMUSIC)

Destroys an instance of the spherical harmonic domain MUSIC implementation.

Parameters
[in]phMUSIC(&) address of the sphMUSIC handle

Definition at line 1332 of file saf_sh.c.

◆ sphPWD_compute()

void sphPWD_compute ( void *const hPWD,
float_complex * Cx,
int nSrcs,
float * P_map,
int * peak_inds )

Computes a power-map based on determining the energy of hyper-cardioid beamformers; optionally, also returning the grid indices corresponding to the N highest peaks (N=nSrcs)

Parameters
[in]hPWDsphPWD handle
[in]CxSignal covariance matrix; FLAT: nSH x nSH
[in]nSrcsNumber of sources (or an estimate of the number of sources), for the optional peak finding (peak_inds)
[in]P_mapPowermap (set to NULL if not wanted); nDirs x 1
[in]peak_indsIndices corresponding to the "nSrcs" highest peaks in the power-map (set to NULL if not wanted); nSrcs x 1

Definition at line 1221 of file saf_sh.c.

◆ sphPWD_create()

void sphPWD_create ( void **const phPWD,
int order,
float * grid_dirs_deg,
int nDirs )

Creates an instance of a spherical harmonic domain implementation of the steer-response power (SRP) approach for computing power-maps, which can then be used for sound-field visualisation/DoA estimation purposes.

Parameters
[in]phPWD(&) address of the sphPWD handle
[in]orderSpherical harmonic input order
[in]grid_dirs_degScanning grid directions; FLAT: nDirs x 2
[in]nDirsNumber of scanning directions

Definition at line 1154 of file saf_sh.c.

◆ sphPWD_destroy()

void sphPWD_destroy ( void **const phPWD)

Destroys an instance of the spherical harmonic domain PWD implementation.

Parameters
[in]phPWD(&) address of the sphPWD handle

Definition at line 1200 of file saf_sh.c.

◆ sphScattererDirModalCoeffs()

void sphScattererDirModalCoeffs ( int order,
double * kr,
double * kR,
int nBands,
double dirCoeff,
double_complex * b_N )

Calculates the modal coefficients for a rigid spherical scatterer with directional sensors.

Assumes all sensors are placed the same distance from the scatterer, w.r.t. the origin

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*array_radius; nBands x 1
[in]kRwavenumber*scatterer_radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[in]dirCoeffDirectivity coefficient, 1: omni, 0.5: card, 0:dipole
[out]b_NModal coefficients per kr and 0:order; FLAT: nBands x (order+1)

Definition at line 2502 of file saf_sh.c.

◆ sphScattererModalCoeffs()

void sphScattererModalCoeffs ( int order,
double * kr,
double * kR,
int nBands,
double_complex * b_N )

Calculates the modal coefficients for a rigid spherical scatterer with omni-directional sensors.

Assumes all sensors are placed the same distance from the scatterer, w.r.t. the origin

Parameters
[in]orderMax order (highest is ~30 given numerical precision)
[in]krwavenumber*array_radius; nBands x 1
[in]kRwavenumber*scatterer_radius; nBands x 1
[in]nBandsNumber of frequency bands/bins
[out]b_NModal coefficients per kr and 0:order; FLAT: nBands x (order+1)

Definition at line 2453 of file saf_sh.c.

◆ unnorm_legendreP()

void unnorm_legendreP ( int n,
double * x,
int lenX,
double * y )

Calculates unnormalised Legendre polynomials up to order N, for all values in vector x [1].

Note
This includes the Condon-Shortley phase term. It is functionally identical to Matlab's legendre function (with default settings ['unnorm']).
Parameters
[in]nOrder of Legendre polynomial
[in]xVector of input values; lenX x 1
[in]lenXNumber of input values
[out]yResulting unnormalised Legendre values for each x value; FLAT: (n+1) x lenX
See also
[1] M, Abramowitz., I.A. Stegun. (1965). "Handbook of Mathematical Functions: Chapter 8", Dover Publications.

Definition at line 53 of file saf_sh.c.

◆ unnorm_legendreP_recur()

void unnorm_legendreP_recur ( int n,
float * x,
int lenX,
float * Pnm_minus1,
float * Pnm_minus2,
float * Pnm )

Calculates unnormalised Legendre polynomials up to order N, for all values in vector x.

It uses a recursive approach, which makes it more suitable for computing the legendre values in a real-time loop.

Note
This does NOT include the Condon-Shortley phase term.
Parameters
[in]nOrder of Legendre polynomial
[in]xVector of input values; lenX x 1
[in]lenXNumber of input values
[in]Pnm_minus1Previous Pnm, (not used for n=1); FLAT: (n+1) x lenX
[in]Pnm_minus2Previous previous Pnm, (not used for n=0); FLAT: (n+1) x lenX
[out]PnmResulting unnormalised Legendre values for each x value; FLAT: (n+1) x lenX

Definition at line 129 of file saf_sh.c.