Chebyshev point

A Chebyshev point (CP) minimizes \(b - A x\) in a pointwise sense, i.e.,

\[\min_x \| A x - b \|_{\infty},\]

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

\[\min_{x,t} \{\; t \; | \; -t \le A x - b \le t \; \},\]

which, in affine conic form, becomes

\[\begin{split}\min_{x,t} \{\; \begin{pmatrix} 0 \\ 1 \end{pmatrix}^T \begin{pmatrix} x \\ t \end{pmatrix} \; | \; \begin{pmatrix} A & -1 \\ -A & -1 \end{pmatrix} \begin{pmatrix} x \\ t \end{pmatrix} \le \begin{pmatrix} b \\ -b \end{pmatrix} \; \}.\end{split}\]

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

Python API

CP(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 CP(const Matrix<Real> &A, const Matrix<Real> &b, Matrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void CP(const ElementalMatrix<Real> &A, const ElementalMatrix<Real> &b, ElementalMatrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void CP(const SparseMatrix<Real> &A, const Matrix<Real> &b, Matrix<Real> &x, const lp::affine::Ctrl<Real> &ctrl = lp::affine::Ctrl<Real>())
void CP(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 ElCP_s(ElConstMatrix_s A, ElConstMatrix_s b, ElMatrix_s x)
ElError ElCPDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElDistMatrix_s x)
ElError ElCPSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElMatrix_s x)
ElError ElCPDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElDistMultiVec_s x)

Double-precision

ElError ElCP_d(ElConstMatrix_d A, ElConstMatrix_d b, ElMatrix_d x)
ElError ElCPDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElDistMatrix_d x)
ElError ElCPSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElMatrix_d x)
ElError ElCPDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElDistMultiVec_d x)

Expert interface

Single-precision

ElError ElCPX_s(ElConstMatrix_s A, ElConstMatrix_s b, ElMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElCPXDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElDistMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElCPXSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElMatrix_s x, ElLPAffineCtrl_s ctrl)
ElError ElCPXDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElDistMultiVec_s x, ElLPAffineCtrl_s ctrl)

Double-precision

ElError ElCPX_d(ElConstMatrix_d A, ElConstMatrix_d b, ElMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElCPXDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElDistMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElCPXSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElMatrix_d x, ElLPAffineCtrl_d ctrl)
ElError ElCPXDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElDistMultiVec_d x, ElLPAffineCtrl_d ctrl)