Clarabel Solver Examples

library(clarabel)

Introduction

These examples are from the original Clarabel documentation.

1. Basic Quadratic Program Example

Suppose that we want to solve the following 2-dimensional quadratic programming problem:

\[ \begin{array}{ll} \text{minimize} & 3x_1^2 + 2x_2^2 - x_1 - 4x_2\\ \text{subject to} & -1 \leq x \leq 1, ~ x_1 = 2x_2 \end{array} \]

We will show how to solve this problem using Clarabel in R.

The first step is to put the problem data into the standard form expected by the solver.

1.1. Objective function

The Clarabel solver’s default configuration expects problem data in the form \(\frac{1}{2}x^\top P x + q^\top x\).
We therefore define the objective function data as

\[ P = 2 \cdot \begin{bmatrix} 3 & 0 \\ 0 & 2\end{bmatrix} \mbox{ and } q = \begin{bmatrix} -1 \\ -4\end{bmatrix}. \]

1.2. Constraints

The solver’s default configuration expects constraints in the form \(Ax + s = b\), where \(s \in \mathcal{K}\) for some composite cone \(\mathcal{K}\). We have 1 equality constraint and 4 inequalities, so we require the first element of \(s\) to be zero (i.e. the first constraint will correspond to the equality) and all other elements \(s_i \ge 0\). Our cone constraint on \(s\) is therefore

\[ s \in \mathcal K = \{0\}^1 \times \mathbb{R}^4_{\ge 0}. \]

Define the constraint data as

\[ A = \begin{bmatrix} 1 & -2 \\ 1 & 0 \\ 0 & 1 \\ -1 & 0 \\ 0 & -1\end{bmatrix} \mbox{ and } b=\begin{bmatrix} 0 \\ 1 \\ 1 \\ 1 \\ 1 \end{bmatrix}. \]

Note that Clarabel expects inputs in Compressed Sparse Column (CSC) format for both \(P\) and \(A\) and will try to convert them if not so.

1.3. Solution

P <- Matrix::Matrix(2 * c(3, 0, 0, 2), nrow = 2, ncol = 2, sparse = TRUE)
P <- as(P, "symmetricMatrix")  # P needs to be a symmetric matrix
q <- c(-1, -4)
A <- Matrix::Matrix(c(1, 1, 0, -1, 0, -2, 0, 1, 0, -1), ncol = 2, sparse = TRUE)
b <- c(0, 1, 1, 1, 1)
cones <- c(z = 1L, l = 4L)  ## 1 equality and 4 inequalities, in order
s <- clarabel(A = A, b = b, q = q, P = P, cones = cones)
cat(sprintf("Solution status, description: = (%d, %s)\n",
            s$status, solver_status_descriptions()[s$status]))
#> Solution status, description: = (2, Solver terminated with a solution.)
cat(sprintf("Solution: (x1, x2) = (%f, %f)\n", s$x[1], s$x[2]))
#> Solution: (x1, x2) = (0.428571, 0.214286)

2. Basic SOCP Example

We want to solve the following 2-dimensional optimization problem:

\[ \begin{array}{ll} \text{minimize} & x_2^2\\[2ex] \text{subject to} & \left\|\begin{pmatrix} 2x_1 \\ x_2 \end{pmatrix} - \begin{pmatrix} 2 \\ 2 \end{pmatrix}\right\|_2 \le 1 \end{array} \]

2.1. Objective function

The Clarabel solver’s default configuration expects problem data in the form \(\frac{1}{2}x^\top P x + q^\top x\).
We therefore define the objective function data as

\[ P = 2 \cdot \begin{bmatrix} 0 & 0 \\ 0 & 1\end{bmatrix} \mbox{ and } q = \begin{bmatrix} 0 \\ 0\end{bmatrix}. \]

2.2. Constraints

The solver’s default configuration expects constraints in the form \(Ax + s = b\), where \(s \in \mathcal{K}\) for some composite cone \(\mathcal{K}\). We have a single constraint on the 2-norm of a vector, so we rewrite

\[ \left\|\begin{pmatrix} 2x_1 \\ x_2 \end{pmatrix} - \begin{pmatrix} 2 \\ 2 \end{pmatrix}\right\|_2 \le 1 \quad \Longleftrightarrow \quad \begin{pmatrix} 1 \\ 2x_1 - 2\\ x_2 - 2 \end{pmatrix} \in \mathcal{K}_{SOC} \] which puts our constraint in the form \(b - Ax \in \mathcal{K}_{SOC}\).

2.3. Solution

P <- Matrix::Matrix(2 * c(0, 0, 0, 1), nrow = 2, ncol = 2, sparse = TRUE)
P <- as(P, "symmetricMatrix") # P needs to be a symmetric matrix
q <- c(0, 0)
A <- Matrix::Matrix(c(0, -2.0, 0, 0, 0, 1.0), nrow = 3, ncol = 2, sparse = TRUE)
b <- c(1, -2, -2)
cones <- c(q = 3L)
s <- clarabel(A = A, b = b, q = q, P = P, cones = cones)
cat(sprintf("Solution status, description: = (%d, %s)\n",
            s$status, solver_status_descriptions()[s$status]))
#> Solution status, description: = (2, Solver terminated with a solution.)
cat(sprintf("Solution (x1, x2) = (%f, %f)\n", s$x[1], s$x[2]))
#> Solution (x1, x2) = (1.000000, -1.000000)

3. Control parameters

Clarabel has a number of parameters that control its behavior, including verbosity, time limits, and tolerances; see help on clarabel_control(). As an example, in the last problem, we can reduce the number of iterations.

P <- Matrix::Matrix(2 * c(0, 0, 0, 1), nrow = 2, ncol = 2, sparse = TRUE)
P <- as(P, "symmetricMatrix") # P needs to be a symmetric matrix
q <- c(0, 0)
A <- Matrix::Matrix(c(0, -2.0, 0, 0, 0, 1.0), nrow = 3, ncol = 2, sparse = TRUE)
b <- c(1, -2, -2)
cones <- c(q = 3L)
s <- clarabel(A = A, b = b, q = q, P = P, cones = cones,
              control = list(max_iter = 3)) ## Reduced number of iterations
cat(sprintf("Solution status, description: = (%d, %s)\n",
            s$status, solver_status_descriptions()[s$status]))
#> Solution status, description: = (5, Solver terminated with a solution (reduced accuracy))
cat(sprintf("Solution (x1, x2) = (%f, %f)\n", s$x[1], s$x[2]))
#> Solution (x1, x2) = (1.000000, -0.999999)

Note the different status, which should always be checked in code.