Getting Started

This section describes the process of creating a Clarabel model directly in Python, populating its settings and problem data, solving the problem and obtaining and understanding results. The complete example from this page can be found here.

The first step is to bring the Clarabel solver and other required packages into scope in your code using:

import clarabel
import numpy as np
from scipy import sparse

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 Python expects matrix data in Compressed Sparse Column format as produced by scipy. We can define our cost data as

P = sparse.csc_matrix(
        [[ 3., 1., -1.],
         [ 1., 4.,  2.],
         [-1., 2.,  5.]])

P = sparse.triu(P).tocsc()

q = np.array([1., 2., -3.])

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

# equality constraint
Aeq = sparse.csc_matrix([1.,1.,-1])

beq = np.array([1.])

# equality constraint
Aineq = sparse.csc_matrix(
        [[0., 1., 0.],
         [0., 0., 1.]])

bineq = np.array([2.,2.])

# SOC constraint
Asoc = -sparse.identity(3)

bsoc = np.array([0.,0.,0.])

# Clarabel constraint data
A = sparse.vstack([Aeq,Aineq,Asoc]).tocsc()

b = np.concatenate([beq,bineq,bsoc])

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

cones = [clarabel.ZeroConeT(1),
         clarabel.NonnegativeConeT(2),
         clarabel.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 PyDefaultSettings object and can be modified by the user. To create a settings object using all defaults you can call the constructor directly:

settings = clarabel.DefaultSettings()

If you want to disable verbose printing and set a 5 second time limit on the solver, you can then just modify the fields:

settings.verbose = False
settings.time_limit = 5.

The Clarabel Python interface set supports the same options as those listed in the Rust API Reference.

Making a Solver

Finally, populate the solver with problem data and solve:

solver = clarabel.DefaultSolver(P,q,A,b,cones,settings)

solution = solver.solve()
solution.x  # primal solution
solution.z  # dual solution
solution.s  # primal slacks

Results

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 solution.solve_time.

CVXPY Interface

The same problem above can also be modelled in CVXPY and solved via Clarabel.

import cvxpy as cp

P = [[3., 1., -1.],
     [1., 4.,  2.],
     [-1., 2.,  5.]]

q = [1., 2., -3.]

# Create optimization variables
x = cp.Variable(3)

# Create two
constraints = [x[0] + x[1] - x[2] == 1,
               x[1] <= 2,
               x[2] <= 2,
               cp.SOC(x[0], x[1:])]

# Form objective.
obj = cp.Minimize(0.5*cp.quad_form(x, P) + q @ x)

# Form and solve problem.
prob = cp.Problem(obj, constraints)
prob.solve(solver='CLARABEL', verbose=True)