Calculate the weights for a spherical Dolph-Chebyshev beamformer.
The design criterion can either be a desired side-lobe attenuation or a
desired main-lobe width. Once one criterion is chosen, the other will
become a dependent property which will be chosen accordingly [1].
Parameters:
n_max (int) – Spherical harmonic order. Must be an integer greater than zero.
design_parameter (float) – This can either be the desired side-lobe attenuation or the width of
the main-lobe in radians.
design_criterion (str) – Can be either 'sidelobe' or 'mainlobe' (see
design_parameter above). The default is 'sidelobe'.
Returns:
weights – An array containing the \((n_\mathrm{max} + 1)^2\) beamformer
weights.
Return type:
array-like, float
References
Examples
Apply Dolph-Chebyshev beamformers with a side-lobe
attenuation of 50 dB to a plane wave from zero degrees azimuth.
>>> importspharpy>>> importpyfaraspf>>> importmatplotlib.pyplotasplt>>> importnumpyasnp>>>>>> # use real valued spherical harmonics of order 7>>> definition=spharpy.SphericalHarmonicDefinition(n_max=7)>>>>>> # define the direction of the incident sound field>>> # (0 degree azimuth)>>> soundfield=spharpy.SphericalHarmonics.from_definition(... definition,pf.Coordinates(1,0,0))>>> a_nm=soundfield.basis>>>>>> # beamforming>>> # a) generate 500 steering vectors>>> steering=spharpy.SphericalHarmonics.from_definition(... definition,pf.Coordinates.from_spherical_elevation(... np.linspace(0,2*np.pi,500),0,1))>>> Y_steering=steering.basis>>>>>> # b) design beamformer weights>>> R=10**(50/20)>>> d_nm=spharpy.beamforming.dolph_chebyshev_weights(... definition.n_max,R,design_criterion='sidelobe')>>>>>> # c) appling beamformer weights yields 500 beamformers>>> beamformer=np.squeeze(Y_steering@np.diag(d_nm))>>>>>> # d) apply beamformers to the soundfield>>> beamformer_output=beamformer@a_nm.T>>>>>> # plot soundfield evaluated with beamformers>>> ax=plt.axes(projection='polar')>>> ax.plot(steering.coordinates.azimuth,... 20*np.log10(np.abs(beamformer_output)))>>> ax.set_rticks([-50,-25,0])>>> ax.set_theta_zero_location('N')>>> ax.set_xlabel('Azimuth in degree')
Compute weights that maximize the front-back ratio of the beam pattern.
This is also often referred to as the super-cardioid beam pattern.
The weights are calculated from an eigenvalue problem [2]. For high
spherical harmonic orders, the eigenvalue problem may not be feasible and
a solution will not be found.
Parameters:
n_max (int) – The spherical harmonic order. Must be an integer greater than zero.
normalize (bool) – If True, the weights will be normalized such that the complex
amplitude of a plane wave is not distorted. The default is True.
Returns:
weights – An array containing the \((n_\mathrm{max} + 1)^2\) beamformer
weights.
Return type:
array-like, float
References
Examples
Apply beamformers maximizing the front-back ratio to a plane wave from
zero degrees azimuth.
>>> importspharpy>>> importpyfaraspf>>> importmatplotlib.pyplotasplt>>> importnumpyasnp>>>>>> # use real valued spherical harmonics of order 5>>> definition=spharpy.SphericalHarmonicDefinition(n_max=5)>>>>>> # define the direction of the incident sound field>>> # (0 degree azimuth)>>> soundfield=spharpy.SphericalHarmonics.from_definition(... definition,pf.Coordinates(1,0,0))>>> a_nm=soundfield.basis>>>>>> # beamforming>>> # a) generate 500 steering vectors>>> steering=spharpy.SphericalHarmonics.from_definition(... definition,pf.Coordinates.from_spherical_elevation(... np.linspace(0,2*np.pi,500),0,1))>>> Y_steering=steering.basis>>>>>> # b) design beamformer weights>>> R=10**(50/20)>>> d_nm=spharpy.beamforming.maximum_front_back_ratio_weights(... definition.n_max)>>>>>> # c) appling beamformer weights yields 500 beamformers>>> beamformer=np.squeeze(Y_steering@np.diag(d_nm))>>>>>> # d) apply beamformers to the soundfield>>> beamformer_output=beamformer@a_nm.T>>>>>> # plot soundfield evaluated with beamformers>>> ax=plt.axes(projection='polar')>>> ax.plot(steering.coordinates.azimuth,... 20*np.log10(np.abs(beamformer_output)))>>> ax.set_rticks([-75,-50,-25,0])>>> ax.set_theta_zero_location('N')>>> ax.set_xlabel('Azimuth in degree')