def check_cond_number_sht(N, dirs, basisType, W=None): """ Computes the condition number for a least-squares SHT. Parameters ---------- N : int Maximum order to be tested for the given set of points. dirs : ndarray Evaluation directions. Dimension = (nDirs, 2). Directions are expected in radians, expressed in pairs [azimuth, inclination]. basisType : str Type of spherical harmonics. Either 'complex' or 'real'. W : ndarray, optional. Weights for each measurement point to condition the inversion, in case of a weighted least-squares. Dimension = (nDirs) Returns ------- cond_N : ndarray Condition number for each order. Dimension = (N+1) Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. Notes ----- Inclination is defined as the angle from zenith: inclination = pi/2-elevation TODO: implement complex basis? TODO: implement W """ _validate_int('N', N, positive=True) _validate_ndarray_2D('dirs', dirs, shape1=C - 1) _validate_string('basisType', basisType, choices=['complex', 'real']) if W is not None: _validate_ndarray_1D('W', W, size=dirs.shape[0]) # Compute the harmonic coefficients Y_N = get_sh(N, dirs, basisType) # Compute condition number for progressively increasing order up to N cond_N = np.zeros(N + 1) for n in range(N + 1): Y_n = Y_N[:, :np.power(n + 1, 2)] if W is None: YY_n = np.dot(Y_n.T, Y_n) else: # YY_n = Y_n.T * np.diag(W) * Y_n raise NotImplementedError cond_N[n] = np.linalg.cond(YY_n) return cond_N
def get_capsule_positions(mic_array_name): """ Retrieve the geometry of a selected set of microphone arrays. Parameters ---------- mic_array_name : str One of: 'eigenmike', 'ambeo' Returns ------- array_sigs : ndarray Capsule positions, in spherical coordinates (radians). Dimension = (nMic, C) Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. """ _validate_string('mic_array_name', mic_array_name, choices=['eigenmike', 'ambeo']) capsule_positions = None if mic_array_name is 'eigenmike': mic_dirs_deg = np.array([[ 0, 32, 0, 328, 0, 45, 69, 45, 0, 315, 291, 315, 91, 90, 90, 89, 180, 212, 180, 148, 180, 225, 249, 225, 180, 135, 111, 135, 269, 270, 270, 271 ], [ 21, 0, -21, 0, 58, 35, 0, -35, -58, -35, 0, 35, 69, 32, -31, -69, 21, 0, -21, 0, 58, 35, 0, -35, -58, -35, 0, 35, 69, 32, -32, -69 ]]) mic_dirs_rad = mic_dirs_deg * np.pi / 180. r = 0.042 mic_dirs_rad = np.row_stack( (mic_dirs_rad, r * np.ones(np.shape(mic_dirs_rad)[1]))) capsule_positions = mic_dirs_rad.T elif mic_array_name is 'ambeo': r = 0.015 capsule_positions = [ [np.pi / 4, np.arcsin(1. / np.sqrt(3)), r], # FLU [7 * np.pi / 4, -1 * np.arcsin(1. / np.sqrt(3)), r], # FRD [3 * np.pi / 4, -1 * np.arcsin(1. / np.sqrt(3)), r], # BLD [5 * np.pi / 4, np.arcsin(1. / np.sqrt(3)), r] ] # BRU return capsule_positions
def cyl_modal_coefs(N, kr, arrayType): """ Modal coefficients for rigid or open cylindrical array Parameters ---------- N : int Maximum spherical harmonic expansion order. kr: ndarray Wavenumber-radius product. Dimension = (l). arrayType: str 'open' or 'rigid'. Returns ------- b_N : ndarray Modal coefficients. Dimension = (l, N+1) 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. """ _validate_int('N', N) _validate_ndarray_1D('kr', kr, positive=True) _validate_string('arrayType', arrayType, choices=['open', 'rigid']) b_N = np.zeros((kr.size, N + 1), dtype='complex') for n in range(N + 1): if arrayType is 'open': b_N[:, n] = np.power(1j, n) * jv(n, kr) elif arrayType is 'rigid': jn = jv(n, kr) jnprime = jvp(n, kr, 1) hn = hankel2(n, kr) hnprime = h2vp(n, kr, 1) temp = np.power(1j, n) * (jn - (jnprime / hnprime) * hn) temp[np.where(kr == 0)] = 1 if n == 0 else 0. b_N[:, n] = temp return b_N
def test_validate_string(): # TypeError: not a list wrong_values = [1, 3.14, True, 3j, np.asarray([2.3]), None, np.nan, np.inf] for wv in wrong_values: with pytest.raises(TypeError, match='must be an instance of str'): _validate_string('vw', wv) # Value error: not a given choice choices = ['foo', 'bar', 'ambisonic'] wrong_values = ['FOO', 'ambisonics', 'random_string'] for wv in wrong_values: with pytest.raises(ValueError, match='must be one of the following'): _validate_string('vw', wv, choices=choices)
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 simulate_cyl_array(N_filt, mic_dirs_rad, src_dirs_rad, arrayType, R, N_order, fs): """ Simulate the impulse responses of a cylindrical 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. Dimension = (N_mic). src_dirs_rad: ndarray Direction of arrival of the indicent plane waves, in radians. Dimension = (N_doa). arrayType: str 'open' or 'rigid'. Target sampling rate R: float Radius of the array cylinder, in meter. N_order: int Maximum cylindrical harmonic expansion order. fs: int Sample rate. 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 cylindrical microphone array for the given directions of incident plane waves. The array type can be either 'open' for omnidirectional microphones in an open setup, or 'rigid' for omnidirectional microphones mounted on a cylinder. """ _validate_int('N_filt', N_filt, positive=True, parity='even') _validate_ndarray_1D('mic_dirs_rad', mic_dirs_rad) _validate_ndarray_1D('src_dirs_rad', src_dirs_rad) _validate_string('arrayType', arrayType, choices=['open', 'rigid']) _validate_float('R', R, positive=True) _validate_int('N_order', N_order, positive=True) _validate_int('fs', fs, positive=True) # 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.cyl_modal_coefs(N_order, kR, arrayType) # 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] 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): angle = mic_dirs_rad - src_dirs_rad[i] C = np.zeros((N_order + 1, N_mic)) for n in range(N_order + 1): # Jacobi-Anger expansion if n == 0: C[n, :] = np.ones(angle.shape) else: C[n, :] = 2 * np.cos(n * angle) h_mic[:, :, i] = np.matmul(b_Nt, C) H_mic[:, :, i] = np.matmul(b_N, C) return h_mic, H_mic
def ims_coreMtx(room, source, receiver, type, typeValue): """ Compute echogram by image source method. Parameters ---------- room : ndarray Room dimensions in cartesian coordinates. Dimension = (3) [x, y, z]. src : ndarray Source position in cartesian coordinates. Dimension = (3) [x, y, z]. rec : ndarray Receiver position in cartesian coordinates. Dimension = (3) [x, y, z]. type : str Restriction type: 'maxTime' or 'maxOrder' typeValue: int or float Value of the chosen restriction. Returns ------- reflections : echogram An Echogram instance. Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. Notes ----- `src` and `rec` positions are specified from the left ground corner of the room, using a left-handed coordinate system. `room` refers to the wall dimensions. Therefore, their values should be positive and smaller than room dimensions. _____ _ | | | | | | x ^ | | l = r[0] | | | | | | o----> - y |-----| w = r[1] """ _validate_ndarray_1D('room', room, size=C, positive=True) _validate_ndarray_1D('source', source, size=C, positive=True, limit=[np.zeros(C),room]) _validate_ndarray_1D('receiver', receiver, size=C, positive=True, limit=[np.zeros(C),room]) _validate_string('type', type, choices=['maxTime', 'maxOrder']) # Room dimensions l, w, h = room # Move source origin to the centrer of the room src = np.empty(C) src[0] = source[0] - l / 2 src[1] = w / 2 - source[1] src[2] = source[2] - h / 2 # Move receiver origin to the centrer of the room rec = np.empty(C) rec[0] = receiver[0] - l / 2 rec[1] = w / 2 - receiver[1] rec[2] = receiver[2] - h / 2 if type is 'maxOrder': maxOrder = typeValue echogram = ims_coreN(room, src, rec, maxOrder) elif type is 'maxTime': maxDelay = typeValue echogram = ims_coreT(room, src, rec, maxDelay) # Sort reflections according to propagation time idx = np.argsort(echogram.time) echogram.time = echogram.time[idx] echogram.value = echogram.value[idx] echogram.order = echogram.order[idx, :] echogram.coords = echogram.coords[idx, :] return echogram
def sph_modal_coefs(N, kr, arrayType, dirCoef=None): """ Modal coefficients for rigid or open spherical array Parameters ---------- N : int Maximum spherical harmonic expansion order. kr: ndarray Wavenumber-radius product. Dimension = (l). arrayType: str 'open', 'rigid' or 'directional'. dirCoef: float, optional Directivity coefficient of the sensor. Default to None. Returns ------- b_N : ndarray Modal coefficients. Dimension = (l, N+1) 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). """ _validate_int('N', N, positive=True) _validate_ndarray_1D('kr', kr, 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) b_N = np.zeros((kr.size, N + 1), dtype='complex') for n in range(N + 1): if arrayType is 'open': b_N[:, n] = 4 * np.pi * np.power(1j, n) * spherical_jn(n, kr) elif arrayType is 'rigid': jn = spherical_jn(n, kr) jnprime = spherical_jn(n, kr, derivative=True) hn = sph_hankel2(n, kr) hnprime = dsph_hankel2(n, kr) temp = 4 * np.pi * np.power(1j, n) * (jn - (jnprime / hnprime) * hn) temp[np.where(kr == 0)] = 4 * np.pi if n == 0 else 0. b_N[:, n] = temp elif arrayType is 'directional': jn = spherical_jn(n, kr) jnprime = spherical_jn(n, kr, derivative=True) temp = 4 * np.pi * np.power(1j, n) * (dirCoef * jn - 1j * (1 - dirCoef) * jnprime) b_N[:, n] = temp return b_N
def get_sh(N, dirs, basisType): """ Get spherical harmonics up to order N evaluated at given angular positions. Parameters ---------- N : int Maximum spherical harmonic expansion order. dirs : ndarray Evaluation directions. Dimension = (nDirs, 2). Directions are expected in radians, expressed in pairs [azimuth, inclination]. basisType : str Type of spherical harmonics. Either 'complex' or 'real'. Returns ------- Y : ndarray Spherical harmonics . Dimension = (nDirs, nHarm). Raises ----- TypeError, ValueError: if method arguments mismatch in type, dimension or value. NotImplementedError: for basisType = 'complex' Notes ----- Ouput dimension is given by: nHarm = (N+1)^2 Inclination is defined as the angle from zenith: inclination = pi/2-elevation TODO: implement complex basis? """ _validate_int('order', N, positive=True) _validate_ndarray_2D('dirs', dirs, shape1=C - 1) _validate_string('basisType', basisType, choices=['complex', 'real']) nDirs = dirs.shape[0] nHarm = np.power(N + 1, 2) Y_N = np.zeros((nDirs, nHarm)) def delta_kronecker(q1, q2): return 1 if q1 == q2 else 0 # TODO # it looks like the output of shs is N3d (max 1, sqrt(3)!3) # so it needs to be scaled as * np.sqrt(4*np.pi) * [1, 1./np.sqrt(3), 1./np.sqrt(3), 1./np.sqrt(3)] def norm(m): """ TODO SN3D FACTOR, REMOVE CONDON SHOTLEY IT MUST BE MULTIPLIED BY sqrt(4PI) to be normalized to 1 :param m: :return: """ return np.power(-1, np.abs(m)) * np.sqrt(2 - delta_kronecker(0, np.abs(m))) if basisType is 'complex': # TODO raise NotImplementedError elif basisType is 'real': harm_idx = 0 for l in range(N + 1): for m in range(-l, 0): Y_N[:, harm_idx] = np.imag( sph_harm(np.abs(m), l, dirs[:, 0], dirs[:, 1])) * norm(m) harm_idx += 1 for m in range(l + 1): Y_N[:, harm_idx] = np.real( sph_harm(m, l, dirs[:, 0], dirs[:, 1])) * norm(m) harm_idx += 1 return Y_N
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