| Type: | Package |
| Title: | Conic Benchmark Format (CBF) Plugin for the R Optimization Infrastructure |
| Version: | 0.1-0 |
| Date: | 2026-03-19 |
| Author: | Benjamin Schwendinger [aut, cre] |
| Maintainer: | Benjamin Schwendinger <benjaminschwe@gmail.com> |
| Description: | Provides read and write support for the Conic Benchmark Format (CBF, version 4) within the R Optimization Infrastructure ('ROI'). Supported cone types include the positive orthant, second-order (SOC), rotated second-order (bridged automatically to standard SOC), exponential (primal and dual), power (primal and dual), and semidefinite (symmetric-vectorised) cones, as well as their mixed-integer variants. The reader translates a .cbf file into an ROI 'OP' object, handling coordinate-convention differences between CBF and ROI transparently; the writer serialises an ROI 'OP' object back to CBF plain-text. |
| License: | GPL-3 |
| URL: | https://gitlab.com/roigrp/tools/roi.format.cbf |
| BugReports: | https://gitlab.com/roigrp/tools/roi.format.cbf/-/issues |
| Imports: | ROI (≥ 1.0-0), slam (≥ 0.1-40) |
| Suggests: | tinytest, ROI.plugin.ecos, ROI.plugin.scs |
| Encoding: | UTF-8 |
| RoxygenNote: | 7.3.2 |
| NeedsCompilation: | no |
| Packaged: | 2026-04-18 14:48:05 UTC; ben |
| Repository: | CRAN |
| Date/Publication: | 2026-04-21 19:42:20 UTC |
Convert variable bounds to a list of CBF VAR cone chunks.
Description
Each "chunk" groups consecutive variables that share the same conic domain. Run-length encoding keeps the VAR section compact.
Usage
.bounds_to_var_chunks(lb, ub)
Arguments
lb |
Numeric vector of lower bounds (length = number of variables). |
ub |
Numeric vector of upper bounds (same length). |
Value
List of list(cone, count) pairs.
Bound -> cone mapping
lb=0, ub=+Inf-
L+(positive orthant) lb=-Inf, ub=0-
L-(negative orthant) lb=0, ub=0-
L=(fixpoint zero) - anything else
-
F.(free)
Reorder a CBF exponential-cone chunk to ROI's coordinate convention.
Description
Reorder a CBF exponential-cone chunk to ROI's coordinate convention.
Usage
.cbf_bridge_exp_cone(A_sub, rhs_sub)
Bridge a rotated quadratic cone chunk to a standard SOC chunk.
Description
Bridge a rotated quadratic cone chunk to a standard SOC chunk.
Usage
.cbf_bridge_rotated_soc(A_sub, rhs_sub)
Convert parsed CBF sections into an ROI OP object.
Description
Convert parsed CBF sections into an ROI OP object.
Usage
.cbf_build_OP(secs)
Build ROI constraint objects from the parsed CON chunks.
Description
Iterates over the cone chunks defined in the CON section. For each chunk it extracts the corresponding rows from the full sparse constraint matrix and RHS vector, then constructs the appropriate ROI constraint object.
Usage
.cbf_build_constraints(con_chunks, A_mat, rhs_vec, powcones, powstar_cones)
Arguments
con_chunks |
List of |
A_mat |
Full sparse constraint matrix (slam STM). |
rhs_vec |
Full RHS vector (already negated: |
powcones |
List of alpha vectors from the POWCONES section. |
powstar_cones |
List of alpha vectors from the POW*CONES section. |
Value
List of ROI constraint objects.
Cone chunk -> ROI constraint mapping
L+(positive orthant)-
Ax + b \geq 0\LeftrightarrowAx \geq -b = \text{rhs}. BecomesL_constraint(A, ">=", rhs). L-(negative orthant)-
Ax + b \leq 0\LeftrightarrowAx \leq \text{rhs}. BecomesL_constraint(A, "<=", rhs). L=(fixpoint zero)-
Ax + b = 0. BecomesL_constraint(A, "==", rhs). F(free)No constraint; chunk is silently skipped.
Q,EXP,EXP*,SVECPSD, geometric/power conesBecome
C_constraint(A, K, rhs).QRIs bridged to a standard SOC by the linear map
(u, v, w) \mapsto (u+v, u-v, \sqrt{2} w), then becomesC_constraint(T A, K_soc, T rhs).
Build ROI conic constraints induced by nonlinear VAR cone chunks.
Description
Build ROI conic constraints induced by nonlinear VAR cone chunks.
Usage
.cbf_build_var_constraints(var_chunks, n_total_var, powcones, powstar_cones)
Map a CBF cone name to the corresponding ROI cone object.
Description
Non-parametric cones are recognised by their literal CBF names (Appendix A
of the spec). Parametric cones carry an @k:POW or @k:POW*
prefix that encodes their 0-based index into the respective lookup table.
Usage
.cbf_cone_to_roi(cone_nm, count, powcones, powstar_cones)
Arguments
cone_nm |
CBF cone token. |
count |
Cone dimension (= number of constraint rows). |
powcones |
Primal power-cone lookup table (list of alpha vectors). |
powstar_cones |
Dual power-cone lookup table. |
Value
ROI cone object, or NULL if unsupported.
ROI cone constructors used
ROI::K_soc(n)Second-order (Lorentz) cone, CBF
Q.CBF QRBridged to
ROI::K_soc(n)via(u, v, w) \mapsto (u+v, u-v, \sqrt{2} w)before cone creation.ROI::K_expp(n)Primal exponential cone, CBF
EXP(n=3).ROI::K_expd(n)Dual exp cone, CBF
EXP*.ROI::K_psd(d)PSD cone;
d= svec dimensiond = n(n+1)/2wherenis the matrix side-length.ROI::K_powp(alpha)Power cone;
alphafrom POWCONES or equal-weight vector forGMEANABS/GMEAN.ROI::K_powd(alpha)Dual power cone.
Normalise (row, col) to lower-triangular form by swapping if row < col.
Description
CBF stores only one triangular half of each symmetric matrix. Most files use the lower triangle (row >= col), but some generators emit upper- triangular entries. Since X[row,col] = X[col,row] for symmetric X, we always swap to the lower triangle before computing the svec index.
Usage
.cbf_lt(row, col)
Arguments
row |
0-based row index. |
col |
0-based column index. |
Value
Integer vector c(row, col) with row >= col.
Resolve a parametric power-cone CBF name to an ROI cone object.
Description
Resolve a parametric power-cone CBF name to an ROI cone object.
Usage
.cbf_parametric_cone(cone_nm, count, powcones, powstar_cones)
Parse a POWCONES (or POW*CONES) section into a list of alpha vectors.
Description
CBF spec Sec.B, POWCONES:
POWCONES N total_params <- header: N cones, total combined param count k_0 <- chunk header: length of alpha for cone 0 alpha_0[0] <- chunk body ... alpha_0[k_0 - 1] k_1 <- chunk header for cone 1 ...
Usage
.cbf_parse_powcones(body)
Arguments
body |
Character vector (body lines including the header line), or
|
Value
A list of numeric vectors (one per parametric cone).
Split stripped CBF lines into a named list of section bodies.
Description
Each element of the returned list is named by the CBF keyword. Its value is a character vector of body lines (the spec's HEADER line, if present, is always the first element of that vector).
Usage
.cbf_parse_sections(lines)
Arguments
lines |
Character vector – trimmed, comment-stripped CBF lines. |
Details
Empty lines between sections and leading/trailing blanks have already been removed by the caller; we drop any remaining blank lines here.
Value
Named list of character vectors.
Split a whitespace-delimited string into an integer vector.
Description
Split a whitespace-delimited string into an integer vector.
Usage
.cbf_split_ints(x)
0-based svec index for lower-triangular (row, col) in an m x m matrix.
Description
CBF uses column-major lower-triangular ordering: column 0: rows 0..(m-1), column 1: rows 1..(m-1), ... The cumulative start of column j is j*m - j*(j-1)/2, so k(row, col) = col*m - col*(col-1)/2 + (row - col).
Usage
.cbf_svec_idx(row, col, m)
Arguments
row |
0-based row index (>= col). |
col |
0-based column index. |
m |
Matrix side-length. |
Value
0-based svec index.
svec scaling factor: 1 on the diagonal, sqrt(2) off-diagonal.
Description
svec(X)[k] = X[i,i] for diagonal entries and sqrt(2)*X[i,j] for off-diagonal entries, ensuring <X,Y> = svec(X)' svec(Y).
Usage
.cbf_svec_scale(row, col)
Arguments
row |
0-based row index. |
col |
0-based column index. |
Value
Numeric scalar: 1 or sqrt(2).
Recover PSD matrix side-length from svec dimension.
Description
svec(PSD^{nxn}) has dimension d = n(n+1)/2.
Solving for n gives n = (-1 + \sqrt{1 + 8d}) / 2.
Usage
.cbf_svec_to_n(d)
Arguments
d |
Integer. Svec (half-vectorisation) dimension. |
Value
Integer matrix side-length n.
Collect unique power-cone alpha vectors from all C_constraints.
Description
Returns a list with elements primal and dual, each a list
of unique alpha vectors for K_powcone / K_powconedual cones
found in the constraints. Must be called before flattening constraints
so that the POWCONES section can be written first.
Usage
.collect_powcone_tables(cons)
Arguments
cons |
ROI constraint object (possibly |
Value
List with primal and dual sub-lists.
Recursively collect all combined cone objects from a constraint object.
Description
Returns a list of "cone" objects (one per C_constraint).
Usage
.extract_all_cones(cons)
Find a power-cone alpha vector in a lookup table, returning its 0-based index.
Description
The table is guaranteed to contain alpha because
.collect_powcone_tables() pre-scans all constraints.
Usage
.find_powcone_idx(alpha, table)
Arguments
alpha |
Numeric vector. Parameter vector to look up. |
table |
List of numeric vectors (the lookup table). |
Value
Integer. 0-based position in table.
Flatten one C_constraint into per-cone chunk records.
Description
Groups rows by cone ID from the combined cones object (ROI >= 1.0
new-style API: class "cone" with $cone, $id, and
optional $params fields). Each unique ID produces one chunk.
Usage
.flatten_C_constraint(cc, pow_tables)
Arguments
cc |
A |
pow_tables |
Power-cone lookup tables (see |
Value
List of chunk records.
Flatten one L_constraint into groups of consecutive equal-direction rows.
Description
Consecutive rows with the same inequality direction (">=", etc.)
are grouped into a single chunk so they map to one CBF cone line.
Usage
.flatten_L_constraint(lc)
Arguments
lc |
An |
Value
List of chunk records (see .flatten_all_constraints).
Flatten all ROI constraints into a uniform list of chunk records.
Description
Supports L_constraint and C_constraint objects, including
combined constraint objects created by c() or rbind().
Usage
.flatten_all_constraints(cons, pow_tables)
Arguments
cons |
ROI constraint object (possibly compound). |
pow_tables |
List with |
Details
Each element of the returned list is a list with fields:
- L
slam
simple_triplet_matrix: the constraint-matrix rows for this chunk.- rhs
Numeric vector: ROI right-hand sides for this chunk.
- cbf_cone
Character: CBF cone token (e.g.
"Q.").
Value
List of chunk records.
Get the total number of rows for a combined ROI cone object.
Description
ROI >= 1.0 uses a unified "cone" class; the dimension is simply
the length of the $cone (or $id) vector.
Usage
.roi_cone_dim(cone)
Arguments
cone |
A |
Value
Integer total row count.
Map an ROI cone type code to its CBF cone name string.
Description
ROI >= 1.0 uses a unified "cone" class; the cone family is encoded
as an integer in $cone. Known codes:
1 = K_zero, 2 = K_lin (nonneg), 3 = K_soc, 4 = K_psd,
5 = K_expp, 6 = K_expd, 7 = K_powp, 8 = K_powd.
Usage
.roi_cone_type_to_cbf_name(type_code, params, pow_tables)
Arguments
type_code |
Integer cone type code from |
params |
Named list of parameters for this cone group
( |
pow_tables |
List with |
Details
For power cones (type_code 7 or 8), params$a holds the
scalar alpha; the 0-based POWCONES table index is embedded in the name.
Value
Character CBF cone token.
Write the ACOORD section: sparse constraint-matrix triplets.
Description
Format (CBF spec Sec.B):
ACOORD count con_idx var_idx value <- 0-based indices ...
Usage
.write_acoord_section(flat, con)
Arguments
flat |
List of chunk records. |
con |
Writable file connection. |
Details
CBF uses 0-based indices; R/slam uses 1-based. The offset -1 is
applied to both the constraint-row index and the variable-column index.
Write the BCOORD section: constraint constant-term offsets.
Description
Format (CBF spec Sec.B):
BCOORD count con_idx value <- 0-based con_idx; only non-zero values are written
Usage
.write_bcoord_section(flat, con)
Arguments
flat |
List of chunk records. |
con |
Writable file connection. |
Sign convention
CBF encodes Ax + b \in K. ROI stores Ax - \text{rhs} \in K.
Therefore b_{\text{CBF}} = -\text{rhs}_{\text{ROI}}.
Internal CBF writer – writes all sections to an open connection.
Description
Internal CBF writer – writes all sections to an open connection.
Usage
.write_cbf_internal(op, con)
Write the CON section: cone chunk declarations for scalar constraints.
Description
Format (CBF spec Sec.B):
CON n_total n_chunks CONE_NAME_1 count_1 ... CONE_NAME_k count_k
Usage
.write_con_section(flat, con)
Arguments
flat |
List of chunk records from |
con |
Writable file connection. |
Write OBJACOORD and OBJBCOORD sections for the objective.
Description
- OBJACOORD
Sparse linear-objective coefficients. Format:
OBJACOORD count var_idx value <- 0-based var_idx; only non-zero values
- OBJBCOORD
Scalar constant. Format:
OBJBCOORD value
Usage
.write_objective_sections(obj, con)
Arguments
obj |
ROI objective object (must be an |
con |
Writable file connection. |
Write a POWCONES or POW*CONES section.
Description
Format (CBF spec Sec.B):
POWCONES N total_params k_0 <- chunk header: length of alpha for cone 0 alpha_0[0] ... alpha_0[k_0-1] k_1 ...
Usage
.write_powcones_section(table, keyword, con)
Arguments
table |
List of numeric alpha vectors (one per parametric cone). |
keyword |
Character. Either |
con |
Writable file connection. |
Download a Single CBLIB Instance
Description
Downloads one Conic Benchmark Library (CBLIB) instance from
https://cblib.zib.de/download/all/. By default the downloaded
.cbf.gz archive is decompressed to a plain .cbf file so it can
be passed directly to read_cbf.
Usage
cblib_download(
instance,
folder,
url = .CBLIB_DOWNLOAD_ALL_URL,
method = NULL,
quiet = TRUE,
force = FALSE,
decompress = TRUE,
keep_archive = !decompress
)
Arguments
instance |
Character string. Instance name, with or without the
|
folder |
Character string. Destination directory. No default is
provided; pass |
url |
Character string. Base URL of the CBLIB download directory. |
method |
Character string passed to |
quiet |
Logical. Passed to |
force |
Logical. If |
decompress |
Logical. If |
keep_archive |
Logical. If |
Details
The CBLIB download README documents recursive retrieval from
https://cblib.zib.de/download/all/ using wget. These helpers
provide the same access pattern directly from R for either a single instance
or the full directory listing.
Value
Character string giving the path to the downloaded file. Returns the
decompressed .cbf path when decompress = TRUE, otherwise the
.cbf.gz archive path.
See Also
Examples
## Not run:
path <- cblib_download("CBLIB_ex1.cbf", folder = tempdir())
op <- read_cbf(path)
## End(Not run)
Download All CBLIB Instances
Description
Downloads every .cbf.gz instance linked from
https://cblib.zib.de/download/all/. By default each archive is
decompressed to a plain .cbf file.
Usage
cblib_download_all(
folder,
url = .CBLIB_DOWNLOAD_ALL_URL,
method = NULL,
quiet = TRUE,
force = FALSE,
decompress = TRUE,
keep_archive = !decompress
)
Arguments
folder |
Character string. Destination directory. No default is
provided; pass |
url |
Character string. Base URL of the CBLIB download directory. |
method |
Character string passed to |
quiet |
Logical. Passed to |
force |
Logical. If |
decompress |
Logical. If |
keep_archive |
Logical. If |
Details
According to the CBLIB download README, the complete collection is intended
to be fetched recursively from the /download/all/ directory. This
helper first reads that index page, extracts all linked .cbf.gz
archives, and then downloads them one by one.
Value
Named character vector of downloaded file paths. Names are the CBLIB archive file names listed in the index.
See Also
Examples
## Not run:
paths <- cblib_download_all(folder = tempdir())
## End(Not run)
Read a CBF File into an ROI OP Object
Description
Parses a file in the Conic Benchmark Format (CBF version 4) and returns
an OP object representing the encoded optimisation
problem.
Usage
read_cbf(file, ...)
Arguments
file |
Character string. Path to the |
... |
Currently unused; present for API compatibility with
|
Details
Supported features
Scalar variables (VAR) with all non-parametric cone types:
F.,L+.,L-.,L=..Scalar constraints (CON) with non-parametric cones:
F.,L+.,L-.,L=.,Q.,QR.,EXP.,EXP*.,SVECPSD.,GMEANABS.,GMEAN., and their duals.Parametric power cones (
@k:POW,@k:POW*) via the POWCONES / POW*CONES lookup tables.Integer variables (INT section).
Linear and constant-term objective coefficients (OBJACOORD, OBJBCOORD).
PSD matrix variables (PSDVAR): each
X_j \in S^{m_j \times m_j}_+is replaced by its symmetric vectorisation\mathrm{svec}(X_j)of lengthm_j(m_j+1)/2, appended as new scalar variables with aK_psd(m_j)cone constraint. Objective terms (OBJFCOORD) and scalar-constraint terms (FCOORD) involvingX_jare translated via the trace inner-product identity\langle F, X_j \rangle = \mathrm{svec}(F)^{\top} \mathrm{svec}(X_j).PSD matrix constraints (PSDCON / HCOORD / DCOORD): each LMI
G_p = \sum_j x_j H_{pj} + D_p \in \mathrm{PSD}^{m_p}is svec-ed into aK_psd(m_p)conic constraint on the scalar variables.
Features that generate a warning() and are skipped:
Hotstart sequences (CHANGE keyword is treated as end-of-file).
Value
An OP object.
See Also
Examples
## round-trip: write a simple LP to a temp file then read it back
op <- ROI::OP(
objective = ROI::L_objective(c(1, 1)),
constraints = ROI::L_constraint(
L = matrix(c(1, 1), nrow = 1),
dir = ">=",
rhs = 1
)
)
f <- file.path(tempdir(), "lp.cbf")
write_cbf(op, f)
op2 <- read_cbf(f)
Write an ROI OP Object to a CBF File
Description
Serialises an ROI optimisation problem (OP) to a plain
text file in the Conic Benchmark Format (CBF version 4).
Usage
write_cbf(op, file, ...)
Arguments
op |
An |
file |
Character string. Destination file path (created or overwritten). |
... |
Currently unused. |
Details
Supported OP components
- Objective
Linear objectives (
L_objective).- Constraints
-
L_constraintis mapped toL+(">="),L-("<="), orL=("==") cone chunks.
C_constraintis mapped to the appropriate CBF cone:Q.,QR.,EXP.,EXP*.,SVECPSD.,@k:POW, or@k:POW*. - Variable types
Continuous (
"C") and integer ("I"); integer variables are written to the INT section.- Variable bounds
Translated to
F.,L+,L-, orL=chunks in the VAR section.
Value
Invisibly returns file.
Sign convention
CBF encodes Ax + b \in K. ROI stores Ax - \text{rhs} \in K.
Therefore write_cbf writes b = -\text{rhs} to the BCOORD
section.
See Also
Examples
## simple LP: min x + y s.t. x + y >= 1, x, y >= 0
op <- ROI::OP(
objective = ROI::L_objective(c(1, 1)),
constraints = ROI::L_constraint(
L = matrix(c(1, 1), nrow = 1),
dir = ">=",
rhs = 1
)
)
f <- file.path(tempdir(), "lp.cbf")
write_cbf(op, f)