Direct Linear Programs

\[\begin{split}\min_x & \{\; c^T x \; | \; A x = b \;\wedge\; x \ge 0 \;\}, \\ \max_{y,z} & \{\; - b^T y \; | \; A^T y - z + c = 0 \;\wedge\; z \ge 0 \;\}\end{split}\]

Python API

LPDirect(A, b, c, x, y, z[, ctrl=None])

C++ API

void LP(const Matrix<Real> &A, const Matrix<Real> &b, const Matrix<Real> &c, Matrix<Real> &x, Matrix<Real> &y, Matrix<Real> &z, const lp::direct::Ctrl<Real> &ctrl = lp::direct::Ctrl<Real>(false))
void LP(const ElementalMatrix<Real> &A, const ElementalMatrix<Real> &b, const ElementalMatrix<Real> &c, ElementalMatrix<Real> &x, ElementalMatrix<Real> &y, ElementalMatrix<Real> &z, const lp::direct::Ctrl<Real> &ctrl = lp::direct::Ctrl<Real>(false))
void LP(const SparseMatrix<Real> &A, const Matrix<Real> &b, const Matrix<Real> &c, Matrix<Real> &x, Matrix<Real> &y, Matrix<Real> &z, const lp::direct::Ctrl<Real> &ctrl = lp::direct::Ctrl<Real>(true))
void LP(const DistSparseMatrix<Real> &A, const DistMultiVec<Real> &b, const DistMultiVec<Real> &c, DistMultiVec<Real> &x, DistMultiVec<Real> &y, DistMultiVec<Real> &z, const lp::direct::Ctrl<Real> &ctrl = lp::direct::Ctrl<Real>(true))

C API

Single-precision

ElError ElLPDirect_s(ElConstMatrix_s A, ElConstMatrix_s b, ElConstMatrix_s c, ElMatrix_s x, ElMatrix_s y, ElMatrix_s z)
ElError ElLPDirectDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElConstDistMatrix_s c, ElDistMatrix_s x, ElDistMatrix_s y, ElDistMatrix_s z)
ElError ElLPDirectSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElConstMatrix_s c, ElMatrix_s x, ElMatrix_s y, ElMatrix_s z)
ElError ElLPDirectDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElConstDistMultiVec_s c, ElDistMultiVec_s x, ElDistMultiVec_s y, ElDistMultiVec_s z)

Double-precision

ElError ElLPDirect_d(ElConstMatrix_d A, ElConstMatrix_d b, ElConstMatrix_d c, ElMatrix_d x, ElMatrix_d y, ElMatrix_d z)
ElError ElLPDirectDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElConstDistMatrix_d c, ElDistMatrix_d x, ElDistMatrix_d y, ElDistMatrix_d z)
ElError ElLPDirectSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElConstMatrix_d c, ElMatrix_d x, ElMatrix_d y, ElMatrix_d z)
ElError ElLPDirectDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElConstDistMultiVec_d c, ElDistMultiVec_d x, ElDistMultiVec_d y, ElDistMultiVec_d z)

Expert interface

Single-precision

ElError ElLPDirectX_s(ElConstMatrix_s A, ElConstMatrix_s b, ElConstMatrix_s c, ElMatrix_s x, ElMatrix_s y, ElMatrix_s z, ElLPDirectCtrl_s ctrl)
ElError ElLPDirectXDist_s(ElConstDistMatrix_s A, ElConstDistMatrix_s b, ElConstDistMatrix_s c, ElDistMatrix_s x, ElDistMatrix_s y, ElDistMatrix_s z, ElLPDirectCtrl_s ctrl)
ElError ElLPDirectXSparse_s(ElConstSparseMatrix_s A, ElConstMatrix_s b, ElConstMatrix_s c, ElMatrix_s x, ElMatrix_s y, ElMatrix_s z, ElLPDirectCtrl_s ctrl)
ElError ElLPDirectXDistSparse_s(ElConstDistSparseMatrix_s A, ElConstDistMultiVec_s b, ElConstDistMultiVec_s c, ElDistMultiVec_s x, ElDistMultiVec_s y, ElDistMultiVec_s z, ElLPDirectCtrl_s ctrl)

Double-precision

ElError ElLPDirectX_d(ElConstMatrix_d A, ElConstMatrix_d b, ElConstMatrix_d c, ElMatrix_d x, ElMatrix_d y, ElMatrix_d z, ElLPDirectCtrl_d ctrl)
ElError ElLPDirectXDist_d(ElConstDistMatrix_d A, ElConstDistMatrix_d b, ElConstDistMatrix_d c, ElDistMatrix_d x, ElDistMatrix_d y, ElDistMatrix_d z, ElLPDirectCtrl_d ctrl)
ElError ElLPDirectXSparse_d(ElConstSparseMatrix_d A, ElConstMatrix_d b, ElConstMatrix_d c, ElMatrix_d x, ElMatrix_d y, ElMatrix_d z, ElLPDirectCtrl_d ctrl)
ElError ElLPDirectXDistSparse_d(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, ElConstDistMultiVec_d c, ElDistMultiVec_d x, ElDistMultiVec_d y, ElDistMultiVec_d z, ElLPDirectCtrl_d ctrl)