Skip to contents

Solves $$arg\min_x 0.5 x'P x + c'x$$ s.t. $$A x = b$$ $$G x \leq h$$ $$x_{lb} \leq x \leq x_{ub}$$ for real matrices P (nxn, positive semidefinite), A (pxn) with m number of equality constraints, and G (mxn) with m number of inequality constraints

Usage

solve_piqp(
  P = NULL,
  c = NULL,
  A = NULL,
  b = NULL,
  G = NULL,
  h = NULL,
  x_lb = NULL,
  x_ub = NULL,
  settings = list(),
  backend = c("auto", "sparse", "dense")
)

Arguments

P

dense or sparse matrix of class dgCMatrix or coercible into such, must be positive semidefinite

c

numeric vector

A

dense or sparse matrix of class dgCMatrix or coercible into such

b

numeric vector

G

dense or sparse matrix of class dgCMatrix or coercible into such

h

numeric vector

x_lb

a numeric vector of lower bounds, default NULL indicating -Inf for all variables, otherwise should be number of variables long

x_ub

a numeric vector of upper bounds, default NULL indicating Inf for all variables, otherwise should be number of variables long

settings

list with optimization parameters, empty by default; see piqp_settings() for a comprehensive list of parameters that may be used

backend

which backend to use, if auto and P, A or G are sparse then sparse backend is used ("auto", "sparse" or "dense") ("auto")

Value

A list with elements solution elements

References

Schwan, R., Jiang, Y., Kuhn, D., Jones, C.N. (2023). ``PIQP: A Proximal Interior-Point Quadratic Programming Solver.'' doi:10.48550/arXiv.2304.00290

See also

piqp(), piqp_settings() and the underlying PIQP documentation: https://predict-epfl.github.io/piqp/

Examples

## example, adapted from PIQP documentation
library(piqp)
library(Matrix)

P <- Matrix(c(6., 0.,
              0., 4.), 2, 2, sparse = TRUE)
c <- c(-1., -4.)
A <- Matrix(c(1., -2.), 1, 2, sparse = TRUE)
b <- c(1.)
G <- Matrix(c(1., 2., -1., 0.), 2, 2, sparse = TRUE)
h <- c(0.2, -1.)
x_lb <- c(-1., -Inf)
x_ub <- c(1., Inf)

settings <- list(verbose = TRUE)

# Solve with PIQP
res <- solve_piqp(P, c, A, b, G, h, x_lb, x_ub, settings)
#> ----------------------------------------------------------
#>                            PIQP                           
#>                     (c) Roland Schwan                     
#>    Ecole Polytechnique Federale de Lausanne (EPFL) 2023   
#> ----------------------------------------------------------
#> sparse backend
#> variables n = 2, nzz(P upper triangular) = 2
#> equality constraints p = 1, nnz(A) = 2
#> inequality constraints m = 2, nnz(G) = 3
#> variable lower bounds n_lb = 1
#> variable upper bounds n_ub = 1
#> 
#> iter  prim_obj       dual_obj       duality_gap   prim_inf      dual_inf      rho         delta       mu          p_step   d_step
#>   0    2.43682e+00   -2.58316e+00   5.01997e+00   2.22782e+00   2.07114e+01   1.000e-06   1.000e-04   1.068e+01   0.0000   0.0000
#>   1    6.54472e+00    3.11191e+00   3.43282e+00   2.20398e-02   2.07118e-01   8.792e-08   8.792e-06   9.387e-01   0.9900   0.9900
#>   2    6.11927e+00    5.84978e+00   2.69498e-01   4.47812e-03   5.52769e-02   6.007e-09   6.007e-07   6.414e-02   0.7962   0.9607
#>   3    6.16036e+00    6.15648e+00   3.87581e-03   4.43795e-05   8.10844e-04   1.000e-10   8.432e-09   9.003e-04   0.9900   0.9870
#>   4    6.16000e+00    6.15996e+00   3.87634e-05   4.43739e-07   8.10934e-06   1.000e-10   1.000e-10   9.007e-06   0.9900   0.9900
#>   5    6.16000e+00    6.16000e+00   3.87671e-07   4.43738e-09   8.10935e-08   1.000e-10   1.000e-10   9.008e-08   0.9900   0.9900
#>   6    6.16000e+00    6.16000e+00   3.87708e-09   4.43737e-11   8.10932e-10   1.000e-10   1.000e-10   9.008e-10   0.9900   0.9900
#> 
#> status:               solved
#> number of iterations: 6
#> objective:            6.16000e+00
res$x
#> [1] -0.6 -0.8