Supported solvers¶
Solvers that are detected as installed on your machine are listed in:
- qpsolvers.available_solvers = ['clarabel', 'cvxopt', 'daqp', 'ecos', 'gurobi', 'highs', 'hpipm', 'mosek', 'osqp', 'piqp', 'proxqp', 'qpalm', 'qpoases', 'qpswift', 'quadprog', 'scs', 'qpax', 'nppro']¶
Built-in mutable sequence.
If no argument is given, the constructor creates a new empty list. The argument must be an iterable if specified.
Clarabel¶
Solver interface for Clarabel.rs.
Clarabel.rs is a Rust implementation of an interior point numerical solver for convex optimization problems using a novel homogeneous embedding. A paper describing the Clarabel solver algorithm and implementation will be forthcoming soon (retrieved: 2023-02-06). Until then, the authors ask that you cite its documentation if you have found Clarabel.rs useful in your work.
- qpsolvers.solvers.clarabel_.clarabel_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using Clarabel.rs.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded as options to Clarabel.rs. For instance, we can call
clarabel_solve_qp(P, q, G, h, u, tol_feas=1e-6)
. Clarabel options include the following:Name
Description
max_iter
Maximum number of iterations.
time_limit
Time limit for solve run in seconds (can be fractional).
tol_gap_abs
absolute duality-gap tolerance
tol_gap_rel
relative duality-gap tolerance
tol_feas
feasibility check tolerance (primal and dual)
Check out the API reference for details.
Lower values for absolute or relative tolerances yield more precise solutions at the cost of computation time. See e.g. [tolerances] for an overview of solver tolerances.
- qpsolvers.solvers.clarabel_.clarabel_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using Clarabel.rs.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using Clarabel.rs.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality matrix.h (
Optional
[ndarray
]) – Linear inequality vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.
CVXOPT¶
Solver interface for CVXOPT.
CVXOPT is a free software package for convex optimization in Python. Its main purpose is to make the development of software for convex optimization applications straightforward by building on Python’s extensive standard library and on the strengths of Python as a high-level programming language. If you are using CVXOPT in some academic work, consider citing the corresponding report [Vandenberghe2010].
- qpsolvers.solvers.cvxopt_.cvxopt_solve_problem(problem, solver=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using CVXOPT.
- Parameters:
problem (
Problem
) – Quadratic program to solve.solver (
Optional
[str
]) – Set to ‘mosek’ to run MOSEK rather than CVXOPT.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError – If the CVXOPT rank assumption is not satisfied.
SolverError – If CVXOPT failed with an error.
Note
Rank assumptions: CVXOPT requires the QP matrices to satisfy the
\[\begin{split}\begin{array}{cc} \mathrm{rank}(A) = p & \mathrm{rank}([P\ A^T\ G^T]) = n \end{array}\end{split}\]where \(p\) is the number of rows of \(A\) and \(n\) is the number of optimization variables. See the “Rank assumptions” paragraph in the report The CVXOPT linear and quadratic cone program solvers for details.
Notes
CVXOPT only considers the lower entries of \(P\), therefore it will use a different cost than the one intended if a non-symmetric matrix is provided.
Keyword arguments are forwarded as options to CVXOPT. For instance, we can call
cvxopt_solve_qp(P, q, G, h, u, abstol=1e-4, reltol=1e-4)
. CVXOPT options include the following:Name
Description
abstol
Absolute tolerance on the duality gap.
feastol
Tolerance on feasibility conditions, that is, on the primal residual.
maxiters
Maximum number of iterations.
refinement
Number of iterative refinement steps when solving KKT equations
reltol
Relative tolerance on the duality gap.
Check out Algorithm Parameters section of the solver documentation for details and default values of all solver parameters. See also [tolerances] for a primer on the duality gap, primal and dual residuals.
- qpsolvers.solvers.cvxopt_.cvxopt_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, solver=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using CVXOPT.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using CVXOPT.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Symmetric cost matrix. Together with \(A\) and \(G\), it should satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\), see the rank assumptions below.q (
ndarray
) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality matrix. Together with \(P\) and \(A\), it should satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\), see the rank assumptions below.h (
Optional
[ndarray
]) – Linear inequality vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix. It needs to be full row rank, and together with \(P\) and \(G\) satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\). See the rank assumptions below.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.solver (
Optional
[str
]) – Set to ‘mosek’ to run MOSEK rather than CVXOPT.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.- Raises:
ProblemError – If the CVXOPT rank assumption is not satisfied.
SolverError – If CVXOPT failed with an error.
DAQP¶
Solver interface for DAQP.
DAQP is a dual active-set algorithm implemented in C [Arnstrom2022]. It has been developed to solve small/medium scale dense problems.
- qpsolvers.solvers.daqp_.daqp_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using DAQP.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded to DAQP. For instance, we can call
daqp_solve_qp(P, q, G, h, u, primal_tol=1e-6, iter_limit=1000)
. DAQP settings include the following:Name
Description
iter_limit
Maximum number of iterations.
primal_tol
Primal feasibility tolerance.
dual_tol
Dual feasibility tolerance.
Check out the DAQP settings documentation for all available settings.
- qpsolvers.solvers.daqp_.daqp_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using DAQP.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using DAQP.
- Parameters:
P (
ndarray
) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded to DAQP. For instance, we can call
daqp_solve_qp(P, q, G, h, u, primal_tol=1e-6, iter_limit=1000)
. DAQP settings include the following:Name
Description
iter_limit
Maximum number of iterations.
primal_tol
Primal feasibility tolerance.
dual_tol
Dual feasibility tolerance.
Check out the DAQP settings documentation for all available settings.
ECOS¶
Solver interface for ECOS.
ECOS is an interior-point solver for convex second-order cone programs (SOCPs). designed specifically for embedded applications. ECOS is written in low footprint, single-threaded, library-free ANSI-C and so runs on most embedded platforms. For small problems, ECOS is faster than most existing SOCP solvers; it is still competitive for medium-sized problems up to tens of thousands of variables. If you are using ECOS in some academic work, consider citing the corresponding paper [Domahidi2013].
- qpsolvers.solvers.ecos_.ecos_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using ECOS.
- Parameters:
P – Primal quadratic cost matrix.
q – Primal quadratic cost vector.
G – Linear inequality constraint matrix.
h – Linear inequality constraint vector.
A – Linear equality constraint matrix.
b – Linear equality constraint vector.
initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError : – If inequality constraints contain infinite values that the solver doesn’t handle.
ValueError : – If the cost matrix is not positive definite.
Notes
All other keyword arguments are forwarded as options to the ECOS solver. For instance, you can call
qpswift_solve_qp(P, q, G, h, abstol=1e-5)
.For a quick overview, the solver accepts the following settings:
Name
Effect
feastol
Tolerance on the primal and dual residual.
abstol
Absolute tolerance on the duality gap.
reltol
Relative tolerance on the duality gap.
feastol_inacc
Tolerance on the primal and dual residual if reduced precisions.
abstol_inacc
Absolute tolerance on the duality gap if reduced precision.
reltolL_inacc
Relative tolerance on the duality gap if reduced precision.
max_iters
Maximum numer of iterations.
nitref
Number of iterative refinement steps.
See the ECOS Python wrapper documentation for more details. You can also check out [tolerances] for a primer on primal-dual residuals or the duality gap.
- qpsolvers.solvers.ecos_.ecos_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using ECOS.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \end{array}\end{split}\end{split}\]It is solved using ECOS.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Primal quadratic cost matrix.q (
ndarray
) – Primal quadratic cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
All other keyword arguments are forwarded as options to the ECOS solver. For instance, you can call
qpswift_solve_qp(P, q, G, h, abstol=1e-5)
.For a quick overview, the solver accepts the following settings:
Name
Effect
feastol
Tolerance on the primal and dual residual.
abstol
Absolute tolerance on the duality gap.
reltol
Relative tolerance on the duality gap.
feastol_inacc
Tolerance on the primal and dual residual if reduced precisions.
abstol_inacc
Absolute tolerance on the duality gap if reduced precision.
reltolL_inacc
Relative tolerance on the duality gap if reduced precision.
max_iters
Maximum numer of iterations.
nitref
Number of iterative refinement steps.
See the ECOS Python wrapper documentation for more details. You can also check out [tolerances] for a primer on primal-dual residuals or the duality gap.
Gurobi¶
Solver interface for Gurobi.
The Gurobi Optimizer suite ships several solvers for mathematical programming, including problems that have linear constraints, bound constraints, integrality constraints, cone constraints, or quadratic constraints. It targets modern CPU architectures and multi-core processors,
See the installation page for additional instructions on installing this solver.
- qpsolvers.solvers.gurobi_.gurobi_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using Gurobi.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
Notes
Keyword arguments are forwarded to Gurobi as parameters. For instance, we can call
gurobi_solve_qp(P, q, G, h, u, FeasibilityTol=1e-8, OptimalityTol=1e-8)
. Gurobi settings include the following:Name
Description
FeasibilityTol
Primal feasibility tolerance.
OptimalityTol
Dual feasibility tolerance.
PSDTol
Positive semi-definite tolerance.
TimeLimit
Run time limit in seconds, 0 to disable.
Check out the Parameter Descriptions documentation for all available Gurobi parameters.
Lower values for primal or dual tolerances yield more precise solutions at the cost of computation time. See e.g. [tolerances] for a primer of solver tolerances.
- qpsolvers.solvers.gurobi_.gurobi_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using Gurobi.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using Gurobi.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Primal quadratic cost matrix.q (
ndarray
) – Primal quadratic cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded to Gurobi as parameters. For instance, we can call
gurobi_solve_qp(P, q, G, h, u, FeasibilityTol=1e-8, OptimalityTol=1e-8)
. Gurobi settings include the following:Name
Description
FeasibilityTol
Primal feasibility tolerance.
OptimalityTol
Dual feasibility tolerance.
PSDTol
Positive semi-definite tolerance.
TimeLimit
Run time limit in seconds, 0 to disable.
Check out the Parameter Descriptions documentation for all available Gurobi parameters.
Lower values for primal or dual tolerances yield more precise solutions at the cost of computation time. See e.g. [tolerances] for a primer of solver tolerances.
HiGHS¶
Solver interface for HiGHS.
HiGHS is an open source, serial and parallel solver for large scale sparse linear programming (LP), mixed-integer programming (MIP), and quadratic programming (QP). It is written mostly in C++11 and available under the MIT licence. HiGHS’s QP solver implements a Nullspace Active Set method. It works best on moderately-sized dense problems. If you are using HiGHS in some academic work, consider citing the corresponding paper [Huangfu2018].
- qpsolvers.solvers.highs_.highs_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using HiGHS.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector for the primal solution.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
Notes
Keyword arguments are forwarded to HiGHS as options. For instance, we can call
highs_solve_qp(P, q, G, h, u, primal_feasibility_tolerance=1e-8, dual_feasibility_tolerance=1e-8)
. HiGHS settings include the following:Name
Description
dual_feasibility_tolerance
Dual feasibility tolerance.
primal_feasibility_tolerance
Primal feasibility tolerance.
time_limit
Run time limit in seconds.
Check out the HiGHS documentation for more information on the solver.
- qpsolvers.solvers.highs_.highs_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using HiGHS.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using HiGHS.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Positive semidefinite cost matrix.q (
ndarray
) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector for the primal solution.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded to HiGHS as options. For instance, we can call
highs_solve_qp(P, q, G, h, u, primal_feasibility_tolerance=1e-8, dual_feasibility_tolerance=1e-8)
. HiGHS settings include the following:Name
Description
dual_feasibility_tolerance
Dual feasibility tolerance.
primal_feasibility_tolerance
Primal feasibility tolerance.
time_limit
Run time limit in seconds.
Check out the HiGHS documentation for more information on the solver.
HPIPM¶
Solver interface for HPIPM.
HPIPM is a high-performance interior point method for solving convex quadratic programs. It is designed to be efficient for small to medium-size problems arising in model predictive control and embedded optmization. If you are using HPIPM in some academic work, consider citing the corresponding paper [Frison2020].
- qpsolvers.solvers.hpipm_.hpipm_solve_problem(problem, initvals=None, mode='balance', verbose=False, **kwargs)¶
Solve a quadratic program using HPIPM.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.mode (
str
) – Solver mode, which provides a set of default solver arguments. Pick one of [“speed_abs”, “speed”, “balance”, “robust”]. These modes are documented in section 4.2 IPM implementation choices of the reference paper [Frison2020]. The default one is “balance”.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
Notes
Keyword arguments are forwarded to HPIPM. For instance, we can call
hpipm_solve_qp(P, q, G, h, u, tol_eq=1e-5)
. HPIPM settings include the following:Name
Description
iter_max
Maximum number of iterations.
tol_eq
Equality constraint tolerance.
tol_ineq
Inequality constraint tolerance.
tol_comp
Complementarity condition tolerance.
tol_stat
Stationarity condition tolerance.
- qpsolvers.solvers.hpipm_.hpipm_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, mode='balance', verbose=False, **kwargs)¶
Solve a quadratic program using HPIPM.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using HPIPM.
- Parameters:
P (
ndarray
) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.mode (
str
) – Solver mode, which provides a set of default solver arguments. Pick one of [“speed_abs”, “speed”, “balance”, “robust”]. Default is “balance”.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
Notes
Keyword arguments are forwarded to HPIPM. For instance, we can call
hpipm_solve_qp(P, q, G, h, u, tol_eq=1e-5)
. HPIPM settings include the following:Name
Description
iter_max
Maximum number of iterations.
tol_eq
Equality constraint tolerance.
tol_ineq
Inequality constraint tolerance.
tol_comp
Complementarity condition tolerance.
tol_stat
Stationarity condition tolerance.
MOSEK¶
Solver interface for MOSEK.
MOSEK is a solver for linear, mixed-integer linear, quadratic, mixed-integer quadratic, quadratically constraint, conic and convex nonlinear mathematical optimization problems. Its interior-point method is geared towards large scale sparse problems, in particular for linear or conic quadratic programs.
- qpsolvers.solvers.mosek_.mosek_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using MOSEK.
- Parameters:
P – Symmetric cost matrix.
q – Cost vector.
G – Linear inequality constraint matrix.
h – Linear inequality constraint vector.
A – Linear equality constraint matrix.
b – Linear equality constraint vector.
lb – Lower bound constraint vector.
ub – Upper bound constraint vector.
initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.
- qpsolvers.solvers.mosek_.mosek_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using MOSEK.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using the MOSEK interface from CVXOPT.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.
OSQP¶
Solver interface for OSQP.
The OSQP solver implements an Operator-Splitting method for large-scale convex quadratic programming. It is designed for both dense and sparse problems, and convexity is the only assumption it makes on problem data (for instance, it does not make any rank assumption contrary to CVXOPT or qpSWIFT). If you are using OSQP in some academic work, consider citing the corresponding paper [Stellato2020].
- qpsolvers.solvers.osqp_.osqp_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using OSQP.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
- Raises:
ValueError – If the problem is clearly non-convex. See this recommendation. Note that OSQP may find the problem unfeasible if the problem is slightly non-convex (in this context, the meaning of “clearly” and “slightly” depends on how close the negative eigenvalues of \(P\) are to zero).
Note
OSQP requires a symmetric P and won’t check for errors otherwise. Check out this point if you get nan values in your solutions.
Notes
Keyword arguments are forwarded to OSQP. For instance, we can call
osqp_solve_qp(P, q, G, h, u, eps_abs=1e-8, eps_rel=0.0)
. OSQP settings include the following:Name
Description
max_iter
Maximum number of iterations.
time_limit
Run time limit in seconds, 0 to disable.
eps_abs
Absolute feasibility tolerance. See Convergence.
eps_rel
Relative feasibility tolerance. See Convergence.
eps_prim_inf
Primal infeasibility tolerance.
eps_dual_inf
Dual infeasibility tolerance.
polish
Perform polishing. See Polishing.
Check out the OSQP settings documentation for all available settings.
Lower values for absolute or relative tolerances yield more precise solutions at the cost of computation time. See e.g. [tolerances] for an overview of solver tolerances.
- qpsolvers.solvers.osqp_.osqp_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using OSQP.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using OSQP.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ValueError –
If the problem is clearly non-convex. See this recommendation. Note that OSQP may find the problem unfeasible if the problem is slightly non-convex (in this context, the meaning of “clearly” and “slightly” depends on how close the negative eigenvalues of \(P\) are to zero).
Note
OSQP requires a symmetric P and won’t check for errors otherwise. Check out this point if you get nan values in your solutions.
Notes
Keyword arguments are forwarded to OSQP. For instance, we can call
osqp_solve_qp(P, q, G, h, u, eps_abs=1e-8, eps_rel=0.0)
. OSQP settings include the following:Name
Description
max_iter
Maximum number of iterations.
time_limit
Run time limit in seconds, 0 to disable.
eps_abs
Absolute feasibility tolerance. See Convergence.
eps_rel
Relative feasibility tolerance. See Convergence.
eps_prim_inf
Primal infeasibility tolerance.
eps_dual_inf
Dual infeasibility tolerance.
polish
Perform polishing. See Polishing.
Check out the OSQP settings documentation for all available settings.
Lower values for absolute or relative tolerances yield more precise solutions at the cost of computation time. See e.g. [tolerances] for an overview of solver tolerances.
PIQP¶
Solver interface for PIQP.
PIQP is a Proximal Interior Point Quadratic Programming solver, which can solve dense and sparse quadratic programs. Combining an infeasible interior point method with the proximal method of multipliers, the algorithm can handle ill-conditioned convex QP problems without the need for linear independence of the constraints. [schwan2023piqp]
- qpsolvers.solvers.piqp_.piqp_solve_problem(problem, initvals=None, verbose=False, backend=None, **kwargs)¶
Solve a quadratic program using PIQP.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).backend (
Optional
[str
]) – PIQP backend to use in[None, "dense", "sparse"]
. IfNone
(default), the backend is selected based on the type ofP
.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP returned by the solver.
Notes
All other keyword arguments are forwarded as options to PIQP. For instance, you can call
piqp_solve_qp(P, q, G, h, eps_abs=1e-6)
. For a quick overview, the solver accepts the following settings:Name
Effect
rho_init
Initial value for the primal proximal penalty parameter rho.
delta_init
Initial value for the augmented lagrangian penalty parameter delta.
eps_abs
Absolute tolerance.
eps_rel
Relative tolerance.
check_duality_gap
Check terminal criterion on duality gap.
eps_duality_gap_abs
Absolute tolerance on duality gap.
eps_duality_gap_rel
Relative tolerance on duality gap.
reg_lower_limit
Lower limit for regularization.
reg_finetune_lower_limit
Fine tune lower limit regularization.
reg_finetune_primal_update_threshold
Threshold of number of no primal updates to transition to fine tune mode.
reg_finetune_dual_update_threshold
Threshold of number of no dual updates to transition to fine tune mode.
max_iter
Maximum number of iterations.
max_factor_retires
Maximum number of factorization retires before failure.
preconditioner_scale_cost
Scale cost in Ruiz preconditioner.
preconditioner_iter
Maximum of preconditioner iterations.
tau
Maximum interior point step length.
iterative_refinement_always_enabled
Always run iterative refinement and not only on factorization failure.
iterative_refinement_eps_abs
Iterative refinement absolute tolerance.
iterative_refinement_eps_rel
Iterative refinement relative tolerance.
iterative_refinement_max_iter
Maximum number of iterations for iterative refinement.
iterative_refinement_min_improvement_rate
Minimum improvement rate for iterative refinement.
iterative_refinement_static_regularization_eps
Static regularization for KKT system for iterative refinement.
iterative_refinement_static_regularization_rel
Static regularization w.r.t. the maximum abs diagonal term of KKT system.
verbose
Verbose printing.
compute_timings
Measure timing information internally.
This list is not exhaustive. Check out the solver documentation for details.
- qpsolvers.solvers.piqp_.piqp_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, backend=None, **kwargs)¶
Solve a quadratic program using PIQP.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{\mbox{minimize}}{x} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using PIQP.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Positive semidefinite cost matrix.q (
Union
[ndarray
,csc_matrix
]) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint vector.lb (
Union
[ndarray
,csc_matrix
,None
]) – Lower bound constraint vector.ub (
Union
[ndarray
,csc_matrix
,None
]) – Upper bound constraint vector.backend (
Optional
[str
]) – PIQP backend to use in[None, "dense", "sparse"]
. IfNone
(default), the backend is selected based on the type ofP
.verbose (
bool
) – Set to True to print out extra information.initvals (
Optional
[ndarray
]) – Warm-start guess vector. Not used.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.
ProxQP¶
Solver interface for ProxQP.
ProxQP is the QP solver from ProxSuite, a collection of open-source solvers rooted in revisited primal-dual proximal algorithms. If you use ProxQP in some academic work, consider citing the corresponding paper [Bambade2022].
- qpsolvers.solvers.proxqp_.proxqp_solve_problem(problem, initvals=None, verbose=False, backend=None, **kwargs)¶
Solve a quadratic program using ProxQP.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.backend (
Optional
[str
]) – ProxQP backend to use in[None, "dense", "sparse"]
. IfNone
(default), the backend is selected based on the type ofP
.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP returned by the solver.
- Raises:
ParamError – If a warm-start value is given both in initvals and the x keyword argument.
Notes
All other keyword arguments are forwarded as options to ProxQP. For instance, you can call
proxqp_solve_qp(P, q, G, h, eps_abs=1e-6)
. For a quick overview, the solver accepts the following settings:Name
Effect
x
Warm start value for the primal variable.
y
Warm start value for the dual Lagrange multiplier for equality constraints.
z
Warm start value for the dual Lagrange multiplier for inequality constraints.
eps_abs
Asbolute stopping criterion of the solver (default: 1e-3, note that this is a laxer default than other solvers). See e.g. [tolerances] for an overview of solver tolerances.
eps_rel
Relative stopping criterion of the solver. See e.g. [tolerances] for an overview of solver tolerances.
mu_eq
Proximal step size wrt equality constraints multiplier.
mu_in
Proximal step size wrt inequality constraints multiplier.
rho
Proximal step size wrt primal variable.
compute_preconditioner
If
True
(default), the preconditioner will be derived.compute_timings
If
True
(default), timings will be computed by the solver (setup time, solving time, and run time = setup time + solving time).max_iter
Maximal number of authorized outer iterations.
initial_guess
Sets the initial guess option for initilizing x, y and z.
This list is not exhaustive. Check out the solver documentation for details.
- qpsolvers.solvers.proxqp_.proxqp_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, backend=None, **kwargs)¶
Solve a quadratic program using ProxQP.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{\mbox{minimize}}{x} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using ProxQP.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Positive semidefinite cost matrix.q (
Union
[ndarray
,csc_matrix
]) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint vector.lb (
Union
[ndarray
,csc_matrix
,None
]) – Lower bound constraint vector.ub (
Union
[ndarray
,csc_matrix
,None
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.backend (
Optional
[str
]) – ProxQP backend to use in[None, "dense", "sparse"]
. IfNone
(default), the backend is selected based on the type ofP
.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.
QPALM¶
Solver interface for QPALM.
C implementation of QPALM, a proximal augmented Lagrangian based solver for (possibly nonconvex) quadratic programs. If you use QPALM in some academic work, consider citing the corresponding paper [Hermans2022].
- qpsolvers.solvers.qpalm_.qpalm_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using QPALM.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP returned by the solver.
- Raises:
ParamError – If a warm-start value is given both in initvals and the x keyword argument.
Note
QPALM internally only uses the upper-triangular part of the cost matrix \(P\).
Notes
Keyword arguments are forwarded as “settings” to QPALM. For instance, we can call
qpalm_solve_qp(P, q, G, h, u, eps_abs=1e-4, eps_rel=1e-4)
.Name
Effect
max_iter
Maximum number of iterations.
eps_abs
Asbolute stopping criterion of the solver. See e.g. [tolerances] for an overview of solver tolerances.
eps_rel
Relative stopping criterion of the solver. See e.g. [tolerances] for an overview of solver tolerances.
rho
Tolerance scaling factor.
theta
Penalty update criterion parameter.
delta
Penalty update factor.
sigma_max
Penalty factor cap.
proximal
Boolean, use proximal method of multipliers or not.
This list is not exhaustive. Check out the solver documentation for details.
- qpsolvers.solvers.qpalm_.qpalm_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using QPALM.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{\mbox{minimize}}{x} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using QPALM.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Positive semidefinite cost matrix.q (
Union
[ndarray
,csc_matrix
]) – Cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint vector.lb (
Union
[ndarray
,csc_matrix
,None
]) – Lower bound constraint vector.ub (
Union
[ndarray
,csc_matrix
,None
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.
qpOASES¶
Solver interface for qpOASES.
qpOASES is an open-source C++ implementation of the online active set strategy, which was inspired by observations from the field of parametric quadratic programming. It has theoretical features that make it suitable to model predictive control. Further numerical modifications have made qpOASES a reliable QP solver, even when tackling semi-definite, ill-posed or degenerated QP problems. If you are using qpOASES in some academic work, consider citing the corresponding paper [Ferreau2014].
See the installation page for additional instructions on installing this solver.
- qpsolvers.solvers.qpoases_.qpoases_solve_problem(problem, initvals=None, verbose=False, max_wsr=1000, time_limit=None, predefined_options=None, **kwargs)¶
Solve a quadratic program using qpOASES.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.max_wsr (
int
) – Maximum number of Working-Set Recalculations given to qpOASES.time_limit (
Optional
[float
]) – Set a run time limit in seconds.predefined_options (
Optional
[str
]) – Set solver options to one of the pre-defined choices provided by qpOASES, to pick in["default", "fast", "mpc", "reliable"]
.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError : – If the problem is ill-formed in some way, for instance if some matrices are not dense.
ValueError : – If
predefined_options
is not a valid choice.
Notes
This function relies on an update to qpOASES to allow empty bounds (lb, ub, lbA or ubA) in Python. This is possible in the C++ API but not by the Python API (as of version 3.2.0). If using them, be sure to update the Cython file (qpoases.pyx) in your distribution of qpOASES to convert
None
to the null pointer. Check out the installation instructions.Keyword arguments are forwarded as options to qpOASES. For instance, we can call
qpoases_solve_qp(P, q, G, h, u, terminationTolerance=1e-14)
. qpOASES options include the following:Name
Description
boundRelaxation
Initial relaxation of bounds to start homotopy and initial value for far bounds.
epsNum
Numerator tolerance for ratio tests.
epsDen
Denominator tolerance for ratio tests.
numRefinementSteps
Maximum number of iterative refinement steps.
terminationTolerance
Relative termination tolerance to stop homotopy.
Check out pages 28 to 30 of qpOASES User’s Manual. for all available options.
- qpsolvers.solvers.qpoases_.qpoases_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, max_wsr=1000, time_limit=None, predefined_options=None, **kwargs)¶
Solve a quadratic program using qpOASES.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using qpOASES.
- Parameters:
P (
ndarray
) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.max_wsr (
int
) – Maximum number of Working-Set Recalculations given to qpOASES.time_limit (
Optional
[float
]) – Set a run time limit in seconds.predefined_options (
Optional
[str
]) – Set solver options to one of the pre-defined choices provided by qpOASES, to pick in["default", "fast", "mpc", "reliable"]
.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError : – If the problem is ill-formed in some way, for instance if some matrices are not dense.
ValueError : – If
predefined_options
is not a valid choice.
Notes
This function relies on an update to qpOASES to allow empty bounds (lb, ub, lbA or ubA) in Python. This is possible in the C++ API but not by the Python API (as of version 3.2.0). If using them, be sure to update the Cython file (qpoases.pyx) in your distribution of qpOASES to convert
None
to the null pointer. Check out the installation instructions.Keyword arguments are forwarded as options to qpOASES. For instance, we can call
qpoases_solve_qp(P, q, G, h, u, terminationTolerance=1e-14)
. qpOASES options include the following:Name
Description
boundRelaxation
Initial relaxation of bounds to start homotopy and initial value for far bounds.
epsNum
Numerator tolerance for ratio tests.
epsDen
Denominator tolerance for ratio tests.
numRefinementSteps
Maximum number of iterative refinement steps.
terminationTolerance
Relative termination tolerance to stop homotopy.
Check out pages 28 to 30 of qpOASES User’s Manual. for all available options.
qpSWIFT¶
Solver interface for qpSWIFT.
qpSWIFT is a light-weight sparse Quadratic Programming solver targeted for embedded and robotic applications. It employs Primal-Dual Interior Point method with Mehrotra Predictor corrector step and Nesterov Todd scaling. For solving the linear system of equations, sparse LDL’ factorization is used along with approximate minimum degree heuristic to minimize fill-in of the factorizations. If you use qpSWIFT in your research, consider citing the corresponding paper [Pandala2019].
- qpsolvers.solvers.qpswift_.qpswift_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using qpSWIFT.
Note
This solver does not handle problems without inequality constraints.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
- Raises:
ProblemError : – If the problem is ill-formed in some way, for instance if some matrices are not dense or the problem has no inequality constraint.
Note
Rank assumptions: qpSWIFT requires the QP matrices to satisfy the
\[\begin{split}\begin{array}{cc} \mathrm{rank}(A) = p & \mathrm{rank}([P\ A^T\ G^T]) = n \end{array}\end{split}\]where \(p\) is the number of rows of \(A\) and \(n\) is the number of optimization variables. This is the same requirement as
cvxopt_solve_qp()
, however qpSWIFT does not perform rank checks as it prioritizes performance. If the solver fails on your problem, try running CVXOPT on it for rank checks.Notes
All other keyword arguments are forwarded as options to the qpSWIFT solver. For instance, you can call
qpswift_solve_qp(P, q, G, h, ABSTOL=1e-5)
. See the solver documentation for details.For a quick overview, the solver accepts the following settings:
Name
Effect
MAXITER
Maximum number of iterations needed.
ABSTOL
Absolute tolerance on the duality gap. See e.g. [tolerances] for a primer on the duality gap and solver tolerances.
RELTOL
Relative tolerance on the residuals \(r_x = P x + G^T z + q\) (dual residual), \(r_y = A x - b\) (primal residual on equality constraints) and \(r_z = h - G x - s\) (primal residual on inequality constraints). See equation (21) in [Pandala2019].
SIGMA
Maximum centering allowed.
If a verbose output shows that the maximum number of iterations is reached, check e.g. (1) the rank of your equality constraint matrix and (2) that your inequality constraint matrix does not have zero rows.
As qpSWIFT does not sanity check its inputs, it should be used with a little more care than the other solvers. For instance, make sure you don’t have zero rows in your input matrices, as it can make the solver numerically unstable.
- qpsolvers.solvers.qpswift_.qpswift_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using qpSWIFT.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using qpSWIFT.
Note
This solver does not handle problems without inequality constraints yet.
- Parameters:
P (
ndarray
) – Symmetric cost matrix. Together with \(A\) and \(G\), it should satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\), see the rank assumptions below.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix. Together with \(P\) and \(A\), it should satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\), see the rank assumptions below.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix. It needs to be full row rank, and together with \(P\) and \(G\) satisfy \(\mathrm{rank}([P\ A^T\ G^T]) = n\). See the rank assumptions below.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector.verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError : – If the problem is ill-formed in some way, for instance if some matrices are not dense or the problem has no inequality constraint.
Note
Rank assumptions: qpSWIFT requires the QP matrices to satisfy the
\[\begin{split}\begin{array}{cc} \mathrm{rank}(A) = p & \mathrm{rank}([P\ A^T\ G^T]) = n \end{array}\end{split}\]where \(p\) is the number of rows of \(A\) and \(n\) is the number of optimization variables. This is the same requirement as
cvxopt_solve_qp()
, however qpSWIFT does not perform rank checks as it prioritizes performance. If the solver fails on your problem, try running CVXOPT on it for rank checks.Notes
All other keyword arguments are forwarded as options to the qpSWIFT solver. For instance, you can call
qpswift_solve_qp(P, q, G, h, ABSTOL=1e-5)
. See the solver documentation for details.For a quick overview, the solver accepts the following settings:
Name
Effect
MAXITER
Maximum number of iterations needed.
ABSTOL
Absolute tolerance on the duality gap. See e.g. [tolerances] for a primer on the duality gap and solver tolerances.
RELTOL
Relative tolerance on the residuals \(r_x = P x + G^T z + q\) (dual residual), \(r_y = A x - b\) (primal residual on equality constraints) and \(r_z = h - G x - s\) (primal residual on inequality constraints). See equation (21) in [Pandala2019].
SIGMA
Maximum centering allowed.
If a verbose output shows that the maximum number of iterations is reached, check e.g. (1) the rank of your equality constraint matrix and (2) that your inequality constraint matrix does not have zero rows.
As qpSWIFT does not sanity check its inputs, it should be used with a little more care than the other solvers. For instance, make sure you don’t have zero rows in your input matrices, as it can make the solver numerically unstable.
quadprog¶
Solver interface for quadprog.
quadprog is a C implementation of the Goldfarb-Idnani dual algorithm [Goldfarb1983]. It works best on well-conditioned dense problems.
- qpsolvers.solvers.quadprog_.quadprog_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using quadprog.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ProblemError : – If the cost matrix of the quadratic program if not positive definite, or if the problem is ill-formed in some way, for instance if some matrices are not dense.
Note
The quadprog solver only considers the lower entries of \(P\), therefore it will use a different cost than the one intended if a non-symmetric matrix is provided.
Notes
All other keyword arguments are forwarded to the quadprog solver. For instance, you can call
quadprog_solve_qp(P, q, G, h, factorized=True)
. See the solver documentation for details.
- qpsolvers.solvers.quadprog_.quadprog_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using quadprog.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using quadprog.
- Parameters:
P (
ndarray
) – Symmetric cost matrix.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.
SCS¶
Solver interface for SCS.
SCS (Splitting Conic Solver) is a numerical optimization package for solving large-scale convex quadratic cone problems, which is a general class of problems that includes quadratic programming. If you use SCS in some academic work, consider citing the corresponding paper [ODonoghue2021].
- qpsolvers.solvers.scs_.scs_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using SCS.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution returned by the solver.
- Raises:
ValueError – If the quadratic program is not unbounded below.
Notes
Keyword arguments are forwarded as is to SCS. For instance, we can call
scs_solve_qp(P, q, G, h, eps_abs=1e-6, eps_rel=1e-4)
. SCS settings include the following:Name
Description
max_iters
Maximum number of iterations to run.
time_limit_secs
Time limit for solve run in seconds (can be fractional). 0 is interpreted as no limit.
eps_abs
Absolute feasibility tolerance. See Termination criteria.
eps_rel
Relative feasibility tolerance. See Termination criteria.
eps_infeas
Infeasibility tolerance (primal and dual), see Certificate of infeasibility.
normalize
Whether to perform heuristic data rescaling. See Data equilibration.
Check out the SCS settings documentation for all available settings.
- qpsolvers.solvers.scs_.scs_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using SCS.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{x}{\mbox{minimize}} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using SCS.
- Parameters:
P (
Union
[ndarray
,csc_matrix
]) – Primal quadratic cost matrix.q (
ndarray
) – Primal quadratic cost vector.G (
Union
[ndarray
,csc_matrix
,None
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Union
[ndarray
,csc_matrix
,None
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
Optional
[ndarray
]- Returns:
Solution to the QP, if found, otherwise
None
.- Raises:
ValueError – If the quadratic program is not unbounded below.
Notes
Keyword arguments are forwarded as is to SCS. For instance, we can call
scs_solve_qp(P, q, G, h, eps_abs=1e-6, eps_rel=1e-4)
. SCS settings include the following:Name
Description
max_iters
Maximum number of iterations to run.
time_limit_secs
Time limit for solve run in seconds (can be fractional). 0 is interpreted as no limit.
eps_abs
Absolute feasibility tolerance. See Termination criteria.
eps_rel
Relative feasibility tolerance. See Termination criteria.
eps_infeas
Infeasibility tolerance (primal and dual), see Certificate of infeasibility.
normalize
Whether to perform heuristic data rescaling. See Data equilibration.
Check out the SCS settings documentation for all available settings.
qpax ===
Solver interface for qpax.
qpax is an open source QP solver that can be combined with JAX’s jit and vmap functionality, as well as differentiated with reverse-mode differentiation. It is based on a primal-dual interior point algorithm. If you are using qpax in some academic work, consider citing the corresponding paper [Tracy2024].
- qpsolvers.solvers.qpax_.qpax_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using qpax.
- Parameters:
problem (
Problem
) – Quadratic program to solve.initvals (
Optional
[ndarray
]) – Warm-start guess vector (not used).verbose (
bool
) – Set to True to print out extra information.
- Return type:
- Returns:
Solution to the QP returned by the solver.
Notes
All other keyword arguments are forwarded as options to qpax. For instance, you can call
qpax_solve_qp(P, q, G, h, solver_tol=1e-5)
. For a quick overview, the solver accepts the following settings:Name
Effect
solver_tol
Tolerance for the solver.
Note that jax by default uses 32-bit floating point numbers, which can lead to numerical instability. If you encounter numerical issues, consider using 64-bit floating point numbers by setting
`python import jax jax.config.update("jax_enable_x64", True) `
- qpsolvers.solvers.qpax_.qpax_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using qpax.
The quadratic program is defined as:
\[\begin{split}\begin{split}\begin{array}{ll} \underset{\mbox{minimize}}{x} & \frac{1}{2} x^T P x + q^T x \\ \mbox{subject to} & G x \leq h \\ & A x = b \\ & lb \leq x \leq ub \end{array}\end{split}\end{split}\]It is solved using qpax. Paper:.
- Parameters:
P (
ndarray
) – Positive semidefinite cost matrix.q (
ndarray
) – Cost vector.G (
Optional
[ndarray
]) – Linear inequality constraint matrix.h (
Optional
[ndarray
]) – Linear inequality constraint vector.A (
Optional
[ndarray
]) – Linear equality constraint matrix.b (
Optional
[ndarray
]) – Linear equality constraint vector.lb (
Optional
[ndarray
]) – Lower bound constraint vector.ub (
Optional
[ndarray
]) – Upper bound constraint vector.verbose (
bool
) – Set to True to print out extra information.initvals (
Optional
[ndarray
]) – Warm-start guess vector. Not used.
- Return type:
Optional
[ndarray
]- Returns:
Primal solution to the QP, if found, otherwise
None
.