PIQP Solver object
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
indicatingInf
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
An R6-object of class "piqp_model" with methods defined which can be further used to solve the problem with updated settings / parameters.
Details
Allows one to solve a parametric
problem with for example warm starts between updates of the parameter, c.f. the examples.
The object returned by piqp
contains several methods which can be used to either update/get details of the
problem, modify the optimization settings or attempt to solve the problem.
Usage
model = piqp(P = NULL, c = NULL, A = NULL, b = NULL, G = NULL, h = NULL, x_lb = NULL, x_ub = NULL, settings = piqp_settings(), backend = c("auto", "sparse", "dense"))
model$solve()
model$update(P = NULL, c = NULL, A = NULL, b = NULL, G = NULL, h = NULL, x_lb = NULL, x_ub = NULL)
model$get_settings()
model$get_dims()
model$update_settings(new_settings = piqp_settings())
print(model)
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)
model <- piqp(P, c, A, b, G, h, x_lb, x_ub, settings)
# Solve
res <- model$solve()
#> ----------------------------------------------------------
#> 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
# Define new data
A_new <- Matrix(c(1., -3.), 1, 2, sparse = TRUE)
h_new <- c(2., 1.)
# Update model and solve again
model$update(A = A_new, h = h_new)
res <- model$solve()
#> ----------------------------------------------------------
#> 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 1.47500e+01 -2.00344e+02 2.15094e+02 8.99959e+00 9.45659e+01 1.000e-10 1.000e-10 2.067e+02 0.0000 0.0000
#> 1 2.19193e+00 -9.18795e+01 9.40714e+01 8.99959e-02 9.45657e-01 1.000e-10 1.000e-10 2.466e+01 0.9900 0.9900
#> 2 1.94403e+00 -2.40809e+00 4.35212e+00 8.99959e-04 5.47442e-02 1.000e-10 1.000e-10 1.089e+00 0.9900 0.9563
#> 3 1.05672e+00 5.35890e-01 5.20827e-01 8.99954e-06 5.48689e-04 1.000e-10 1.000e-10 1.302e-01 0.9900 0.9900
#> 4 9.63579e-01 9.18538e-01 4.50412e-02 8.99741e-08 4.99746e-06 1.000e-10 1.000e-10 1.126e-02 0.9900 0.9900
#> 5 9.57033e-01 9.53468e-01 3.56551e-03 8.91871e-10 1.69943e-07 1.000e-10 1.000e-10 8.914e-04 0.9900 0.9900
#> 6 9.56897e-01 9.56828e-01 6.90282e-05 7.56337e-12 4.20226e-09 1.000e-10 1.000e-10 1.726e-05 0.9900 0.9900
#> 7 9.56897e-01 9.56896e-01 6.89596e-07 5.49193e-14 2.38823e-10 1.000e-10 1.000e-10 1.724e-07 0.9900 0.9900
#> 8 9.56897e-01 9.56897e-01 6.89545e-09 4.54040e-16 5.93214e-12 1.000e-10 1.000e-10 1.724e-09 0.9900 0.9900
#>
#> status: solved
#> number of iterations: 8
#> objective: 9.56897e-01
res$x
#> [1] 0.4310345 -0.1896552