piqp {piqp}R Documentation

PIQP Solver object

Description

PIQP Solver object

Usage

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")

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.

Value

An R6-object of class "piqp_model" with methods defined which can be further used to solve the problem with updated settings / parameters.

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)

See Also

solve_piqp(), piqp_settings()

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()
res$x

# 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()
res$x


[Package piqp version 0.2.2 Index]