Least Absolute Value regression

Least Absolute Value (LAV) regression minimizes a residual \(b - A x\) in the one-norm, i.e.,

\[\min_x \| A x - b \|_1,\]

and is known to be the solution of a linear program when \(A\) and \(b\) are real. In particular, we may solve

\[\begin{split}\min_{x,u,v} \{\; 1^T \begin{pmatrix} u \\ v \end{pmatrix} \; | \; \begin{pmatrix} A & I & -I \end{pmatrix} \begin{pmatrix} x \\ u \\ v \end{pmatrix} = b \; \wedge \; u,v \ge 0 \; \}.\end{split}\]

By default, Elemental solves this linear program via a Mehrotra Predictor-Corrector primal-dual Interior Point Method.

Python API

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

  • b – dense right-hand side vector (with type compatible to A)

  • ctrl – (optional) LPAffineCtrl instance

Return type

dense solution vector (with type matching that of b)

C++ API

void LAV(const Matrix<Real> &A, const Matrix<Real> &b, Matrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void LAV(const ElementalMatrix<Real> &A, const ElementalMatrix<Real> &b, ElementalMatrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void LAV(const SparseMatrix<Real> &A, const Matrix<Real> &b, Matrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void LAV(const DistSparseMatrix<Real> &A, const DistMultiVec<Real> &b, DistMultiVec<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())

C API

Single-precision

ElError ElLAV_s(ElConstMatrix_s A, ElConstMatrix_s b, ElMatrix_s x)
ElError ElLAVDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElDistMatrix_s x)
ElError ElLAVSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElMatrix_s x)
ElError ElLAVDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElDistMultiVec_s x)

Double-precision

ElError ElLAV_d(ElConstMatrix_d A, ElConstMatrix_d b, ElMatrix_d x)
ElError ElLAVDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElDistMatrix_d x)
ElError ElLAVSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElMatrix_d x)
ElError ElLAVDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElDistMultiVec_d x)

Expert interface

Single-precision

ElError ElLAVX_s(ElConstMatrix_s A, ElConstMatrix_s b, ElMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElLAVXDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElDistMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElLAVXSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElLAVXDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElDistMultiVec_s x, ElLPAffineCtrl_s ctrl)

Double-precision

ElError ElLAVX_d(ElConstMatrix_d A, ElConstMatrix_d b, ElMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElLAVXDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElDistMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElLAVXSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElLAVXDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElDistMultiVec_d x, ElLPAffineCtrl_d ctrl)