Maximum-likelihood based estimators¶
API references¶
-
class
doatools.estimation.ml.CovarianceBasedMLEstimator(array, wavelength)[source]¶ Bases:
abc.ABCAbstract base class for covariance based maximum-likelihood estimators.
Parameters: - array (ArrayDesign) – Sensor array design.
- wavelength (float) – Wavelength of the carrier wave.
-
get_max_resolvable_sources()[source]¶ Returns the maximum number of sources resolvable.
This default implementation returns (array size - 1), which is suitable for most ML based estimators because the projection matrix of the steering matrix is not well-defined when the number of sources is greater than or equal to the number of sensors.
-
estimate(R, sources0, **kwargs)[source]¶ Solves the ML problem for the given inputs.
Parameters: - R (ndarray) – The sample covariance matrix.
- sources0 (SourcePlacement) –
The initial guess of source locations. Its type determines the source type and its size determines the number of sources.
Because the log-likelihood function is highly non-convex, the initial guess of source locations will greatly affect the final estimates. It is recommended to use the output of another estimator (e.g. conventional beamformer) as the initial guess.
- **kwargs – Additional keyword arguments for the solver.
Notes
In general, ML estimates are computationally expensive to obtain and sensitive to initialization. They are generally used in theoretical performance analyses.
Returns: A tuple containing the following elements: - resolved (
bool):Trueif the optimizer exited successfully. This flag does not guarantee that the estimated source locations are correct. The estimated source locations may be completely wrong! If resolved is False,estimateswill beNone. - estimates (
SourcePlacement): ASourcePlacementinstance of the same type as that ofsources0, represeting the estimated source locations. Will beNoneif resolved isFalse.
Return type: tuple
-
class
doatools.estimation.ml.AMLEstimator(array, wavelength)[source]¶ Bases:
doatools.estimation.ml.CovarianceBasedMLEstimatorAsymptotic maximum-likelihood (AML) estimator.
The AML estimator maximizes the following log-likelihood function:
\[- \log\det \mathbf{S} - \mathrm{tr}(\mathbf{S}^{-1} \mathbf{R})\]where \(\mathbf{S} = \mathbf{A}\mathbf{P}\mathbf{A}^H + \sigma^2\mathbf{I}\), \(\mathbf{A}\) is the steering matrix, \(\mathbf{P}\) is the source covariance matrix, \(\sigma^2\) is the noise variance, and \(\hat{\mathbf{R}}\) is the sample covariance matrix.
Here the unknown parameters include the source locations, \(\mathbf{P}\), and \(\sigma^2\). The MLE of \(\mathbf{P}\) and \(\sigma^2\) can be analytically obtained in terms of the source locations. The final optimization problem only involves the source locations, \(\mathbf{\theta}\), as unknown variables:
\[\min_{\mathbf{\theta}} \log\det\left\lbrack \mathbf{P}_{\mathbf{A}} \hat{\mathbf{R}} \mathbf{P}_{\mathbf{A}} + \frac{ \mathrm{tr}(\mathbf{P}^\perp_{\mathbf{A}} \hat{\mathbf{R}}) \mathbf{P}^\perp_{\mathbf{A}} }{N-D} \right\rbrack.\]References
[1] H. L. Van Trees, Optimum array processing. New York: Wiley, 2002.
-
class
doatools.estimation.ml.CMLEstimator(array, wavelength)[source]¶ Bases:
doatools.estimation.ml.CovarianceBasedMLEstimatorConditional maximum-likelihood (CML) estimator.
Given the conditional observation model (the source signals are assumed to be deterministic unknown):
\[\mathbf{y}(t) = \mathbf{A}(\mathbf{\theta})\mathbf{x}(t) + \mathbf{n}(t), t = 1,2,...,T,\]the CML estimator maximizes the following log-likelihood function:
\[- TM\log\sigma^2 - \sigma^{-2} \sum_{t=1}^T \| \mathbf{y}(t) - \mathbf{A}\mathbf{x}(t) \|^2,\]where \(M\) is the number of sensors, \(T\) is the number of snapshots, \(\mathbf{A}\) is the steering matrix, \(\sigma^2\) is the noise variance.
Here the unknown parameters include the source locations, \(\mathbf{\theta}\), as well as \(\mathbf{x}(t)\) and \(\sigma^2\). With further computations, it can be shown that the final optimization problem only involves the source locations:
\[\mathrm{tr}(\mathbf{P}^\perp_{\mathbf{A}} \hat{\mathbf{R}}),\]where \(\hat{\mathbf{R}} = 1/T \sum_{t=1}^T \mathbf{x}(t)\mathbf{x}^H(t)\).
References
[1] H. L. Van Trees, Optimum array processing. New York: Wiley, 2002.
-
class
doatools.estimation.ml.WSFEstimator(array, wavelength)[source]¶ Bases:
doatools.estimation.ml.CovarianceBasedMLEstimatorWeighted subspace fitting (WSF) estimator.
WSF is based on the CML estimator, with the objective function given by
\[\mathrm{tr}(\mathbf{P}^\perp_{\mathbf{A}} \hat{\mathbf{U}}_\mathrm{s} \hat{\mathbf{W}} \hat{\mathbf{U}}_\mathrm{s}^H),\]where \(\hat{\mathbf{U}}_\mathrm{s}\) consists of the eigenvectors of the signal subspace of \(\hat{\mathbf{R}}\), and \(\hat{\mathbf{W}}\) is a diagonal matrix consists of asymptotically optimal weights.
References
[1] M. Viberg and B. Ottersten, “Sensor array processing based on subspace fitting,” IEEE Transactions on Signal Processing, vol. 39, no. 5, pp. 1110-1121, May 1991.
[2] H. L. Van Trees, Optimum array processing. New York: Wiley, 2002.
[3] P. Stoica and K. Sharman, “Maximum likelihood methods for direction-of-arrival estimation,” IEEE Trans. Acoust., Speech, Signal Process., vol. 38, pp. 1132-1143, July 1990.