1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163
|
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/ecos.R
\name{ECOS_csolve}
\alias{ECOS_csolve}
\title{Solve a conic optimization problem}
\usage{
ECOS_csolve(
c = numeric(0),
G = NULL,
h = numeric(0),
dims = list(l = integer(0), q = NULL, e = integer(0)),
A = NULL,
b = numeric(0),
bool_vars = integer(0),
int_vars = integer(0),
control = ecos.control()
)
}
\arguments{
\item{c}{the coefficients of the objective function; the length of
this determines the number of variables \eqn{n} in the problem.}
\item{G}{the inequality constraint matrix in one of three forms: a
plain matrix, simple triplet matrix, or compressed column
format, e.g. \link[Matrix]{dgCMatrix-class}. Can also be
\code{NULL}}
\item{h}{the right hand side of the inequality constraint. Can be
empty numeric vector.}
\item{dims}{is a list of three named elements: \code{dims['l']} an
integer specifying the dimension of positive orthant cone,
\code{dims['q']} an integer vector specifying dimensions of
second-order cones, \code{dims['e']} an integer specifying the
number of exponential cones}
\item{A}{the optional equality constraint matrix in one of three
forms: a plain matrix, simple triplet matrix, or compressed
column format, e.g. \link[Matrix]{dgCMatrix-class}. Can be
\code{NULL}}
\item{b}{the right hand side of the equality constraint, must be
specified if \eqn{A} is. Can be empty numeric vector.}
\item{bool_vars}{the indices of the variables, 1 through \eqn{n},
that are boolean; that is, they are either present or absent in
the solution}
\item{int_vars}{the indices of the variables, 1 through \eqn{n},
that are integers}
\item{control}{is a named list that controls various optimization
parameters; see \link[ECOSolveR]{ecos.control}.}
}
\value{
a list of 8 named items
\describe{
\item{x}{primal variables}
\item{y}{dual variables for equality constraints}
\item{s}{slacks for \eqn{Gx + s <= h}, \eqn{s \in K}}
\item{z}{dual variables for inequality constraints \eqn{s \in K}}
\item{infostring}{gives information about the status of solution}
\item{retcodes}{a named integer vector containing four elements
\describe{
\item{exitflag}{0=\code{ECOS_OPTIMAL}, 1=\code{ECOS_PINF},
2=\code{ECOS_DINF}, 10=\code{ECOS_INACC_OFFSET}, -1=\code{ECOS_MAXIT},
-2=\code{ECOS_NUMERICS}, -3=\code{ECOS_OUTCONE}, -4=\code{ECOS_SIGINT},
-7=\code{ECOS_FATAL}. See \link[ECOSolveR]{ECOS_exitcodes}}.
\item{iter}{the number of iterations used}
\item{mi_iter}{the number of iterations for mixed integer problems}
\item{numerr}{a non-zero number if a numeric error occurred}
}
}
\item{summary}{a named numeric vector containing
\describe{
\item{pcost}{value of primal objective}
\item{dcost}{value of dual objective}
\item{pres}{primal residual on inequalities and equalities}
\item{dres}{dual residual}
\item{pinf}{primal infeasibility measure}
\item{dinf}{dual infeasibility measure}
\item{pinfres}{primal infeasibility residual}
\item{dinfres}{dual infeasibility residual}
\item{gap}{duality gap}
\item{relgap}{relative duality gap}
\item{r0}{Unknown at the moment to this R package maintainer.}
}
}
\item{timing}{a named numeric vector of timing information consisting of
\describe{
\item{runtime}{the total runtime in ecos}
\item{tsetup}{the time for setup of the problem}
\item{tsolve}{the time to solve the problem}
}
}
}
}
\description{
The function \code{ECOS_csolve} is a wrapper around the ecos
\code{csolve} C function. Conic constraints are specified using the
\eqn{G} and \eqn{h} parameters and can be \code{NULL} and zero
length vector respectively indicating an absence of conic
constraints. Similarly, equality constraints are specified via
\eqn{A} and \eqn{b} parameters with \code{NULL} and empty vector
values representing a lack of such constraints. At most one of the
pair \eqn{(G , h)} or \eqn{(A, b)} is allowed to be absent.
}
\section{Details}{
A call to this function will solve the problem:
minimize \eqn{c^Tx}, subject to \eqn{Ax = b}, and \eqn{h - G*x \in K}.
Variables can be constrained to be boolean (1 or 0) or integers. This is indicated
by specifying parameters \code{bool_vars} and/or \code{int_vars} respectively. If so
indicated, the solutions will be found using a branch and bound algorithm.
}
\examples{
## githubIssue98
cat("Basic matrix interface\n")
Gmat <- matrix(c(0.416757847405471, 2.13619609566845, 1.79343558519486, 0, 0,
0, 0, -1, 0, 0, 0, 0.056266827226329, -1.64027080840499, 0.841747365656204,
0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0.416757847405471, 2.13619609566845,
1.79343558519486, 0, 0, 0, -1, 0, 0, 0, 0, 0.056266827226329, -1.64027080840499,
0.841747365656204, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0), ncol = 5L)
c <- as.numeric(c(0, 0, 0, 0, 1))
h <- as.numeric(c(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0))
dims <- list(l = 6L, q = 5L, e = 0L)
ECOS_csolve(c = c, G = Gmat, h = h,
dims = dims,
A = NULL, b = numeric(0))
cat("Simple Triplet Matrix interface, if you have package slam\n")
if (requireNamespace("slam")) {
ECOS_csolve(c = c, G = slam::as.simple_triplet_matrix(Gmat), h = h,
dims = dims,
A = NULL, b = numeric(0))
}
if (requireNamespace("Matrix")) {
ECOS_csolve(c = c, G = Matrix::Matrix(Gmat), h = h,
dims = dims,
A = NULL, b = numeric(0))
}
## Larger problems using saved data can be found in the test suite.
## Here is one
if (requireNamespace("Matrix")) {
MPC01 <- readRDS(system.file("testdata", "MPC01_1.RDS", package = "ECOSolveR"))
G <- Matrix::sparseMatrix(x = MPC01$Gpr, i = MPC01$Gir, p = MPC01$Gjc,
dims = c(MPC01$m, MPC01$n), index1 = FALSE)
h <- MPC01$h
dims <- lapply(list(l = MPC01$l, q=MPC01$q, e=MPC01$e), as.integer)
retval <- ECOS_csolve(c = MPC01$c, G=G, h = h, dims = dims, A = NULL, b = NULL,
control = ecos.control(verbose=1L))
retval$retcodes
retval$infostring
retval$summary
}
}
|