def simulate_sph_array(N_filt, mic_dirs_rad, src_dirs_rad, arrayType, R, N_order, fs, dirCoef=None): """ Simulate the impulse responses of a spherical array. Parameters ---------- N_filt : int Number of frequencies where to compute the response. It must be even. mic_dirs_rad: ndarray Directions of microphone capsules, in radians. Expressed in [azi, ele] pairs. Dimension = (N_mic, C-1). src_dirs_rad: ndarray Direction of arrival of the indicent plane waves, in radians. Expressed in [azi, ele] pairs. Dimension = (N_doa, C-1). arrayType: str 'open', 'rigid' or 'directional'. Target sampling rate R: float Radius of the array sphere, in meter. N_order: int Maximum spherical harmonic expansion order. fs: int Sample rate. dirCoef: float, optional Directivity coefficient of the sensors. Default to None. Returns ------- h_mic: ndarray Computed IRs in time-domain. Dimension = (N_filt, N_mic, N_doa). H_mic: ndarray, dtype='complex' Frequency responses of the computed IRs. Dimension = (N_filt//2+1, N_mic, N_doa). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. Notes ----- This method computes the impulse responses of the microphones of a spherical microphone array for the given directions of incident plane waves. The array type can be either 'open' for omnidirectional microphones in an open setup, 'rigid' for omnidirectional microphones mounted on a sphere, or 'directional' for an open array of first-order directional microphones determined by `dirCoef`. """ _validate_int('N_filt', N_filt, positive=True, parity='even') _validate_ndarray_2D('mic_dirs_rad', mic_dirs_rad, shape1=masp.C - 1) _validate_ndarray_2D('src_dirs_rad', src_dirs_rad, shape1=masp.C - 1) _validate_string('arrayType', arrayType, choices=['open', 'rigid', 'directional']) _validate_float('R', R, positive=True) _validate_int('N_order', N_order, positive=True) _validate_int('fs', fs, positive=True) if arrayType is 'directional': if dirCoef is None: raise ValueError( 'dirCoef must be defined in the directional case.') _validate_float('dirCoef', dirCoef) # Compute the frequency-dependent part of the microphone responses (radial dependence) f = np.arange(N_filt // 2 + 1) * fs / N_filt kR = 2 * np.pi * f * R / masp.c b_N = asr.sph_modal_coefs(N_order, kR, arrayType, dirCoef) # Handle Nyquist for real impulse response temp = b_N.copy() temp[-1, :] = np.real(temp[-1, :]) # Create the symmetric conjugate negative frequency response for a real time-domain signal b_Nt = np.real( np.fft.fftshift(np.fft.ifft(np.append(temp, np.conj(temp[-2:0:-1, :]), axis=0), axis=0), axes=0)) # Compute angular-dependent part of the microphone responses # Unit vectors of DOAs and microphones N_doa = src_dirs_rad.shape[0] N_mic = mic_dirs_rad.shape[0] U_doa = masp.sph2cart( np.column_stack((src_dirs_rad[:, 0], src_dirs_rad[:, 1], np.ones(N_doa)))) U_mic = masp.sph2cart( np.column_stack((mic_dirs_rad[:, 0], mic_dirs_rad[:, 1], np.ones(N_mic)))) h_mic = np.zeros((N_filt, N_mic, N_doa)) H_mic = np.zeros((N_filt // 2 + 1, N_mic, N_doa), dtype='complex') for i in range(N_doa): cosangle = np.dot(U_mic, U_doa[i, :]) P = np.zeros((N_order + 1, N_mic)) for n in range(N_order + 1): for nm in range(N_mic): # The Legendre polynomial gives the angular dependency Pn = scipy.special.lpmn(n, n, cosangle[nm])[0][0, -1] P[n, nm] = (2 * n + 1) / (4 * np.pi) * Pn h_mic[:, :, i] = np.matmul(b_Nt, P) H_mic[:, :, i] = np.matmul(b_N, P) return h_mic, H_mic
def array_sht_filters_theory_radInverse(R, nMic, order_sht, Lfilt, fs, amp_threshold): """ Generate SHT filters based on theoretical responses (regularized radial inversion) Parameters ---------- R : float Microphone array radius, in meter. nMic : int Number of microphone capsules. order_sht : int Spherical harmonic transform order. Lfilt : int Number of FFT points for the output filters. It must be even. fs : int Sample rate for the output filters. amp_threshold : float Max allowed amplification for filters, in dB. Returns ------- h_filt : ndarray Impulse responses of the filters. Dimension = (Lfilt, order_sht+1). H_filt: ndarray Frequency-domain filters. Dimension = (Lfilt//2+1, order_sht+1). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. UserWarning: if `nMic` not big enough for the required sht order. Notes ----- Generate the filters to convert microphone signals from a spherical microphone array to SH signals, based on an ideal theoretical model of the array. The filters are generated from an inversion of the radial components of the response, neglecting spatial aliasing effects and non-ideal arrangements of the microphones. One filter is shared by all SH signals of the same order in this case. Here this single channel inversion problem is done with a constraint on max allowed amplification regularization, using Tikhonov regularization, similar e.g. to: 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. """ _validate_float('R', R, positive=True) _validate_int('nMic', nMic, positive=True) _validate_int('order_sht', order_sht, positive=True) _validate_int('Lfilt', Lfilt, positive=True, parity='even') _validate_int('fs', fs, positive=True) _validate_float('amp_threshold', amp_threshold) f = np.arange(Lfilt // 2 + 1) * fs / Lfilt # Adequate sht order to the number of microphones if order_sht > np.sqrt(nMic) - 1: order_sht = int(np.floor(np.sqrt(nMic) - 1)) warnings.warn( "Set order too high for the number of microphones, should be N<=np.sqrt(Q)-1. Auto set to " + str(order_sht), UserWarning) # Modal responses kR = 2 * np.pi * f * R / c bN = sph_modal_coefs(order_sht, kR, 'rigid') / (4 * np.pi) # Encoding matrix with regularization a_dB = amp_threshold alpha = complex(np.sqrt(nMic) * np.power( 10, a_dB / 20)) # Explicit casting to allow negative sqrt (a_dB < 0) beta = np.sqrt( (1 - np.sqrt(1 - 1 / np.power(alpha, 2))) / (1 + np.sqrt(1 - 1 / np.power(alpha, 2)))) # Moreau & Daniel # Regularized single channel equalization filters per order H_filt = (bN).conj() / (np.power(np.abs(bN), 2) + np.power(beta, 2) * np.ones( (Lfilt // 2 + 1, order_sht + 1))) # Time domain filters temp = H_filt.copy() temp[-1, :] = np.real(temp[-1, :]) # Create the symmetric conjugate negative frequency response for a real time-domain signal h_filt = np.real( np.fft.fftshift(np.fft.ifft(np.append(temp, np.conj(temp[-2:0:-1, :]), axis=0), axis=0), axes=0)) # TODO: check return ordering return h_filt, H_filt
def array_sht_filters_theory_regLS(R, mic_dirsAziElev, order_sht, Lfilt, fs, amp_threshold): """ Generate SHT filters based on theoretical responses (regularized least-squares) Parameters ---------- R : float Microphone array radius, in meter. mic_dirsAziElev : ndarray Evaluation directions. Dimension = (nDirs, 2). Directions are expected in radians, expressed in pairs [azimuth, elevation]. order_sht : int Spherical harmonic transform order. Lfilt : int Number of FFT points for the output filters. It must be even. fs : int Sample rate for the output filters. amp_threshold : float Max allowed amplification for filters, in dB. Returns ------- h_filt : ndarray Impulse responses of the filters. Dimension = ( nMic, (order_sht+1)^2, Lfilt ). H_filt: ndarray Frequency-domain filters. Dimension = ( nMic, (order_sht+1)^2, Lfilt//2+1 ). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. UserWarning: if `nMic` not big enough for the required sht order. TODO: ValueError: if the array order is not big enough. Notes ----- Generate the filters to convert microphone signals from a spherical microphone array to SH signals, based on an ideal theoretical model of the array. The filters are generated as a least-squares solution with a constraint on filter amplification, using Tikhonov regularization. The method formulates the LS problem in the spherical harmonic domain, by expressing the array response to an order-limited series of SH coefficients, similar to 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. """ _validate_float('R', R, positive=True) _validate_ndarray_2D('mic_dirsAziElev', mic_dirsAziElev, shape1=C - 1) _validate_int('order_sht', order_sht, positive=True) _validate_int('Lfilt', Lfilt, positive=True, parity='even') _validate_int('fs', fs, positive=True) _validate_float('amp_threshold', amp_threshold) f = np.arange(Lfilt // 2 + 1) * fs / Lfilt num_f = f.size kR = 2 * np.pi * f * R / c kR_max = kR[-1] nMic = mic_dirsAziElev.shape[0] # Adequate sht order to the number of microphones if order_sht > np.sqrt(nMic) - 1: order_sht = int(np.floor(np.sqrt(nMic) - 1)) warnings.warn( "Set order too high for the number of microphones, should be N<=np.sqrt(Q)-1. Auto set to " + str(order_sht), UserWarning) mic_dirsAziIncl = elev2incl(mic_dirsAziElev) order_array = int(np.floor(2 * kR_max)) # TODO: check validity of the approach if order_array <= 1: raise ValueError("Order array <= 1. Consider increasing R or fs.") Y_array = np.sqrt(4 * np.pi) * get_sh(order_array, mic_dirsAziIncl, 'real') # Modal responses bN = sph_modal_coefs(order_array, kR, 'rigid') / ( 4 * np.pi ) # due to modified SHs, the 4pi term disappears from the plane wave expansion # Array response in the SHD H_array = np.zeros((nMic, np.power(order_array + 1, 2), num_f), dtype='complex') for kk in range(num_f): temp_b = bN[kk, :].T B = np.diag(replicate_per_order(temp_b)) H_array[:, :, kk] = np.matmul(Y_array, B) a_dB = amp_threshold alpha = complex(np.power( 10, a_dB / 20)) # Explicit casting to allow negative sqrt (a_dB < 0) beta = 1 / (2 * alpha) H_filt = np.zeros((np.power(order_sht + 1, 2), nMic, num_f), dtype='complex') for kk in range(num_f): tempH_N = H_array[:, :, kk] tempH_N_trunc = tempH_N[:, :np.power(order_sht + 1, 2)] H_filt[:, :, kk] = np.matmul( tempH_N_trunc.T.conj(), np.linalg.inv( np.matmul(tempH_N, tempH_N.T.conj()) + np.power(beta, 2) * np.eye(nMic))) # Time domain filters h_filt = H_filt.copy() h_filt[:, :, -1] = np.abs(h_filt[:, :, -1]) h_filt = np.concatenate((h_filt, np.conj(h_filt[:, :, -2:0:-1])), axis=2) h_filt = np.real(np.fft.ifft(h_filt, axis=2)) h_filt = np.fft.fftshift(h_filt, axes=2) # TODO: check return ordering return h_filt, H_filt
def array_sht_filters_theory_softLim(R, nMic, order_sht, Lfilt, fs, amp_threshold): """ Generate SHT filters based on theoretical responses (soft-limiting) Parameters ---------- R : float Microphone array radius, in meter. nMic : int Number of microphone capsules. order_sht : int Spherical harmonic transform order. Lfilt : int Number of FFT points for the output filters. It must be even. fs : int Sample rate for the output filters. amp_threshold : float Max allowed amplification for filters, in dB. Returns ------- h_filt : ndarray Impulse responses of the filters. Dimension = (Lfilt, order_sht+1). H_filt: ndarray Frequency-domain filters. Dimension = (Lfilt//2+1, order_sht+1). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. UserWarning: if `nMic` not big enough for the required sht order. Notes ----- Generate the filters to convert microphone signals from a spherical microphone array to SH signals, based on an ideal theoretical model of the array. The filters are generated from an inversion of the radial components of the response, neglecting spatial aliasing effects and non-ideal arrangements of the microphones. One filter is shared by all SH signals of the same order in this case. Here this single channel inversion problem is done through a thresholding approach on the inversion, limited to a max allowed amplification. The limiting follows the approach of Bernschutz, B., Porschmann, C., Spors, S., Weinzierl, S., Versterkung, B., 2011. Soft-limiting der modalen amplitudenverst?rkung bei sph?rischen mikrofonarrays im plane wave decomposition verfahren. Proceedings of the 37. Deutsche Jahrestagung fur Akustik (DAGA 2011) """ _validate_float('R', R, positive=True) _validate_int('nMic', nMic, positive=True) _validate_int('order_sht', order_sht, positive=True) _validate_int('Lfilt', Lfilt, positive=True, parity='even') _validate_int('fs', fs, positive=True) _validate_float('amp_threshold', amp_threshold) f = np.arange(Lfilt // 2 + 1) * fs / Lfilt # Adequate sht order to the number of microphones if order_sht > np.sqrt(nMic) - 1: order_sht = int(np.floor(np.sqrt(nMic) - 1)) warnings.warn( "Set order too high for the number of microphones, should be N<=np.sqrt(Q)-1. Auto set to " + str(order_sht), UserWarning) # Modal responses kR = 2 * np.pi * f * R / c bN = sph_modal_coefs(order_sht, kR, 'rigid') / ( 4 * np.pi ) # due to modified SHs, the 4pi term disappears from the plane wave expansion # Single channel equalization filters per order inv_bN = 1. / bN inv_bN[0, 1:] = 0. # Encoding matrix with regularization a_dB = amp_threshold alpha = complex(np.sqrt(nMic) * np.power( 10, a_dB / 20)) # Explicit casting to allow negative sqrt (a_dB < 0) # Regularized single channel equalization filters per order H_filt = (2 * alpha / np.pi) * (np.abs(bN) * inv_bN) * np.arctan( (np.pi / (2 * alpha)) * np.abs(inv_bN)) # Time domain filters temp = H_filt.copy() temp[-1, :] = np.real(temp[-1, :]) # Create the symmetric conjugate negative frequency response for a real time-domain signal h_filt = np.real( np.fft.fftshift(np.fft.ifft(np.append(temp, np.conj(temp[-2:0:-1, :]), axis=0), axis=0), axes=0)) # TODO: check return order return h_filt, H_filt
def sph_array_noise(R, nMic, maxN, arrayType, f): """ Returns noise amplification curves of a spherical mic. array Parameters ---------- R : float Microphone array radius, in meter. nMic : int Number of microphone capsules. maxN : int Maximum spherical harmonic expansion order. arrayType : str 'open' or 'rigid'. f : ndarray Frequencies where to perform estimation. Dimension = (l). Returns ------- g2 : ndarray Noise amplification curve. Dimension = (l, maxN) g2_lin: ndarray Linear log-log interpolation of noise. Dimension = (l, maxN). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. Notes ----- The `arrayType` options are: - 'open' for open array of omnidirectional sensors, or - 'rigid' for sensors mounted on a rigid baffle. `g2_lin` is an approximation of the curves at low frequencies showing the linear behaviour in the log-log axis, with a 6n dB slope. """ _validate_float('R', R, positive=True) _validate_int('nMic', nMic, positive=True) _validate_int('maxN', maxN, positive=True) _validate_string('arrayType', arrayType, choices=['open', 'rigid']) _validate_ndarray_1D('f', f, positive=True) # Frequency axis kR = 2 * np.pi * f * R / c # Modal responses bN = sph_modal_coefs(maxN, kR, arrayType) / (4 * np.pi) # Noise power response g2 = 1. / (nMic * np.power(np.abs(bN), 2)) # Approximate linearly p = -(6 / 10) * np.arange(1, maxN + 1) / np.log10(2) bN_lim0 = (sph_modal_coefs(maxN, np.asarray([1]), arrayType) / (4 * np.pi)).squeeze() a = 1. / (nMic * np.power(np.abs(bN_lim0[1:]), 2)) g2_lin = np.zeros((kR.size, maxN)) for n in range(maxN): g2_lin[:, n] = a[n] * np.power(kR, p[n]) return g2, g2_lin
def sph_array_noise_threshold(R, nMic, maxG_db, maxN, arrayType, dirCoef=None): """ Returns frequency limits for noise amplification of a spherical mic. array Parameters ---------- R : float Microphone array radius, in meter. nMic : int Number of microphone capsules. maxG_db : float max allowed amplification for the noise level. maxG_db = 20*log10(maxG) maxN : int Maximum spherical harmonic expansion order. arrayType: str 'open', 'rigid' or 'directional' dirCoef: float, optional Directivity coefficient of the sensor. Default to None. Returns ------- f_lim : ndarray Frequency points where threhsold is reached, for orders 1:maxN. Dimension = (maxN). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. Notes ----- The `arrayType` options are: - 'open' for open array of omnidirectional sensors, - 'rigid' for sensors mounted on a rigid baffle, - 'directional' for an open array of first-order directional microphones determined by `dirCoef`. `dirCoef` is relevant (and required) only in the 'directional' type. `dirCoef` ranges from 0 (omni) to 1 (dipole), where for example 0.5 is a cardioid sensor. In the 0 case it is equivalent to an open array of omnis. The first order directivity function is defined as d(theta) = dirCoeff + (1-dirCoeff)*cos(theta). The method returns the frequencies that the noise in the output channels of a SMA, after performing the SHT and equalization 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 Sector-based Parametric Sound Field Reproduction in the Spherical Harmonic Domain A Politis, J Vilkamo, V Pulkki IEEE Journal of Selected Topics in Signal Processing 9 (5), 852 - 866 """ _validate_float('R', R, positive=True) _validate_int('nMic', nMic, positive=True) _validate_float('maxG_db', maxG_db) _validate_int('maxN', maxN, positive=True) _validate_string('arrayType', arrayType, choices=['open', 'rigid', 'directional']) if arrayType is 'directional': if dirCoef is None: raise ValueError( 'dirCoef must be defined in the directional case.') _validate_float('dirCoef', dirCoef) f_lim = np.zeros(maxN) for n in range(1, maxN + 1): bn = sph_modal_coefs(n, np.asarray([1]), arrayType, dirCoef) / (4 * np.pi) bn = bn.flatten()[-1] maxG = np.power(10, (maxG_db / 10)) kR_lim = np.power((maxG * nMic * np.power(np.abs(bn), 2)), (-10 * np.log10(2) / (6 * n))) f_lim[n - 1] = kR_lim * c / (2 * np.pi * R) return f_lim