Non-negative least squares

Non-negative least squares minimizes the residual \(b - A x\) under the constraint that \(x\) must have non-negative entries, i.e.,

\[\min_x \| A x - b \|_2 \text{ such that } x \ge 0,\]

and real instances of this problem is well-known to be expressable as the quadratic program

\[\begin{split}& \min_x \frac{1}{2} x^T (A^T A) x + (-A^T b)^T x \\ & \text{s.t. } x \ge 0.\end{split}\]

Elemental defaults to solving this QP via a Mehrotra Predictor-Corrector primal-dual Interior Point Method, but a (prototype) batched version of QP ADMM is also available.

Python API

NNLS(A, B[, ctrl=None])
Parameters
  • A – dense or sparse, sequential or distributed matrix

  • B – dense right-hand side matrix (with type compatible to A)

  • ctrl – (optional) NNLSCtrl_s or NNLSCtrl_d instance, depending upon whether data is single-precision or double-precision

Return type

dense solution vector (with type matching that of b)

C++ API

void NNLS(const Matrix<Real> &A, const Matrix<Real> &B, Matrix<Real> &X, const NNLSCtrl<Real> &ctrl = NNLSCtrl<Real>())
void NNLS(const ElementalMatrix<Real> &A, const ElementalMatrix<Real> &B, ElementalMatrix<Real> &X, const NNLSCtrl<Real> &ctrl = NNLSCtrl<Real>())
void NNLS(const SparseMatrix<Real> &A, const Matrix<Real> &B, Matrix<Real> &X, const NNLSCtrl<Real> &ctrl = NNLSCtrl<Real>())
void NNLS(const DistSparseMatrix<Real> &A, const DistMultiVec<Real> &B, DistMultiVec<Real> &X, const NNLSCtrl<Real> &ctrl = NNLSCtrl<Real>())

C API

Single-precision

ElError ElNNLS_s(ElConstMatrix_s A, ElConstMatrix_s B, ElMatrix_s X)
ElError ElNNLSDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s B, ElDistMatrix_s X)
ElError ElNNLSSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s B, ElMatrix_s X)
ElError ElNNLSDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s B, ElDistMultiVec_s X)

Double-precision

ElError ElNNLS_d(ElConstMatrix_d A, ElConstMatrix_d B, ElMatrix_d X)
ElError ElNNLSDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d B, ElDistMatrix_d X)
ElError ElNNLSSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d B, ElMatrix_d X)
ElError ElNNLSDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d B, ElDistMultiVec_d X)

Expert interface

Single-precision

ElError ElNNLSX_s(ElConstMatrix_s A, ElConstMatrix_s B, ElMatrix_s X, ElNNLSCtrl_s ctrl)
ElError ElNNLSXDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s B, ElDistMatrix_s X, ElNNLSCtrl_s ctrl)
ElError ElNNLSXSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s B, ElMatrix_s X, ElNNLSCtrl_s ctrl)
ElError ElNNLSXDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s B, ElDistMultiVec_s X, ElNNLSCtrl_s ctrl)

Double-precision

ElError ElNNLSX_d(ElConstMatrix_d A, ElConstMatrix_d B, ElMatrix_d X, ElNNLSCtrl_d ctrl)
ElError ElNNLSXDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d B, ElDistMatrix_d X, ElNNLSCtrl_d ctrl)
ElError ElNNLSXSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d B, ElMatrix_d X, ElNNLSCtrl_d ctrl)
ElError ElNNLSXDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d B, ElDistMultiVec_d X, ElNNLSCtrl_d ctrl)