Least Squares¶
Given \(A \in \mathbb{F}^{m \times n}\) and a right-hand side \(b \in \mathbb{F}^m\), a least-squares method solves \(A x \approx b\) differently depending upon whether \(m \ge n\).
When \(m \ge n\), there are at least as many constraints as degrees of freedom, and so a solution is sought for
This problem is solved through the use of QR()
.
When \(m < n\), the problem is under-constrained and a solution is sought for the problem
This problem is solved through the use of LQ()
.
The above optimization problems can be readily generalized to multiple right-hand sides by switching to Frobenius norms.
Note
If orientation is set to NORMAL
, then solve \(AX=B\), otherwise
orientation must be equal to ADJOINT
and \(A^H X=B\) will
be solved. Upon completion, \(A\) is overwritten with its QR or LQ
factorization, and each column of \(X\) is overwritten with a solution
vector.
C++ API¶
-
void
LeastSquares
(Orientation orientation, Matrix<F> &A, const Matrix<F> &B, Matrix<F> &X)¶
-
void
LeastSquares
(Orientation orientation, AbstractDistMatrix<F> &A, const AbstractDistMatrix<F> &B, AbstractDistMatrix<F> &X)¶
C API¶
-
ElError
ElLeastSquares_s
(ElOrientation orientation, ElMatrix_s A, ElConstMatrix_s B, ElMatrix_s X)¶
-
ElError
ElLeastSquares_d
(ElOrientation orientation, ElMatrix_d A, ElConstMatrix_d B, ElMatrix_d X)¶
-
ElError
ElLeastSquares_c
(ElOrientation orientation, ElMatrix_c A, ElConstMatrix_c B, ElMatrix_c X)¶
-
ElError
ElLeastSquares_z
(ElOrientation orientation, ElMatrix_z A, ElConstMatrix_z B, ElMatrix_z X)¶
-
ElError
ElLeastSquaresDist_s
(ElOrientation orientation, ElDistMatrix_s A, ElConstDistMatrix_s B, ElDistMatrix_s X)¶
-
ElError
ElLeastSquaresDist_d
(ElOrientation orientation, ElDistMatrix_d A, ElConstDistMatrix_d B, ElDistMatrix_d X)¶
-
ElError
ElLeastSquaresDist_c
(ElOrientation orientation, ElDistMatrix_c A, ElConstDistMatrix_c B, ElDistMatrix_c X)¶
-
ElError
ElLeastSquaresDist_z
(ElOrientation orientation, ElDistMatrix_z A, ElConstDistMatrix_z B, ElDistMatrix_z X)¶