Least Squares¶
A Least Squares problem involves the solution of
and is applicable for overdetermined (as well as square) linear system. Whenever the system is underdetermined, it is more appropriate to solve a Minimum Length problem,
Dense algorithm¶
Elemental solves dense instances of these problems in a straight-forward manner using QR and LQ factorizations, respectively.
Sparse-direct algorithm¶
Sparse instances are solved via applying a priori regularization to the symmetric quasi-semidefinite augmented systems
and
where \(\alpha\) should ideally be chosen near \(\sigma_{\text{min}}(A)\) in order to minimize the condition number (relative to both the augmented system and the solution for \(X\) [Bjorck92] [Bjorck96]). The augmented systems are of interest because they have nearby quasi-definite [Vanderbei95] matrices which be factored with a Cholesky-like sparse-direct solver (and can therefore be effectively preconditioned) [Saunders96].
Lastly, Elemental in fact allows for slight generalizations of the above problems: \(A^T\) or \(A^H\) may also be used in the above equations rather than only \(A\).
C++ API¶
-
void
LeastSquares
(Orientation orientation, const Matrix<F> &A, const Matrix<F> &B, Matrix<F> &X)¶
-
void
LeastSquares
(Orientation orientation, const AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &X)¶
-
void
LeastSquares
(Orientation orientation, const SparseMatrix<F> &A, const Matrix<F> &B, Matrix<F> &X, const LeastSquaresCtrl<Base<F>> &ctrl = LeastSquaresCtrl<Base<F>>())¶
-
void
LeastSquares
(Orientation orientation, const DistSparseMatrix<F> &A, const DistMultiVec<F> &B, DistMultiVec<F> &X, const LeastSquaresCtrl<Base<F>> &ctrl = LeastSquaresCtrl<Base<F>>())¶
Dense versions which overwrite some of the input¶
-
void
ls
::
Overwrite
(Orientation orientation, Matrix<F> &A, const Matrix<F> &B, Matrix<F> &X)¶
-
void
ls
::
Overwrite
(Orientation orientation, AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &X)¶
C API¶
Standard interface¶
Single-precision¶
-
ElError
ElLeastSquares_s
(ElOrientation orientation, ElConstMatrix_s A, ElConstMatrix_s B, ElMatrix_s X)¶
-
ElError
ElLeastSquaresDist_s
(ElOrientation orientation, ElConstDistMatrix_s A, ElConstDistMatrix_s B, ElDistMatrix_s X)¶
-
ElError
ElLeastSquaresSparse_s
(ElOrientation orientation, ElConstSparseMatrix_s A, ElConstMatrix_s B, ElMatrix_s X)¶
-
ElError
ElLeastSquaresDistSparse_s
(ElOrientation orientation, ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s B, ElDistMultiVec_s X)¶
Double-precision¶
-
ElError
ElLeastSquares_d
(ElOrientation orientation, ElConstMatrix_d A, ElConstMatrix_d B, ElMatrix_d X)¶
-
ElError
ElLeastSquaresDist_d
(ElOrientation orientation, ElConstDistMatrix_d A, ElConstDistMatrix_d B, ElDistMatrix_d X)¶
-
ElError
ElLeastSquaresSparse_d
(ElOrientation orientation, ElConstSparseMatrix_d A, ElConstMatrix_d B, ElMatrix_d X)¶
-
ElError
ElLeastSquaresDistSparse_d
(ElOrientation orientation, ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d B, ElDistMultiVec_d X)¶
Single-precision complex¶
-
ElError
ElLeastSquares_c
(ElOrientation orientation, ElConstMatrix_c A, ElConstMatrix_c B, ElMatrix_c X)¶
-
ElError
ElLeastSquaresDist_c
(ElOrientation orientation, ElConstDistMatrix_c A, ElConstDistMatrix_c B, ElDistMatrix_c X)¶
-
ElError
ElLeastSquaresSparse_c
(ElOrientation orientation, ElConstSparseMatrix_c A, ElConstMatrix_c B, ElMatrix_c X)¶
-
ElError
ElLeastSquaresDistSparse_c
(ElOrientation orientation, ElConstDistSparseMatrix_c A, ElConstDistMultiVec_c B, ElDistMultiVec_c X)¶
Double-precision complex¶
-
ElError
ElLeastSquares_z
(ElOrientation orientation, ElConstMatrix_z A, ElConstMatrix_z B, ElMatrix_z X)¶
-
ElError
ElLeastSquaresDist_z
(ElOrientation orientation, ElConstDistMatrix_z A, ElConstDistMatrix_z B, ElDistMatrix_z X)¶
-
ElError
ElLeastSquaresSparse_z
(ElOrientation orientation, ElConstSparseMatrix_z A, ElConstMatrix_z B, ElMatrix_z X)¶
-
ElError
ElLeastSquaresDistSparse_z
(ElOrientation orientation, ElConstDistSparseMatrix_z A, ElConstDistMultiVec_z B, ElDistMultiVec_z X)¶
Expert versions¶
Single-precision¶
-
ElError
ElLeastSquaresXSparse_s
(ElOrientation orientation, ElConstSparseMatrix_s A, ElConstMatrix_s B, ElMatrix_s X, ElLeastSquaresCtrl_s ctrl)¶
-
ElError
ElLeastSquaresXDistSparse_s
(ElOrientation orientation, ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s B, ElDistMultiVec_s X, ElLeastSquaresCtrl_s ctrl)¶
Double-precision¶
-
ElError
ElLeastSquaresXSparse_d
(ElOrientation orientation, ElConstSparseMatrix_d A, ElConstMatrix_d B, ElMatrix_d X, ElLeastSquaresCtrl_d ctrl)¶
-
ElError
ElLeastSquaresXDistSparse_d
(ElOrientation orientation, ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d B, ElDistMultiVec_d X, ElLeastSquaresCtrl_d ctrl)¶
Single-precision complex¶
-
ElError
ElLeastSquaresXSparse_c
(ElOrientation orientation, ElConstSparseMatrix_c A, ElConstMatrix_c B, ElMatrix_c X, ElLeastSquaresCtrl_s ctrl)¶
-
ElError
ElLeastSquaresXDistSparse_c
(ElOrientation orientation, ElConstDistSparseMatrix_c A, ElConstDistMultiVec_c B, ElDistMultiVec_c X, ElLeastSquaresCtrl_s ctrl)¶
Double-precision complex¶
-
ElError
ElLeastSquaresXSparse_z
(ElOrientation orientation, ElConstSparseMatrix_z A, ElConstMatrix_z B, ElMatrix_z X, ElLeastSquaresCtrl_d ctrl)¶
-
ElError
ElLeastSquaresXDistSparse_z
(ElOrientation orientation, ElConstDistSparseMatrix_z A, ElConstDistMultiVec_z B, ElDistMultiVec_z X, ElLeastSquaresCtrl_d ctrl)¶
- Bjorck92
Ake Bjorck, Pivoting and stability in the augmented system method. In D.F. Griffiths and G.A. Watson (eds.), Proc. 14th Dundee Conf., Pitman Research Notes in Math., pp. 1–16, 1992.
- Bjorck96
Ake Bjorck, Numerical methods for least squares problems, SIAM, Philadelphia, 1996. DOI
- Saunders96
Michael Saunders, Chapter 8, Cholesky-based Methods for Sparse Least Squares: The Benefits of Regularization, in L. Adams and J.L. Nazareth (eds.), Linear and Nonlinear Conjugate Gradient-Related Methods, SIAM, Philadelphia, pp. 92–100, 1996. Currently available here
- Vanderbei95
R.J. Vanderbei, Symmetric quasi-definite matrices, SIAM J. Optim., 5(1), pp. 100–113, 1995. Preprint available here