fbeamform#
Implements beamformers in the frequency domain.
Basic class for implementing steering vectors with monopole source transfer models. |
|
Beamforming using the basic delay-and-sum algorithm in the frequency domain. |
|
Functional beamforming algorithm. |
|
Beamforming using the Capon (Mininimum Variance) algorithm. |
|
Beamforming using eigenvalue and eigenvector techniques. |
|
Beamforming using the MUSIC algorithm. |
|
CLEAN deconvolution algorithm. |
|
DAMAS deconvolution algorithm. |
|
DAMAS deconvolution [2] for solving the system of equations. |
|
Orthogonal deconvolution algorithm. |
|
CLEAN-SC deconvolution algorithm. |
|
Covariance Matrix Fitting algorithm. |
|
Source directivity modeling in the cross-spectral matrix (SODIX) algorithm. |
|
Beamforming GIB methods with different normalizations. |
|
Abstract base class for array methods without predefined grid. |
|
Orthogonal beamforming without predefined grid. |
|
The point spread function. |
|
|
Calculates the sound pressure level from the squared sound pressure. |
|
Integrates a sound pressure map over a given sector. |
- class acoular.fbeamform.SteeringVector#
Bases:
HasStrictTraitsBasic class for implementing steering vectors with monopole source transfer models.
Handles four different steering vector formulations. See [3] for details.
- steer_type = Enum('true level', 'true location', 'classic', 'inverse')#
Type of steering vectors, see also [3]. Defaults to âtrue levelâ.
- env = Instance(Environment(), Environment)#
Environmentor derived object, which provides information about the sound propagation in the medium. Defaults to standardEnvironmentobject.
- r0 = Property()#
Sound travel distances from microphone array center to grid : points or reference position (readonly). Feature may change.
- rm = Property()#
Sound travel distances from array microphones to grid : points (readonly). Feature may change.
- ref = Property()#
Reference position or distance at which to evaluate the sound pressure of a grid point. If set to a scalar, this is used as reference distance to the grid points. If set to a vector, this is interpreted as x,y,z coordinates of the reference position. Defaults to [0.,0.,0.].
- digest = Property(depends_on=['steer_type', 'env.digest', 'grid.digest', 'mics.digest', '_ref'])#
A unique identifier for the steering vector, based on its properties. (read-only)
- inv_digest = Property(depends_on=['env.digest', 'grid.digest', 'mics.digest', '_ref'])#
A unique identifier for the grid, excluding
steer_type. Use for inverse methods. (read-only)
- transfer(f, ind=None)#
Calculates the transfer matrix for one frequency.
- Parameters:
- ffloat
Frequency for which to calculate the transfer matrix
- ind(optional) array of ints
If set, only the transfer function of the gridpoints addressed by the given indices will be calculated. Useful for algorithms like CLEAN-SC, where not the full transfer matrix is needed
- Returns:
- array of complex128
array of shape (ngridpts, nmics) containing the transfer matrix for the given frequency
- steer_vector(f, ind=None)#
Calculates the steering vectors based on the transfer function.
See also [3].
- Parameters:
- ffloat
Frequency for which to calculate the transfer matrix
- ind(optional) array of ints
If set, only the steering vectors of the gridpoints addressed by the given indices will be calculated. Useful for algorithms like CLEAN-SC, where not the full transfer matrix is needed
- Returns:
- array of complex128
array of shape (ngridpts, nmics) containing the steering vectors for the given frequency
- class acoular.fbeamform.BeamformerBase#
Bases:
HasStrictTraitsBeamforming using the basic delay-and-sum algorithm in the frequency domain.
- steer = Instance(SteeringVector, args=())#
Instance of
SteeringVectoror its derived classes that contains information about the steering vector. This is a private trait. Do not set this directly, use steer trait instead.
- freq_data = Instance(PowerSpectra)#
PowerSpectraobject that provides the cross spectral matrix and eigenvalues
- r_diag = Bool(True)#
Boolean flag, if âTrueâ (default), the main diagonal is removed before beamforming.
- r_diag_norm = Float(0.0)#
If diagonal of the CSM is removed, some signal energy is lost. This is handled via this normalization factor. Internally, the default is: num_mics / (num_mics - 1).
If r_diag==True: if r_diag_norm==0.0, the default normalization = num_mics/(num_mics-1) is used. If r_diag_norm !=0.0, the user input is used instead. If r_diag==False, the normalization is 1.0 either way.
- precision = Enum('float64', 'float32')#
Floating point precision of property result. Corresponding to numpy dtypes. Default = 64 Bit.
- cached = Bool(True)#
Boolean flag, if âTrueâ (default), the result is cached in h5 files.
- h5f = Instance(H5CacheFileBase, transient=True)#
hdf5 cache file
- result = Property()#
The beamforming result as squared sound pressure values at all grid point locations (readonly). Returns a (number of frequencies, number of gridpoints) array-like of floats. Values can only be accessed via the index operator [].
- digest = Property(depends_on=BEAMFORMER_BASE_DIGEST_DEPENDENCIES)#
A unique identifier for the beamformer, based on its properties. (read-only)
- sig_loss_norm()#
If the diagonal of the CSM is removed one has to handle the loss of signal energy.
Done via a normalization factor.
- synthetic(f, num=0)#
Evaluates the beamforming result for an arbitrary frequency band.
- Parameters:
- f: float
Band center frequency.
- numinteger
Controls the width of the frequency bands considered; defaults to 0 (single frequency line).
num
frequency band width
0
single frequency line
1
octave band
3
third-octave band
n
1/n-octave band
- Returns:
- array of floats
The synthesized frequency band values of the beamforming result at each grid point . Note that the frequency resolution and therefore the bandwidth represented by a single frequency line depends on the
sampling frequencyand usedFFT block size.
- integrate(sector, frange=None, num=0)#
Integrates result map over a given sector.
- Parameters:
- sector: array of floats or :class:`~acoular.grids.Sector`
either an array, tuple or list with arguments for the âindicesâ method of a
Grid-derived class (e.g.RectGrid.indicesorRectGrid3D.indices). Possible sectors would be array([xmin, ymin, xmax, ymax]) or array([x, y, radius]) or an instance of aSector-derived class- frange: tuple or None
a tuple of (fmin,fmax) frequencies to include in the result if num*==0, or band center frequency/frequencies for which to return the results if *num>0; if None, then the frequency range is determined from the settings of the
PowerSpectra.ind_lowandPowerSpectra.ind_highoffreq_data- numinteger
Controls the width of the frequency bands considered; defaults to 0 (single frequency line). Only considered if frange is not None.
num
frequency band width
0
single frequency line
1
octave band
3
third-octave band
n
1/n-octave band
- Returns:
- res or (f, res): array of floats or tuple(array of floats, array of floats)
If frange*==None or *num>0, the spectrum (all calculated frequency bands) for the integrated sector is returned as res. The dimension of this array is the number of frequencies given by
freq_dataand entries not computed are zero. If frange!=None and num*==0, then (f, res) is returned where *f are the (band) frequencies and the dimension of both arrays is determined from frange
- class acoular.fbeamform.BeamformerFunctional#
Bases:
BeamformerBaseFunctional beamforming algorithm.
See [4] for details.
- gamma = Float(1)#
Functional exponent, defaults to 1 (= Classic Beamforming).
- r_diag = Enum(False)#
Functional Beamforming is only well defined for full CSM
- r_diag_norm = Enum(1.0)#
Normalization factor in case of CSM diagonal removal. Defaults to 1.0 since Functional Beamforming is only well defined for full CSM.
- digest = Property(depends_on=BEAMFORMER_BASE_DIGEST_DEPENDENCIES + ['gamma'])#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerCapon#
Bases:
BeamformerBaseBeamforming using the Capon (Mininimum Variance) algorithm.
See [5] for details.
- r_diag = Enum(False)#
Boolean flag, if
True, the main diagonal is removed before beamforming; for Capon beamformingr_diagis set toFalse.
- r_diag_norm = Enum(1.0)#
Normalization factor in case of CSM diagonal removal. Defaults to
1.0since Beamformer Capon is only well defined for full CSM.
- class acoular.fbeamform.BeamformerEig#
Bases:
BeamformerBaseBeamforming using eigenvalue and eigenvector techniques.
See [6] for details.
- n = Int(-1)#
Number of component to calculate:
0(smallest) âŠnum_channels-1; defaults to-1, i.e.num_channels-1.
- na = Property()#
No. of eigenvalue
- digest = Property(depends_on=BEAMFORMER_BASE_DIGEST_DEPENDENCIES + ['n'])#
A uniquernal identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerMusic#
Bases:
BeamformerEigBeamforming using the MUSIC algorithm.
See [7] for details.
- r_diag = Enum(False)#
Boolean flag, if
True, the main diagonal is removed before beamforming; for MUSIC beamformingr_diagis set toFalse.
- r_diag_norm = Enum(1.0)#
Normalization factor in case of CSM diagonal removal. Defaults to 1.0 since BeamformerMusic is only well defined for full CSM.
- n = Int(1)#
assumed number of sources
- class acoular.fbeamform.PointSpreadFunction#
Bases:
HasStrictTraitsThe point spread function.
This class provides tools to calculate the PSF depending on the used microphone geometry, focus grid, flow environment, etc. The PSF is needed by several deconvolution algorithms to correct the aberrations when using simple delay-and-sum beamforming.
- steer = Instance(SteeringVector, args=())#
Instance of
SteeringVectoror its derived classes that contains information about the steering vector. This is a private trait. Do not set this directly, usesteertrait instead.
- grid_indices = CArray( âŠ#
Indices of grid points to calculate the PSF for.
- calcmode = Enum('single', 'block', 'full', 'readonly')#
Flag that defines how to calculate and store the point spread function defaults to âsingleâ.
- âfullâ: Calculate the full PSF (for all grid points) in one go (should be used if the PSF
at all grid points is needed, as with
DAMAS)
- âsingleâ: Calculate the PSF for the grid points defined by
grid_indices, one by one (useful if not all PSFs are needed, as with
CLEAN)
- âsingleâ: Calculate the PSF for the grid points defined by
- âblockâ: Calculate the PSF for the grid points defined by
grid_indices, in one go (useful if not all PSFs are needed, as with
CLEAN)
- âblockâ: Calculate the PSF for the grid points defined by
- âreadonlyâ: Do not attempt to calculate the PSF since it should already be cached (useful
if multiple processes have to access the cache file)
- precision = Enum('float64', 'float32')#
Floating point precision of property psf. Corresponding to numpy dtypes. Default = 64 Bit.
- psf = Property()#
The actual point spread function.
- freq = Float(1.0)#
Frequency to evaluate the PSF for; defaults to 1.0.
- digest = Property(depends_on=['steer.digest', 'precision'], cached=True)#
A unique identifier for the object, based on its properties. (read-only)
- calc_psf(ac, gp)#
point-spread function calculation.
- class acoular.fbeamform.BeamformerDamas#
Bases:
BeamformerBaseDAMAS deconvolution algorithm.
See [2] for details.
- psf_precision = Enum('float64', 'float32')#
The floating-number-precision of the PSFs. Default is 64 bit.
- n_iter = Int(100)#
Number of iterations, defaults to 100.
- damp = Float(1.0)#
Damping factor in modified gauss-seidel
- calcmode = Enum('full', 'single', 'block', 'readonly')#
Flag that defines how to calculate and store the point spread function defaults to âsingleâ.
- âfullâ: Calculate the full PSF (for all grid points) in one go (should be used if the PSF
at all grid points is needed, as with
DAMAS)
- âsingleâ: Calculate the PSF for the grid points defined by
grid_indices, one by one (useful if not all PSFs are needed, as with
CLEAN)
- âsingleâ: Calculate the PSF for the grid points defined by
- âblockâ: Calculate the PSF for the grid points defined by
grid_indices, in one go (useful if not all PSFs are needed, as with
CLEAN)
- âblockâ: Calculate the PSF for the grid points defined by
- âreadonlyâ: Do not attempt to calculate the PSF since it should already be cached (useful
if multiple processes have to access the cache file)
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerDamasPlus#
Bases:
BeamformerDamasDAMAS deconvolution [2] for solving the system of equations.
Instead of the original Gauss-Seidel iterations, this class employs the NNLS or linear programming solvers from scipy.optimize or one of several optimization algorithms from the scikit-learn module. Needs a-priori delay-and-sum beamforming (
BeamformerBase).- method = Enum('NNLS', 'LP', 'LassoLars', 'OMPCV')#
Type of fit method to be used (âLassoLarsâ, âOMPCVâ, âLPâ, or âNNLSâ, defaults to âNNLSâ). These methods are implemented in the scikit-learn module or within scipy.optimize respectively.
- n_iter = Int(500)#
Maximum number of iterations, tradeoff between speed and precision; defaults to 500
- unit_mult = Float(1e9)#
Unit multiplier for evaluating, e.g., nPa instead of Pa. Values are converted back before returning. Temporary conversion may be necessary to not reach machine epsilon within fitting method algorithms. Defaults to 1e9.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerOrth#
Bases:
BeamformerBaseOrthogonal deconvolution algorithm.
See [8] for details. New faster implementation without explicit (
BeamformerEig).- eva_list = CArray(dtype=int, value=np.array([-1]))#
List of components to consider, use this to directly set the eigenvalues used in the beamformer. Alternatively, set
n.
- n = Int(1)#
Number of components to consider, defaults to
1. If set,eva_listwill contain the indices of the n largest eigenvalues. Settingeva_listafterwards will override this value.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerCleansc#
Bases:
BeamformerBaseCLEAN-SC deconvolution algorithm.
See [9] for details. Classic delay-and-sum beamforming is already included.
- n_iter = Int(0)#
no of CLEAN-SC iterations defaults to 0, i.e. automatic (max 2*num_channels)
- damp = Range(0.01, 1.0, 0.6)#
iteration damping factor defaults to 0.6
- stopn = Int(3)#
iteration stop criterion for automatic detection iteration stops if power[i]>power[i-stopn] defaults to 3
- digest = Property(depends_on=BEAMFORMER_BASE_DIGEST_DEPENDENCIES + ['n_iter', 'damp', 'stopn'])#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerClean#
Bases:
BeamformerBaseCLEAN deconvolution algorithm.
See [10] for details.
- psf_precision = Enum('float64', 'float32')#
The floating-number-precision of the PSFs. Default is 64 bit.
- damp = Range(0.01, 1.0, 0.6)#
damping factor
- n_iter = Int(100)#
maximum number of iterations
- calcmode = Enum('block', 'full', 'single', 'readonly')#
Flag that defines how to calculate and store the point spread function defaults to âsingleâ.
- âfullâ: Calculate the full PSF (for all grid points) in one go (should be used if the PSF
at all grid points is needed, as with
DAMAS)
- âsingleâ: Calculate the PSF for the grid points defined by
grid_indices, one by one (useful if not all PSFs are needed, as with
CLEAN)
- âsingleâ: Calculate the PSF for the grid points defined by
- âblockâ: Calculate the PSF for the grid points defined by
grid_indices, in one go (useful if not all PSFs are needed, as with
CLEAN)
- âblockâ: Calculate the PSF for the grid points defined by
- âreadonlyâ: Do not attempt to calculate the PSF since it should already be cached (useful
if multiple processes have to access the cache file)
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerCMF#
Bases:
BeamformerBaseCovariance Matrix Fitting algorithm.
This is not really a beamformer, but an inverse method. See [11] for details.
- method = Enum( âŠ#
Type of fit method to be used (âLassoLarsâ, âLassoLarsBICâ, âOMPCVâ or âNNLSâ, defaults to âLassoLarsâ). These methods are implemented in the scikit-learn module.
- alpha = Range(0.0, 1.0, 0.0)#
Weight factor for LassoLars method, defaults to 0.0. (Use values in the order of 10^â»9 for good results.)
- n_iter = Int(500)#
Total or maximum number of iterations (depending on
method), tradeoff between speed and precision; defaults to 500
- unit_mult = Float(1e9)#
Unit multiplier for evaluating, e.g., nPa instead of Pa. Values are converted back before returning. Temporary conversion may be necessary to not reach machine epsilon within fitting method algorithms. Defaults to 1e9.
- show = Bool(False)#
If True, shows the status of the PyLops solver. Only relevant in case of FISTA or Split_Bregman
- r_diag_norm = Enum(None)#
Energy normalization in case of diagonal removal not implemented for inverse methods.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerSODIX#
Bases:
BeamformerBaseSource directivity modeling in the cross-spectral matrix (SODIX) algorithm.
See [12] and [13] for details.
- method = Enum('fmin_l_bfgs_b')#
Type of fit method to be used (âfmin_l_bfgs_bâ). These methods are implemented in the scipy module.
- n_iter = Int(200)#
Maximum number of iterations, tradeoff between speed and precision; defaults to 200
- alpha = Range(0.0, 1.0, 0.0)#
Weight factor for regularization, defaults to 0.0.
- unit_mult = Float(1e9)#
Unit multiplier for evaluating, e.g., nPa instead of Pa. Values are converted back before returning. Temporary conversion may be necessary to not reach machine epsilon within fitting method algorithms. Defaults to 1e9.
- r_diag_norm = Enum(None)#
Energy normalization in case of diagonal removal not implemented for inverse methods.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerGIB#
Bases:
BeamformerEigBeamforming GIB methods with different normalizations.
See [14] for details.
- unit_mult = Float(1e9)#
Unit multiplier for evaluating, e.g., nPa instead of Pa. Values are converted back before returning. Temporary conversion may be necessary to not reach machine epsilon within fitting method algorithms. Defaults to 1e9.
- n_iter = Int(10)#
Total or maximum number of iterations (depending on
method), tradeoff between speed and precision; defaults to 10
- method = Enum( âŠ#
Type of fit method to be used (âSuzukiâ, âLassoLarsâ, âLassoLarsCVâ, âLassoLarsBICâ, âOMPCVâ or âNNLSâ, defaults to âSuzukiâ). These methods are implemented in the scikit-learn module.
- alpha = Range(0.0, 1.0, 0.0)#
Weight factor for LassoLars method, defaults to 0.0.
- pnorm = Float(1)#
Norm to consider for the regularization in InverseIRLS and Suzuki methods defaults to L-1 Norm
- beta = Float(0.9)#
Beta - Fraction of sources maintained after each iteration defaults to 0.9
- eps_perc = Float(0.05)#
eps - Regularization parameter for Suzuki algorithm defaults to 0.05.
- m = Int(0)#
First eigenvalue to consider. Defaults to 0.
- r_diag_norm = Enum(None)#
Energy normalization in case of diagonal removal not implemented for inverse methods.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- class acoular.fbeamform.BeamformerAdaptiveGrid#
Bases:
BeamformerBase,GridAbstract base class for array methods without predefined grid.
- integrate(sector)#
Integrates result map over a given sector.
- Parameters:
- sector: :class:`~acoular.grids.Sector` or derived
Gives the sector over which to integrate
- Returns:
- array of floats
The spectrum (all calculated frequency bands) for the integrated sector.
- class acoular.fbeamform.BeamformerGridlessOrth#
Bases:
BeamformerAdaptiveGridOrthogonal beamforming without predefined grid.
See [15] for details.
- eva_list = CArray(dtype=int, value=np.array([-1]))#
List of components to consider, use this to directly set the eigenvalues used in the beamformer. Alternatively, set
n.
- n = Int(1)#
Number of components to consider, defaults to 1. If set,
eva_listwill contain the indices of the n largest eigenvalues. Settingeva_listafterwards will override this value.
- bounds = List(Tuple(Float, Float), minlen=3, maxlen=3, value=[(-1.0, 1.0), (-1.0, 1.0), (0.01, 1.0)])#
Geometrical bounds of the search domain to consider.
boundis a list that contains exactly three tuple of (min,max) for each of the coordinates x, y, z. Defaults to [(-1.,1.),(-1.,1.),(0.01,1.)]
- shgo = Dict#
options dictionary for the SHGO solver, see scipy docs. Default is Sobol sampling Nelder-Mead local minimizer, 256 initial sampling points and 1 iteration
- r_diag_norm = Enum(1.0)#
If diagonal of the csm is removed, some signal energy is lost. This is handled via this normalization factor. For this class, normalization is not implemented. Defaults to 1.0.
- digest = Property( âŠ#
A unique identifier for the beamformer, based on its properties. (read-only)
- acoular.fbeamform.L_p(x)#
Calculates the sound pressure level from the squared sound pressure.
\(L_p = 10 \lg ( x / 4\cdot 10^{-10})\)
- Parameters:
- x: array of floats
The squared sound pressure values
- Returns:
- array of floats
The corresponding sound pressure levels in dB. If x<0, -350.0 dB is returned.
- acoular.fbeamform.integrate(data, grid, sector)#
Integrates a sound pressure map over a given sector.
This function can be applied on beamforming results to quantitatively analyze the sound pressure in a given sector. If used with
Beamformer.result(), the output is identical to the result of the intrinsicBeamformer.integratemethod. It can, however, also be used with theBeamformer.syntheticoutput.- Parameters:
- data: array of floats
Contains the calculated squared sound pressure values in Pa**2. If data has the same number of entries than the number of grid points only one value is returned. In case of a 2-D array with the second dimension identical to the number of grid points an array containing as many entries as the first dimension is returned.
- grid: Grid object
Object of a
Grid-derived class that provides the grid locations.- sector: array of floats or :class:`~acoular.grids.Sector`-derived object
Tuple with arguments for the indices method of a
Grid-derived class (e.g.RectGrid.indicesorRectGrid3D.indices). Possible sectors would be np.array([xmin, ymin, xmax, ymax]) or np.array([x, y, radius]). Alternatively, aSector-derived object can be used.
- Returns:
- array of floats
The spectrum (all calculated frequency bands) for the integrated sector.