def test_validate_float():

    # TypeError: not a float (NOTE: np.nan, np.inf are actually floats!)
    wrong_values = ['1', True, 1, 3j, [2.3], None, np.asarray([0.5])]
    for wv in wrong_values:
        with pytest.raises(TypeError, match='must be an instance of float'):
            _validate_float('vw', wv)

    # ValueError: not positive
    wrong_values = [-1.1, -3.14, -1e-100]
    for wv in wrong_values:
        with pytest.raises(ValueError, match='must be positive'):
            _validate_float('vw', wv, positive=True)

    # ValueError: greater than limit
    limit = 5.
    wrong_values = [5 + 1e-5, 6., np.inf]
    for wv in wrong_values:
        with pytest.raises(ValueError, match='must be smaller than'):
            _validate_float('vw', wv, limit=limit)
Esempio n. 2
0
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
Esempio n. 3
0
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
Esempio n. 4
0
def spherical_scatterer(mic_dirs_rad, src_dirs_rad, R, N_order, N_filt, fs):
    """
    Compute the pressure due to a spherical scatterer

    The function computes the impulse responses of the pressure measured
    at some points in the field with a spherical rigid scatterer centered
    at the origin and due to incident plane waves.

    Parameters
    ----------
    mic_dirs_rad: ndarray
        Position of microphone capsules. Dimension = (N_mic, C).
        Positions are expected in radians, expressed in triplets [azimuth, elevation, distance].
    src_dirs_rad: ndarray
        Direction of arrival of the indicent plane waves. Dimension = (N_doa, C-1).
        Directions are expected in radians, expressed in pairs [azimuth, elevation].
    R: float
        Radius of the array sphere, in meter.
    N_order: int
        Maximum cylindrical harmonic expansion order.
    N_filt : int
        Number of frequencies where to compute the response. It must be even.
    fs: int
        Sample rate.

    Returns
    -------
    h_mic: ndarray
        Computed IRs in time-domain. Dimension = (N_filt, Nmic, Ndoa)
    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.

    """

    _validate_ndarray_2D('mic_dirs_rad', mic_dirs_rad, shape1=masp.C)
    _validate_ndarray_2D('src_dirs_rad', src_dirs_rad, shape1=masp.C - 1)
    _validate_float('R', R, positive=True)
    _validate_int('N_order', N_order, positive=True)
    _validate_int('N_filt', N_filt, positive=True, parity='even')
    _validate_int('fs', fs, positive=True)
    if np.any(mic_dirs_rad[:, 2] < R):
        raise ValueError(
            'mic_dirs_rad: The distance of the measurement point cannot be less than the radius:'
            + str(R))

    # Compute the frequency-dependent part of the microphone responses (radial dependence)
    K = N_filt // 2 + 1
    f = np.arange(K) * fs / N_filt
    kR = 2 * np.pi * f * R / masp.c
    N_mic = mic_dirs_rad.shape[0]
    N_doa = src_dirs_rad.shape[0]

    # Check if all microphones at same radius
    same_radius = np.sum(mic_dirs_rad[1:, 2] - mic_dirs_rad[:-1, 2]) == 0
    if same_radius:
        # Spherical modal coefs for rigid sphere
        b_N = np.zeros((K, N_order + 1), dtype='complex')
        r = mic_dirs_rad[0, 2]
        kr = 2 * np.pi * f * r / masp.c

        # Similar to the sph_modal_coefs for the rigid case
        for n in range(N_order + 1):
            jn = spherical_jn(n, kr)
            jnprime = spherical_jn(n, kR, derivative=True)
            hn = asr.sph_hankel2(n, kr)
            hnprime = asr.dsph_hankel2(n, kR)
            b_N[:,
                n] = (2 * n + 1) * np.power(1j, n) * (jn -
                                                      (jnprime / hnprime) * hn)
    else:
        # Spherical modal coefs for rigid sphere, but at different distances
        b_N = np.zeros((K, N_order + 1, N_mic), dtype='complex')
        for nm in range(N_mic):
            r = mic_dirs_rad[nm, 2]
            kr = 2 * np.pi * f * r / masp.c

            # Similar to the sph_modal_coefs for the rigid case
            for n in range(N_order + 1):
                jn = spherical_jn(n, kr)
                jnprime = spherical_jn(n, kR, derivative=True)
                hn = asr.sph_hankel2(n, kr)
                hnprime = asr.dsph_hankel2(n, kR)
                b_N[:, n, nm] = (2 * n + 1) * np.power(
                    1j, n) * (jn - (jnprime / hnprime) * hn)

    # Avoid NaNs for very high orders, instead of (very) very small values
    b_N[np.isnan(b_N)] = 0.

    # Compute angular-dependent part of the microphone responses
    H_mic = np.zeros((K, N_mic, N_doa), dtype='complex')
    for nd in range(N_doa):
        # Unit vectors of DOAs and microphones
        azi0 = src_dirs_rad[nd, 0]
        elev0 = src_dirs_rad[nd, 1]
        azi = mic_dirs_rad[:, 0]
        elev = mic_dirs_rad[:, 1]
        cosAlpha = np.sin(elev) * np.sin(elev0) + np.cos(elev) * np.cos(
            elev0) * np.cos(azi - azi0)

        P_N = np.zeros((N_order + 1, N_mic))
        for n in range(N_order + 1):
            for nm in range(N_mic):
                P_N[n, nm] = lpmn(n, n, cosAlpha[nm])[0][0, -1]

        # Accumulate across orders
        if same_radius:
            H_mic[:, :, nd] = np.matmul(b_N, P_N)
        else:
            for nm in range(N_mic):
                H_mic[:, nm, nd] = np.matmul(b_N[:, :, nm], P_N[:, nm])

    # Handle Nyquist for real impulse response
    tempH_mic = H_mic.copy()
    # TODO: in `simulate_sph_array()` it was real, not abs. Why?
    tempH_mic[-1, :] = np.abs(tempH_mic[-1, :])
    # Conjugate ifft and fftshift for causal IR
    h_mic = np.real(
        np.fft.fftshift(np.fft.ifft(np.append(tempH_mic,
                                              np.conj(tempH_mic[-2:0:-1, :]),
                                              axis=0),
                                    axis=0),
                        axes=0))

    return h_mic, H_mic
Esempio n. 5
0
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
Esempio n. 6
0
def render_quantised(qechogram, endtime, fs, fractional):
    """
    Render a quantised echogram array into a quantised impulse response matrix.

    Parameters
    ----------
    qechograms : ndarray, dtype = QuantisedEchogram
        Target quantised echograms. Dimension = (nDirs).
    endtime : float
        Maximum time of rendered reflections, in seconds.
    fs : int
        Target sampling rate
    fractional : bool, optional
        Use fractional or integer (round) delay. Default to True.

    Returns
    -------
    qIR : ndarray
        Rendered quantised echograms. Dimension = (ceil(endtime * fs), nChannels)
    idx_nonzero : 1D ndarray
        Indices of non-zero elements.

    Raises
    -----
    TypeError, ValueError: if method arguments mismatch in type, dimension or value.

    TODO: expose fractional as parameter?
    """

    _validate_quantised_echogram_array(qechogram)
    _validate_float('edntime', endtime, positive=True)
    _validate_int('fs', fs, positive=True)
    _validate_boolean('fractional', fractional)

    # Number of grid directions of the quantization
    nDirs = qechogram.size

    # Render echograms
    L_rir = int(np.ceil(endtime * fs))
    qIR = np.zeros((L_rir, nDirs))
    idx_nonzero = []
    for nq in range(nDirs):
        tempgram = Echogram(
            time=qechogram[nq].time,
            value=qechogram[nq].value,
            order=np.zeros((qechogram[nq].time.size, 3),
                           dtype=int),  # whatever to pass the size validation
            coords=np.zeros((qechogram[nq].time.size,
                             3)))  # whatever to pass the size validation
        # Omit if there are no echoes in the specific one
        if qechogram[nq].isActive:
            idx_nonzero.append(nq)
            # Number of reflections inside the time limit')
            idx_limit = tempgram.time[tempgram.time < endtime].size
            tempgram.time = tempgram.time[:idx_limit + 1]
            tempgram.value = tempgram.value[:idx_limit + 1]
            tempgram.order = tempgram.order[:idx_limit +
                                            1]  # whatever to pass the size validation
            tempgram.coords = tempgram.coords[:idx_limit +
                                              1]  # whatever to pass the size validation

            qIR[:, nq] = render_rirs(tempgram, endtime, fs,
                                     fractional).squeeze()

    return qIR, np.asarray(idx_nonzero)
Esempio n. 7
0
def render_rirs(echogram, endtime, fs, fractional=True):
    """
    Render an echogram into an impulse response.

    Parameters
    ----------
    echogram : Echogram
        Target Echogram.
    endtime : float
        Maximum time of rendered reflections, in seconds.
    fs : int
        Target sampling rate
    fractional : bool, optional
        Use fractional or integer (round) delay. Default to True.

    Returns
    -------
    ir : ndarray
        Rendered echogram. Dimension = (ceil(endtime * fs), nChannels)

    Raises
    -----
    TypeError, ValueError: if method arguments mismatch in type, dimension or value.

    Notes
    -----
    TODO: expose filter order as parameter?
    """

    _validate_echogram(echogram)
    _validate_float('endtime', endtime, positive=True)
    _validate_int('fs', fs, positive=True)
    _validate_boolean('fractional', fractional)

    nChannels = 1 if np.ndim(echogram.value) <= 1 else np.shape(
        echogram.value)[1]
    L_ir = int(np.ceil(endtime * fs))
    ir = np.zeros((L_ir, nChannels))
    # Number of reflections inside the time limit
    idx_trans = echogram.time[echogram.time < endtime].size

    if fractional:
        # Get lagrange interpolating filter of order 100 (filter length 101)
        order = 100
        h_offset = 50
        h_idx = np.arange(-(order / 2), (order / 2) + 1).astype(
            int)  # only valid for even orders

        # Make a filter table for quick access for quantized fractional samples
        fractions = np.linspace(0, 1, 101)
        H_frac = lagrange(order, 50 + fractions)

        # Initialise array
        tmp_ir = np.zeros((int(L_ir + (2 * h_offset)), nChannels))

        for i in range(idx_trans):
            refl_idx = int(np.floor(echogram.time[i] * fs) + 1)
            refl_frac = np.remainder(echogram.time[i] * fs, 1)
            filter_idx = np.argmin(np.abs(refl_frac - fractions))
            h_frac = H_frac[:, filter_idx]

            tmp_ir[h_offset + refl_idx + h_idx -
                   1, :] += h_frac[:, np.newaxis] * echogram.value[i]

        ir = tmp_ir[h_offset:-h_offset, :]

    else:
        refl_idx = (np.round(echogram.time[:idx_trans] * fs)).astype(int)
        # Filter out exceeding indices
        refl_idx = refl_idx[refl_idx < L_ir]
        ir[refl_idx, :] = echogram.value[:refl_idx.size]

    return ir
Esempio n. 8
0
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
Esempio n. 9
0
def array_sht_filters_measure_regLSHD(H_array,
                                      order_sht,
                                      grid_dirs_rad,
                                      w_grid=None,
                                      nFFT=1024,
                                      amp_threshold=10.):
    """
    Generate SHT filters based on measured responses  (regularized least-squares in the SHD)

    Parameters
    ----------
    H_array : ndarray, dtype=complex
        Frequency domain measured array responses. Dimension = ( nFFT//2+1, nMics, nGrids )
    order_sht : int
        Spherical harmonic transform order.
    grid_dirs_rad: ndarray,
        Grid positions in [azimuth, elevation] pairs (in radians). Dimension = (nGrid, 2).
    w_grid : ndarray, optional
        Weights for weighted-least square solution, based on the importance or area
        around each measurement point (leave empty if not known, or not important)
        Dimension = ( nGrid ).
    nFFT : int, optional
        Number of points for the FFT.
    amp_threshold : float
        Max allowed amplification for filters, in dB.

    Returns
    -------
    h_filt : ndarray
        Impulse responses of the filters. Dimension = ( (order_sht+1)^2, nMics, nFFT ).
    H_filt: ndarray
        Frequency-domain filters. Dimension = ( (order_sht+1)^2, nMics, nFFT//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 micorphone signals from a spherical
    microphone array to SH signals, based on a least-squares solution with
    a constraint on noise amplification, using Tikhonov regularization. The
    method formulates the LS problem in the space domain, using the
    directional measurements of the array response, similar to e.g.

        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.

    # TODO: nFFT argument is redundant!!!
    Due to the matrix nature of computations,
    the minimum valid value for `nMics` and `nGrid` is 2 and 16 respectively.
    """

    _validate_int('nFFT', nFFT, positive=True, parity='even')
    nBins = nFFT // 2 + 1

    _validate_ndarray_3D('H_array', H_array, shape0=nBins)
    nMic = H_array.shape[1]
    _validate_number('nMic', nMic, limit=[2, np.inf])
    nGrid = H_array.shape[2]
    _validate_number('nGrid', nGrid, limit=[16, np.inf])

    _validate_int('order_sht', order_sht, positive=True)
    _validate_ndarray_2D('grid_dirs_rad', grid_dirs_rad, shape1=C - 1)

    if w_grid is None:
        w_grid = 1 / nGrid * np.ones(nGrid)
    _validate_ndarray_1D('w_grid', w_grid, size=nGrid)

    _validate_float('amp_threshold', amp_threshold)

    # 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)

    order_array = int(np.floor(np.sqrt(nGrid) / 2 - 1))
    # TODO: check validity of the approach
    # order_array must be greater or equal than requested order_sht ( nGrid > (2*(order_sht+1))^2 )
    if order_array < order_sht:
        raise ValueError("Order array < Order SHT. Consider increasing nGrid")

    # SH matrix at grid directions
    Y_grid = np.sqrt(
        4 * np.pi) * get_sh(order_array, elev2incl(grid_dirs_rad),
                            'real').T  # SH matrix for grid directions

    # Compute inverse matrix
    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)
    W_grid = np.diag(w_grid)
    H_nm = np.zeros((nBins, nMic, np.power(order_array + 1, 2)),
                    dtype='complex')
    for kk in range(nBins):
        tempH = H_array[kk, :, :]
        H_nm[kk, :, :] = np.matmul(
            np.matmul(np.matmul(tempH, W_grid), Y_grid.T),
            np.linalg.inv(np.matmul(np.matmul(Y_grid, W_grid), Y_grid.T)))

    # Compute the inverse matrix in the SHD with regularization
    H_filt = np.zeros((np.power(order_sht + 1, 2), nMic, nBins),
                      dtype='complex')
    for kk in range(nBins):
        tempH_N = H_nm[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
Esempio n. 10
0
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
Esempio n. 11
0
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
Esempio n. 12
0
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
Esempio n. 13
0
def sph_array_alias_lim(R, nMic, maxN, mic_dirs_rad, mic_weights=None):
    """
    Get estimates of the aliasing limit of a spherical array, in three different ways.

    Parameters
    ----------
    R : float
        Microphone array radius, in meter.
    nMic : int
        Number of microphone capsules.
    maxN : int
        Maximum spherical harmonic expansion order.
    mic_dirs_rad : ndarray
        Evaluation directions. Dimension = (nDirs, 2).
        Directions are expected in radians, expressed in pairs [azimuth, elevation].
    dirCoef: ndarray, optional
       Vector of weights used to improve orthogonality of the SH transform. Dimension = (nDirs)

    Returns
    -------
    f_alias : ndarray
       the alising limit estimates. Dimension = (3).

    Raises
    -----
    TypeError, ValueError: if method arguments mismatch in type, dimension or value.

    Notes
    -----
    First estimate takes into account only the radius and a nominal order
    that the array is expected to support, it is the simplest one and it
    expresses the kR = maxN rule.
    The second estimate is based on the number of microphones, and it can
    be more relaxed than the first, if the nominal supported order is less
    than maxN<floor(sqrt(Nmic)-1).
    The third takes into account microphone numbers and directions, and it
    is based on the orthogonality of the SH matrix for the microphone
    positions expressed through the condition number.

    # todo check RETURN VALUES (need for very big rtol if cond_N returned)
    """

    _validate_float('R', R, positive=True)
    _validate_int('nMic', nMic, positive=True)
    _validate_int('maxN', maxN, positive=True)
    _validate_ndarray_2D('mic_dirs_rad', mic_dirs_rad, shape1=C - 1)
    if mic_weights is not None:
        _validate_ndarray_1D('mic_weights',
                             mic_weights,
                             size=mic_dirs_rad.shape[0])

    f_alias = np.zeros(3)

    # Conventional kR = N assumption
    f_alias[0] = c * maxN / (2 * np.pi * R)

    # Based on the floor of the number of microphones, uniform arrangement
    f_alias[1] = c * np.floor(np.sqrt(nMic) - 1) / (2 * np.pi * R)

    # Based on condition number of the SHT matrix
    maxOrder = int(np.ceil(np.sqrt(nMic) - 1))

    cond_N = check_cond_number_sht(maxOrder, elev2incl(mic_dirs_rad), 'real',
                                   mic_weights)
    trueMaxOrder = np.flatnonzero(
        cond_N < np.power(10, 4))[-1]  # biggest element passing the condition
    f_alias[2] = c * trueMaxOrder / (2 * np.pi * R)

    return f_alias
Esempio n. 14
0
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