Unsupported solvers¶
Unsupported solvers will be made available if they are detected on your system, but their performance is not guaranteed as they are not part of continuous integration (typically because they are not open source).
PDHCG¶
Solver interface for PDHCG.
PDHCG (Primal-Dual Hybrid Conjugate Gradient) is a high-performance, GPU-accelerated solver designed for large-scale convex Quadratic Programming (QP). It is particularly efficient for huge-scale problems by fully leveraging NVIDIA CUDA architectures.
Note
To use this solver, you need an NVIDIA GPU and the pdhcg package
installed via pip install pdhcg. For advanced installation (e.g.,
custom CUDA paths), please refer to the
official PDHCG-II repository.
References
- qpsolvers.solvers.pdhcg_.pdhcg_solve_problem(problem, initvals=None, verbose=False, **kwargs)¶
Solve a quadratic program using PDHCG.
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}\]- 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 to the QP, if found, otherwise
None.
Notes
Keyword arguments are forwarded to PDHCG as solver parameters. For instance, you can call
pdhcg_solve_qp(..., TimeLimit=60). Common PDHCG parameters include:Name
Description
TimeLimitMaximum wall-clock time in seconds (default: 3600.0).
IterationLimitMaximum number of iterations.
OptimalityTolRelative tolerance for optimality gap (default: 1e-4).
FeasibilityTolRelative feasibility tolerance for residuals (default: 1e-4).
OutputFlagEnable (True) or disable (False) console logging output.
For advanced parameters, please refer to the PDHCG Documentation.
- qpsolvers.solvers.pdhcg_.pdhcg_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 PDHCG.
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 PDHCG.
- 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 PDHCG as solver parameters. For instance, you can call
pdhcg_solve_qp(..., TimeLimit=60). Common PDHCG parameters include:Name
Description
TimeLimitMaximum wall-clock time in seconds (default: 3600.0).
IterationLimitMaximum number of iterations.
OptimalityTolRelative tolerance for optimality gap (default: 1e-4).
FeasibilityTolRelative feasibility tolerance for residuals (default: 1e-4).
OutputFlagEnable (True) or disable (False) console logging output.
For advanced parameters, please refer to the PDHCG Documentation.
NPPro¶
Solver interface for NPPro.
The NPPro solver implements an enhanced Newton Projection with Proportioning method for strictly convex quadratic programming. Currently, it is designed for dense problems only.
- qpsolvers.solvers.nppro_.nppro_solve_problem(problem, initvals=None, **kwargs)¶
Solve a quadratic program using NPPro.
- Parameters:
problem (
Problem) – Quadratic program to solve.initvals (
Optional[ndarray]) – Warm-start guess vector.
- Return type:
- Returns:
Solution returned by the solver.
Notes
All other keyword arguments are forwarded as options to NPPro. For instance, you can call
nppro_solve_qp(P, q, G, h, MaxIter=15). For a quick overview, the solver accepts the following settings:Name
Effect
MaxIterMaximum number of iterations.
SkipPreprocessingSkip preprocessing phase or not.
SkipPhaseOneSkip feasible starting point finding or not.
InfValValues are assumed to be infinite above this threshold.
HessianUpdatesEnable Hessian updates or not.
- qpsolvers.solvers.nppro_.nppro_solve_qp(P, q, G=None, h=None, A=None, b=None, lb=None, ub=None, initvals=None, **kwargs)¶
Solve a quadratic program using NPPro.
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 NPPro.
- Parameters:
P (
ndarray) – Positive definite 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.
- Return type:
Optional[ndarray]- Returns:
Solution to the QP, if found, otherwise
None.
Notes
See the Notes section in
nppro_solve_problem().