Modal Beamforming¶
Submodules for modal beamforming
Angular¶
-
micarray.modal.angular.
sht_matrix
(N, azi, elev, weights=None)[source]¶ Matrix of spherical harmonics up to order N for given angles.
Computes a matrix of spherical harmonics up to order \(N\) for the given angles/grid.
\[\begin{split}\mathbf{Y} = \left[ \begin{array}{cccccc} Y_0^0(\theta[0], \phi[0]) & Y_1^{-1}(\theta[0], \phi[0]) & Y_1^0(\theta[0], \phi[0]) & Y_1^1(\theta[0], \phi[0]) & \dots & Y_N^N(\theta[0], \phi[0]) \\ Y_0^0(\theta[1], \phi[1]) & Y_1^{-1}(\theta[1], \phi[1]) & Y_1^0(\theta[1], \phi[1]) & Y_1^1(\theta[1], \phi[1]) & \dots & Y_N^N(\theta[1], \phi[1]) \\ \vdots & \vdots & \vdots & \vdots & \vdots & \vdots \\ Y_0^0(\theta[Q-1], \phi[Q-1]) & Y_1^{-1}(\theta[Q-1], \phi[Q-1]) & Y_1^0(\theta[Q-1], \phi[Q-1]) & Y_1^1(\theta[Q-1], \phi[Q-1]) & \dots & Y_N^N(\theta[Q-1], \phi[Q-1]) \end{array} \right]\end{split}\]where
\[Y_n^m(\theta, \phi) = \sqrt{\frac{2n + 1}{4 \pi} \frac{(n-m)!}{(n+m)!}} P_n^m(\cos \theta) e^{i m \phi}\]Parameters: - N (int) – Maximum order.
- azi ((Q,) array_like) – Azimuth.
- elev ((Q,) array_like) – Elevation.
- weights ((Q,) array_like, optional) – Quadrature weights.
Returns: Ymn (((N+1)**2, Q) numpy.ndarray) – Matrix of spherical harmonics.
-
micarray.modal.angular.
Legendre_matrix
(N, ctheta)[source]¶ Matrix of weighted Legendre Polynominals.
Computes a matrix of weighted Legendre polynominals up to order N for the given angles
\[L_n(\theta) = \frac{2n+1}{4 \pi} P_n(\theta)\]Parameters: - N (int) – Maximum order.
- ctheta ((Q,) array_like) – Angles.
Returns: Lmn (((N+1), Q) numpy.ndarray) – Matrix containing Legendre polynominals.
Radial¶
-
micarray.modal.radial.
spherical_pw
(N, k, r, setup)[source]¶ Radial coefficients for a plane wave.
Computes the radial component of the spherical harmonics expansion of a plane wave impinging on a spherical array.
\[\mathring{P}_n(k) = 4 \pi i^n b_n(kr)\]Parameters: - N (int) – Maximum order.
- k ((M,) array_like) – Wavenumber.
- r (float) – Radius of microphone array.
- setup ({‘open’, ‘card’, ‘rigid’}) – Array configuration (open, cardioids, rigid).
Returns: bn ((M, N+1) numpy.ndarray) – Radial weights for all orders up to N and the given wavenumbers.
-
micarray.modal.radial.
spherical_ps
(N, k, r, rs, setup)[source]¶ Radial coefficients for a point source.
Computes the radial component of the spherical harmonics expansion of a point source impinging on a spherical array.
\[\mathring{P}_n(k) = 4 \pi (-i) k h_n^{(2)}(k r_s) b_n(kr)\]Parameters: - N (int) – Maximum order.
- k ((M,) array_like) – Wavenumber.
- r (float) – Radius of microphone array.
- rs (float) – Distance of source.
- setup ({‘open’, ‘card’, ‘rigid’}) – Array configuration (open, cardioids, rigid).
Returns: bn ((M, N+1) numpy.ndarray) – Radial weights for all orders up to N and the given wavenumbers.
-
micarray.modal.radial.
weights
(N, kr, setup)[source]¶ Radial weighing functions.
Computes the radial weighting functions for diferent array types (cf. eq.(2.62), Rafaely 2015).
For instance for an rigid array
\[b_n(kr) = j_n(kr) - \frac{j_n^\prime(kr)}{h_n^{(2)\prime}(kr)}h_n^{(2)}(kr)\]Parameters: - N (int) – Maximum order.
- kr ((M,) array_like) – Wavenumber * radius.
- setup ({‘open’, ‘card’, ‘rigid’}) – Array configuration (open, cardioids, rigid).
Returns: bn ((M, N+1) numpy.ndarray) – Radial weights for all orders up to N and the given wavenumbers.
-
micarray.modal.radial.
regularize
(dn, a0, method)[source]¶ Regularization (amplitude limitation) of radial filters.
Amplitude limitation of radial filter coefficients, methods according to (cf. Rettberg, Spors : DAGA 2014)
Parameters: - dn (numpy.ndarray) – Values to be regularized
- a0 (float) – Parameter for regularization (not required for all methods)
- method ({‘none’, ‘discard’, ‘softclip’, ‘Tikh’, ‘wng’}) – Method used for regularization/amplitude limitation (none, discard, hardclip, Tikhonov, White Noise Gain).
Returns: - dn (numpy.ndarray) – Regularized values.
- hn (array_like)
-
micarray.modal.radial.
diagonal_mode_mat
(bk)[source]¶ Diagonal matrix of radial coefficients for all modes/wavenumbers.
Parameters: bk ((M, N+1) numpy.ndarray) – Vector containing values for all wavenumbers \(M\) and modes up to order \(N\) Returns: Bk ((M, (N+1)**2, (N+1)**2) numpy.ndarray) – Multidimensional array containing diagnonal matrices with input vector on main diagonal.