Getting Started

This section describes the process of creating a Clarabel.rs model, populating its settings and problem data, solving the problem and obtaining and understanding results. It is assumed here that you are building your project using cargo.

Full documentation for the Rust API is available in the API Reference.

The first step is to make the Clarabel solver a dependency in your project by adding:

[dependencies]
clarabel = {version = "0"}

to your project's Cargo.toml file. Then bring the solver into scope in your source files:

use clarabel::algebra::*;
use clarabel::solver::*;

The algebra module defines the CscMatrix type for defining matrices in compressed sparse column format. It also contains some basic utilities for creating and manipulating sparse matrices.

Problem Format

Clarabel solves optimisation problems in the format:

\[\begin{array}{ll} \text{minimize} & \textstyle{\frac{1}{2}}x^\top Px + q^\top x\\ \text{subject to} & Ax + s = b \\ & s \in \mathcal{K}, \end{array}\]

with decision variables $x \in \mathbb{R}^n$, $s \in \mathbb{R}^m$ and data matrices $P=P^\top \succeq 0$, $q \in \mathbb{R}^n$, $A \in \mathbb{R}^{m \times n}$, and $b \in \mathbb{R}^m$. The convex cone $\mathcal{K}$ is a composition of smaller convex cones $\mathcal{K} = \mathcal{K}_1 \times \mathcal{K}_2 \dots \mathcal{K}_p$. Equality conditions can be modelled in this format using the solver's ZeroCone type.

To initialize the solver with an optimisation problem we require three things:

  • The objective function, i.e. the matrix P and the vector q in $\frac{1}{2}x^\top P x + q^\top x$.
  • The data matrix A and vector b, along with a description of the composite cone $\mathcal{K}$ and the dimensions of its constituent pieces.
  • A settings object that specifies how Clarabel solves the problem.

Objective Function

To set the objective function of your optimisation problem simply define the square positive semidefinite matrix $P \in \mathrm{R}^{n\times n}$ and the vector $q \in \mathrm{R}^{n}$.

Suppose that we have a problem with decision variable $x \in \mathbb{R}^3$ and our objective function is:

\[\begin{equation*} \min ~ \frac{1}{2} \left[ \begin{array}{l} x_1 \\ x_2 \\x_3 \end{array} \right] ^T \left[ \begin{array}{rrr} 3.0 & 1.0 & -1.0 \\ 1.0 & 4.0 & 2.0 \\ -1.0 & 2.0 & 5.0 \end{array} \right] \left[ \begin{array}{l} x_1 \\ x_2 \\x_3 \end{array} \right] + \left[ \begin{array}{r} 1 \\ 2 \\-3 \end{array} \right]^T \left[ \begin{array}{l} x_1 \\ x_2 \\x_3 \end{array} \right] \end{equation*}\]

Clarabel expects the P matrix to be supplied in Compressed Sparse Column format. P is assumed by the solver to be symmetric and only values in the upper triangular part of P are needed by the solver, i.e. you only need to provide

\[\begin{equation*} P = \left[ \begin{array}{rrr} 3.0 & 1.0 & -1.0 \\ ⋅ & 4.0 & 2.0 \\ ⋅ & ⋅ & 5.0 \end{array} \right] \end{equation*}\]

The Clarabel default implementation in Rust expects matrix data as a CscMatrix object and provides a set of basic utilities for sparse matrix construction. We can define our cost data as

let P = CscMatrix::new(
    3,                             // m
    3,                             // n
    vec![0, 1, 3, 6],              // colptr
    vec![0, 0, 1, 0, 1, 2],        // rowval
    vec![3., 1., 4., -1., 2., 5.], // nzval
);

let q = vec![1., 2., -3.];
Tip

To specify P = I, you can use

let P = CscMatrix::identity(2);

where in this case we have had to be specific about the floating point data type we want. To use a zero matrix (e.g. if solving an LP), you can use

let P = CscMatrix::spalloc((2,2),0);

to construct a sparse matrix with no entries.

The solver will not conduct any check on the internal correctness of matrices passed in CscMatrix format. You can do this externally using the check_format method, e.g.:

assert!(P.check_format().is_ok());

Constraints

The Clarabel interface expects constraints to be presented in the single vectorized form $Ax + s = b, s \in \mathcal{K}$, where $\mathcal{K} = \mathcal{K}_1 \times \dots \times \mathcal{K}_p$ and each $\mathcal{K}_i$ is one of the solver's supported cone types.

Suppose that we have a problem with decision variable $x \in \mathbb{R}^3$ and our constraints are:

  • A single equality constraint $x_1 + x_2 - x_3 = 1$.
  • A pair of inequalities such that $x_2$ and $x_3$ are each less than 2.
  • A second order cone constraint on the 3-dimensional vector $x$.

For the three constraints above, we have

\[ \begin{align*} A_{eq} &= \left[ \begin{array}{lll} 1 & 1 & -1 \end{array} \right], \quad & b_{eq} &= \left[ \begin{array}{l} 1 \end{array} \right], \\[4ex] A_{ineq} &= \left[ \begin{array}{lll} 0 & 1 & 0 \\ 0 & 0 & 1 \end{array} \right], \quad & b_{ineq} &= \left[ \begin{array}{l} 2\\2 \end{array} \right], \\[4ex] A_{soc} &= \left[ \begin{array}{rrr} -1 & 0 & 0 \\ 0 & -1 & 0 \\ 0 & 0 & -1 \end{array} \right], \quad & b_{soc} &= \left[ \begin{array}{l} 0 \\0 \\0 \end{array} \right] \end{align*} \]

We can define our constraint data as

let Aeq = CscMatrix::new(
        1,                 // m
        3,                 // n
        vec![0, 1, 2, 3],  // colptr
        vec![0, 0, 0],     // rowval
        vec![1., 1., -1.], // nzval
    );

    let Aineq = CscMatrix::new(
        2,                // m
        3,                // n
        vec![0, 0, 1, 2], // colptr
        vec![0, 1],       // rowval
        vec![1., 1.],     // nzval
    );

    let mut Asoc = CscMatrix::identity(3);
    Asoc.negate();

    let A = CscMatrix::vcat(&Aeq, &Aineq);
    let A = CscMatrix::vcat(&A, &Asoc);

    let b = vec![1., 2., 2., 0., 0., 0.];

    // optional correctness check 
    assert!(A.check_format().is_ok());

Clarabel.rs expects to receive a vector of cone specifications. For the above constraints we should also define

# Clarabel.jl cone specification
let cones = [ZeroConeT(1), NonnegativeConeT(2), SecondOrderConeT(3)];

There is no restriction on the ordering of the cones that appear in cones, nor on the number of instances of each type that appear. Your input vector b should be compatible with the sum of the cone dimensions.

Note

Note carefully the signs in the above example. The inequality condition is $A_{ineq} x \le b_{ineq}$, which is equivalent to $A_{ineq} x + s = b_{ineq}$ with $s \ge 0$, i.e. $s$ in the Nonnegative cone. The SOC condition is $x \in \mathcal{K}_{SOC}$, or equivalently $-x + s = 0$ with $s \in \mathcal{K}_{SOC}$.

Solver Settings

Solver settings for the Clarabel's default implementation in Rust are stored in a DefaultSettings object and can be modified by the user. To create a settings object using all defaults you can call the constructor directly:

let settings = DefaultSettings::default();

Alternatively, you can use the DefaultSettingsBuilder to specify custom settings. For example, if you want to disable verbose printing and set a 5 second time limit on the solver, you can use:

let settings = DefaultSettingsBuilder::default()
    .verbose(false)
    .time_limit(1.)
    .build()
    .unwrap();

The full set of user configurable solver settings are listed in the Rust API Reference.

Making a Solver

Finally populate the solver with problem data and solve:

let mut solver = DefaultSolver::new(&P, &q, &A, &b, &cones, settings);

solver.solve();

Results

Once the solver algorithm terminates you can inspect the solution using the solution field of the solver.

println!("Solution status = {:?}", solver.solution.status);
println!("Primal solution = {:?}", solver.solution.x);
println!("Dual solution   = {:?}", solver.solution.z);
println!("Primal slacks   = {:?}", solver.solution.s);

The outcome of the solve is specified in solver.solution.status and will be one of the following :

Status CodeDescription
UnsolvedDefault value, only occurs prior to calling solve
SolvedSolution found
PrimalInfeasibleProblem is primal infeasible
DualInfeasibleProblem is dual infeasible
AlmostSolvedSolution found (reduced accuracy)
AlmostPrimalInfeasibleProblem is primal infeasible (reduced accuracy)
AlmostDualInfeasibleProblem is dual infeasible (reduced accuracy)
MaxIterationsSolver halted after reaching iteration limit
MaxTimeSolver halted after reaching time limit
NumericalErrorSolver terminated with a numerical error
InsufficientProgressSolver terminated due to lack of progress

The total solution time is available in solver.solution.solve_time.

Warning

Be careful to retrieve solver solutions from the solution that is returned by the solver, or directly from a solver object from the solver.solution field. Do not use the solver.variables, since these have both homogenization and equilibration scaling applied and therefore do not solve the optimization problem posed to the solver.