SAF
Loading...
Searching...
No Matches
saf_sh.h
Go to the documentation of this file.
1/*
2 * Copyright 2016-2018 Leo McCormack
3 *
4 * Permission to use, copy, modify, and/or distribute this software for any
5 * purpose with or without fee is hereby granted, provided that the above
6 * copyright notice and this permission notice appear in all copies.
7 *
8 * THE SOFTWARE IS PROVIDED "AS IS" AND THE AUTHOR DISCLAIMS ALL WARRANTIES WITH
9 * REGARD TO THIS SOFTWARE INCLUDING ALL IMPLIED WARRANTIES OF MERCHANTABILITY
10 * AND FITNESS. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY SPECIAL, DIRECT,
11 * INDIRECT, OR CONSEQUENTIAL DAMAGES OR ANY DAMAGES WHATSOEVER RESULTING FROM
12 * LOSS OF USE, DATA OR PROFITS, WHETHER IN AN ACTION OF CONTRACT, NEGLIGENCE OR
13 * OTHER TORTIOUS ACTION, ARISING OUT OF OR IN CONNECTION WITH THE USE OR
14 * PERFORMANCE OF THIS SOFTWARE.
15 */
16
39#ifndef __SAF_SH_H_INCLUDED__
40#define __SAF_SH_H_INCLUDED__
41
42#ifdef __cplusplus
43extern "C" {
44#endif /* __cplusplus */
45
47
51#define ORDER2NSH(order) ((order+1)*(order+1))
55#define NSH2ORDER(nSH) ( (int)(sqrt((double)nSH)-0.999) )
56
57/* ========================================================================== */
58/* Enums */
59/* ========================================================================== */
60
73
89
105
106
107/* ========================================================================== */
108/* Misc. Functions */
109/* ========================================================================== */
110
128void unnorm_legendreP(/* Input Arguments */
129 int n,
130 double* x,
131 int lenX,
132 /* Output Arguments */
133 double* y);
134
153void unnorm_legendreP_recur(/* Input Arguments */
154 int n,
155 float* x,
156 int lenX,
157 float* Pnm_minus1,
158 float* Pnm_minus2,
159 /* Output Arguments */
160 float* Pnm);
161
162
163/* ========================================================================== */
164/* SH and Beamforming related Functions */
165/* ========================================================================== */
166
192void getSHreal(/* Input Arguments */
193 int order,
194 float* dirs_rad,
195 int nDirs,
196 /* Output Arguments */
197 float* Y);
198
227void getSHreal_recur(/* Input Arguments */
228 int order,
229 float* dirs_rad,
230 int nDirs,
231 /* Output Arguments */
232 float* Y);
233
234void getSHreal_part
235(
236 int order_start,
237 int order_end,
238 float* dirs_rad,
239 int nDirs,
240 float* Y /* the SH weights: (order_end+1)^2 x nDirs */
241);
242
265void getSHcomplex(/* Input Arguments */
266 int order,
267 float* dirs_rad,
268 int nDirs,
269 /* Output Arguments */
270 float_complex* Y);
271
286void complex2realSHMtx(/* Input Arguments */
287 int order,
288 /* Output Arguments */
289 float_complex* T_c2r);
290
305void real2complexSHMtx(/* Input Arguments */
306 int order,
307 /* Output Arguments */
308 float_complex* T_r2c);
309
318void complex2realCoeffs(/* Input Arguments */
319 int order,
320 float_complex* C_N,
321 int K,
322 /* Output Arguments */
323 float* R_N);
324
351void getSHrotMtxReal(float R[3][3],
352 float* RotMtx,
353 int L);
354
373void computeVelCoeffsMtx(/* Input Arguments */
374 int sectorOrder,
375 /* Output Arguments */
376 float_complex* A_xyz);
377
418float computeSectorCoeffsEP(/* Input Arguments */
419 int orderSec,
420 float_complex* A_xyz,
421 SECTOR_PATTERNS pattern,
422 float* sec_dirs_deg,
423 int nSecDirs,
424 /* Output Arguments */
425 float* sectorCoeffs);
426
465float computeSectorCoeffsAP(/* Input Arguments */
466 int orderSec,
467 float_complex* A_xyz,
468 SECTOR_PATTERNS pattern,
469 float* sec_dirs_deg,
470 int nSecDirs,
471 /* Output Arguments */
472 float* sectorCoeffs);
473
485void beamWeightsCardioid2Spherical(/* Input Arguments */
486 int N,
487 /* Output Arguments */
488 float* b_n);
489
509void beamWeightsDolphChebyshev2Spherical(/* Input Arguments */
510 int N,
511 int paramType,
512 float arrayParam,
513 /* Output Arguments */
514 float* b_n);
515
528void beamWeightsHypercardioid2Spherical(/* Input Arguments */
529 int N,
530 /* Output Arguments */
531 float* b_n);
532
553void beamWeightsMaxEV(/* Input Arguments */
554 int N,
555 /* Output Arguments */
556 float* b_n);
557
581void beamWeightsVelocityPatternsReal(/* Input Arguments */
582 int order,
583 float* b_n,
584 float azi_rad,
585 float elev_rad,
586 float_complex* A_xyz,
587 /* Output Arguments */
588 float* velCoeffs);
589
613void beamWeightsVelocityPatternsComplex(/* Input Arguments */
614 int order,
615 float* b_n,
616 float azi_rad,
617 float elev_rad,
618 float_complex* A_xyz,
619 /* Output Arguments */
620 float_complex* velCoeffs);
621
634void rotateAxisCoeffsReal(/* Input arguments */
635 int order,
636 float* c_n,
637 float theta_0,
638 float phi_0,
639 /* Output arguments */
640 float* c_nm);
641
654void rotateAxisCoeffsComplex(/* Input arguments */
655 int order,
656 float* c_n,
657 float theta_0,
658 float phi_0,
659 /* Output arguments */
660 float_complex* c_nm);
661
674void checkCondNumberSHTReal(/* Input arguments */
675 int order,
676 float* dirs_rad,
677 int nDirs,
678 float* w,
679 /* Output arguments */
680 float* cond_N);
681
694int calculateGridWeights(/* Input arguments */
695 float* dirs_rad,
696 int nDirs,
697 int order,
698 /* Output arguments */
699 float* w);
700
701
702/* ========================================================================== */
703/* Localisation Functions in the SHD */
704/* ========================================================================== */
705
716void sphPWD_create(void ** const phPWD,
717 int order,
718 float* grid_dirs_deg,
719 int nDirs);
720
726void sphPWD_destroy(void ** const phPWD);
727
742void sphPWD_compute(/* Input arguments */
743 void* const hPWD,
744 float_complex *Cx,
745 int nSrcs,
746 /* Output arguments */
747 float* P_map,
748 int* peak_inds);
749
766void sphMUSIC_create(void ** const phMUSIC,
767 int order,
768 float* grid_dirs_deg,
769 int nDirs);
770
776void sphMUSIC_destroy(void ** const phMUSIC);
777
794void sphMUSIC_compute(/* Input arguments */
795 void* const hMUSIC,
796 float_complex *Vn,
797 int nSrcs,
798 /* Output arguments */
799 float* P_music,
800 int* peak_inds);
801
823void sphESPRIT_create(void ** const phESPRIT,
824 int order);
825
832void sphESPRIT_destroy(void ** const phESPRIT);
833
848void sphESPRIT_estimateDirs(/* Input arguments */
849 void * const hESPRIT,
850 float_complex* Us,
851 int K,
852 /* Output arguments */
853 float* src_dirs_rad);
854
867void generatePWDmap(/* Input arguments */
868 int order,
869 float_complex* Cx,
870 float_complex* Y_grid,
871 int nGrid_dirs,
872 /* Output arguments */
873 float* pmap);
874
890void generateMVDRmap(/* Input arguments */
891 int order,
892 float_complex* Cx,
893 float_complex* Y_grid,
894 int nGrid_dirs,
895 float regPar,
896 /* Output arguments */
897 float* pmap,
898 float_complex* w_MVDR);
899
929void generateCroPaCLCMVmap(/* Input arguments */
930 int order,
931 float_complex* Cx,
932 float_complex* Y_grid,
933 int nGrid_dirs,
934 float regPar,
935 float lambda,
936 /* Output arguments */
937 float* pmap);
938
953void generateMUSICmap(/* Input arguments */
954 int order,
955 float_complex* Cx,
956 float_complex* Y_grid,
957 int nSources,
958 int nGrid_dirs,
959 int logScaleFlag,
960 /* Output arguments */
961 float* pmap);
962
977void generateMinNormMap(/* Input arguments */
978 int order,
979 float_complex* Cx,
980 float_complex* Y_grid,
981 int nSources,
982 int nGrid_dirs,
983 int logScaleFlag,
984 /* Output arguments */
985 float* pmap);
986
987
988/* ========================================================================== */
989/* Microphone/Hydrophone array processing functions */
990/* ========================================================================== */
991
1007void arraySHTmatrices(/* Input arguments */
1008 ARRAY_SHT_OPTIONS method,
1009 int order,
1010 float amp_thresh_dB,
1011 float_complex* H_array,
1012 float* grid_dirs_deg,
1013 int nBins,
1014 int nMics,
1015 int nGrid,
1016 float* w_grid,
1017 /* Output arguments */
1018 float_complex* H_sht);
1019
1035void arraySHTfilters(/* Input arguments */
1036 ARRAY_SHT_OPTIONS method,
1037 int order,
1038 float amp_thresh_dB,
1039 float_complex* H_array,
1040 float* grid_dirs_deg,
1041 int nFFT,
1042 int nMics,
1043 int nGrid,
1044 float* w_grid,
1045 /* Output arguments */
1046 float* h_sht);
1047
1061void arraySHTmatricesDiffEQ(/* Input arguments */
1062 float_complex* H_sht,
1063 float_complex* DCM,
1064 float* freqVector,
1065 float alias_freq_hz,
1066 int nBins,
1067 int order,
1068 int nMics,
1069 /* Output arguments */
1070 float_complex* H_sht_eq);
1071
1082void cylModalCoeffs(/* Input arguments */
1083 int order,
1084 double* kr,
1085 int nBands,
1086 ARRAY_CONSTRUCTION_TYPES arrayType,
1087 /* Output arguments */
1088 double_complex* b_N);
1089
1098float sphArrayAliasLim(/* Input arguments */
1099 float r,
1100 float c,
1101 int maxN);
1102
1131void sphArrayNoiseThreshold(/* Input arguments */
1132 int maxN,
1133 int Nsensors,
1134 float r,
1135 float c,
1136 ARRAY_CONSTRUCTION_TYPES arrayType,
1137 double dirCoeff,
1138 float maxG_db,
1139 /* Output arguments */
1140 float* f_lim);
1141
1154void sphModalCoeffs(/* Input arguments */
1155 int order,
1156 double* kr,
1157 int nBands,
1158 ARRAY_CONSTRUCTION_TYPES arrayType,
1159 double dirCoeff,
1160 /* Output arguments */
1161 double_complex* b_N);
1162
1177void sphScattererModalCoeffs(/* Input arguments */
1178 int order,
1179 double* kr,
1180 double* kR,
1181 int nBands,
1182 /* Output arguments */
1183 double_complex* b_N);
1184
1200void sphScattererDirModalCoeffs(/* Input arguments */
1201 int order,
1202 double* kr,
1203 double* kR,
1204 int nBands,
1205 double dirCoeff,
1206 /* Output arguments */
1207 double_complex* b_N);
1208
1225void sphDiffCohMtxTheory(/* Input arguments */
1226 int order,
1227 float* sensor_dirs_rad,
1228 int N_sensors,
1229 ARRAY_CONSTRUCTION_TYPES arrayType,
1230 double dirCoeff,
1231 double* kr,
1232 int nBands,
1233 /* Output arguments */
1234 double* M_diffcoh);
1235
1247void diffCohMtxMeas(/* Input arguments */
1248 float_complex* H_array,
1249 int nBins,
1250 int N_sensors,
1251 int nGrid,
1252 float* w_grid,
1253 /* Output arguments */
1254 float_complex* M_diffcoh);
1255
1266void diffCohMtxMeasReal(/* Input arguments */
1267 float* H_array,
1268 int N_sensors,
1269 int nGrid,
1270 float* w_grid,
1271 /* Output arguments */
1272 float* M_diffcoh);
1273
1292void simulateCylArray(/* Input arguments */
1293 int order,
1294 double* kr,
1295 int nBands,
1296 float* sensor_dirs_rad,
1297 int N_sensors,
1298 float* src_dirs_deg,
1299 int N_srcs,
1300 ARRAY_CONSTRUCTION_TYPES arrayType,
1301 /* Output arguments */
1302 float_complex* H_array);
1303
1326void simulateSphArray(/* Input arguments */
1327 int order,
1328 double* kr,
1329 double* kR,
1330 int nBands,
1331 float* sensor_dirs_rad,
1332 int N_sensors,
1333 float* src_dirs_deg,
1334 int N_srcs,
1335 ARRAY_CONSTRUCTION_TYPES arrayType,
1336 double dirCoeff,
1337 /* Output arguments */
1338 float_complex* H_array);
1339
1372void evaluateSHTfilters(/* Input arguments */
1373 int order,
1374 float_complex* M_array2SH,
1375 int nSensors,
1376 int nBands,
1377 float_complex* H_array,
1378 int nDirs,
1379 float_complex* Y_grid,
1380 /* Output arguments */
1381 float* cSH,
1382 float* lSH);
1383
1384
1385#ifdef __cplusplus
1386} /* extern "C" */
1387#endif /* __cplusplus */
1388
1389#endif /* __SAF_SH_H_INCLUDED__ */
1390
1391 /* doxygen addtogroup SH */
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.
Definition saf_sh.c:1917
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)...
Definition saf_sh.c:1699
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.
Definition saf_sh.c:2196
void sphMUSIC_destroy(void **const phMUSIC)
Destroys an instance of the spherical harmonic domain MUSIC implementation.
Definition saf_sh.c:1333
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].
Definition saf_sh.c:54
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)
Definition saf_sh.c:966
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) a...
Definition saf_sh.c:1155
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.
Definition saf_sh.c:2267
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.
Definition saf_sh.c:2847
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 ...
Definition saf_sh.c:1356
void sphESPRIT_create(void **const phESPRIT, int order)
Creates an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.
Definition saf_sh.c:1421
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,...
Definition saf_sh.c:704
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.
Definition saf_sh.c:130
int calculateGridWeights(float *dirs_rad, int nDirs, int order, float *w)
Computes approximation of quadrature/integration weights for a given grid.
Definition saf_sh.c:1066
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...
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.
Definition saf_sh.c:2647
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.
Definition saf_sh.c:440
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)
Definition saf_sh.c:885
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 comp...
Definition saf_sh.c:1285
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)
Definition saf_sh.c:906
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 harmon...
Definition saf_sh.c:1546
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.
Definition saf_sh.c:2370
void beamWeightsCardioid2Spherical(int N, float *b_n)
Generates spherical coefficients for generating cardioid beampatterns.
Definition saf_sh.c:823
void beamWeightsMaxEV(int N, float *b_n)
Generates beamweights in the SHD for maximum energy-vector beampatterns.
Definition saf_sh.c:858
void sphPWD_destroy(void **const phPWD)
Destroys an instance of the spherical harmonic domain PWD implementation.
Definition saf_sh.c:1201
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)
Definition saf_sh.c:946
float sphArrayAliasLim(float r, float c, int maxN)
Returns a simple estimate of the spatial aliasing limit (the kR = maxN rule)
Definition saf_sh.c:2332
void complex2realSHMtx(int order, float_complex *T_c2r)
Computes a complex to real spherical harmonic transform matrix.
Definition saf_sh.c:491
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,...
Definition saf_sh.c:1222
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.
Definition saf_sh.c:2684
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 maximu...
Definition saf_sh.c:2342
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.
Definition saf_sh.c:1869
void beamWeightsHypercardioid2Spherical(int N, float *b_n)
Generates beamweights in the SHD for hypercardioid beampatterns.
Definition saf_sh.c:840
void complex2realCoeffs(int order, float_complex *C_N, int K, float *R_N)
Converts SH coefficients from the complex to real basis.
Definition saf_sh.c:555
ARRAY_CONSTRUCTION_TYPES
Microphone/Hydrophone array construction types.
Definition saf_sh.h:64
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) sour...
Definition saf_sh.c:2717
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-doma...
Definition saf_sh.c:1978
void sphESPRIT_destroy(void **const phESPRIT)
Destroys an instance of the spherical harmonic domain ESPRIT-based direction of arrival estimator.
Definition saf_sh.c:1496
SECTOR_PATTERNS
Sector pattern designs for directionally-constraining sound-fields [1].
Definition saf_sh.h:81
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...
Definition saf_sh.c:2769
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 qu...
Definition saf_sh.c:586
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.
Definition saf_sh.c:254
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)
Definition saf_sh.c:2155
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.
Definition saf_sh.c:1657
void computeVelCoeffsMtx(int sectorOrder, float_complex *A_xyz)
Computes the matrices which generate the coefficients of a beampattern of order (sectorOrder+1) that ...
Definition saf_sh.c:672
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,...
Definition saf_sh.c:770
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.
Definition saf_sh.c:191
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.
Definition saf_sh.c:2570
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.
Definition saf_sh.c:2454
void checkCondNumberSHTReal(int order, float *dirs_rad, int nDirs, float *w, float *cond_N)
Computes the condition numbers for a least-squares SHT.
Definition saf_sh.c:991
void real2complexSHMtx(int order, float_complex *T_r2c)
Computes a real to complex spherical harmonic transform matrix.
Definition saf_sh.c:523
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.
Definition saf_sh.c:2503
ARRAY_SHT_OPTIONS
Microphone array to spherical harmonic domain conversion options.
Definition saf_sh.h:98
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]
Definition saf_sh.c:1765
@ ARRAY_CONSTRUCTION_RIGID_DIRECTIONAL
Rigid baffle, directional sensors.
Definition saf_sh.h:70
@ ARRAY_CONSTRUCTION_RIGID
Rigid baffle, omni-directional sensors.
Definition saf_sh.h:68
@ ARRAY_CONSTRUCTION_OPEN_DIRECTIONAL
Open array, directional sensors.
Definition saf_sh.h:67
@ ARRAY_CONSTRUCTION_OPEN
Open array, omni-directional sensors.
Definition saf_sh.h:65
@ SECTOR_PATTERN_PWD
Plane-wave decomposition/Hyper-cardioid.
Definition saf_sh.h:82
@ SECTOR_PATTERN_MAXRE
Spatially tapered hyper-cardioid, such that it has maximum energy concentrated in the look- direction...
Definition saf_sh.h:83
@ SECTOR_PATTERN_CARDIOID
Cardioid pattern.
Definition saf_sh.h:86
@ ARRAY_SHT_DEFAULT
The default SHT filters are ARRAY_SHT_REG_LS.
Definition saf_sh.h:99
@ ARRAY_SHT_REG_LS
Regularised least-squares (LS)
Definition saf_sh.h:100
@ ARRAY_SHT_REG_LSHD
Regularised least-squares (LS) in the SH domain (similar as in [1])
Definition saf_sh.h:101
Contains wrappers for handling complex numbers across both C99-compliant compilers and Microsoft Visu...