C

Main solver API

The main C API is imported from the header osqp.h. It is divided into the following headers: * Constants in osqp_api_constants.h * Functions in osqp_api_functions.h * Types in osqp_api_types.h

Main solver data types

struct OSQPSolver

Main OSQP solver structure that holds all information.

Public Members

OSQPSettings *settings

Problem settings.

OSQPSolution *solution

Computed solution.

OSQPInfo *info

Solver information.

OSQPWorkspace *work

Internal solver workspace (contents not public)

struct OSQPInfo

Information about the solution process.

Public Members

char status[32]

Status string, e.g. ‘solved’.

OSQPInt status_val

Status as OSQPInt, defined in osqp_api_constants.h.

OSQPInt status_polish

Polishing status: successful (1), unperformed (0), unsuccessful (-1)

OSQPFloat obj_val

Primal objective value.

OSQPFloat dual_obj_val

Dual objective value.

OSQPFloat prim_res

Norm of primal residual.

OSQPFloat dual_res

Norm of dual residual.

OSQPFloat duality_gap

Duality gap (Primal obj - Dual obj)

OSQPInt iter

Number of iterations taken.

OSQPInt rho_updates

Number of rho updates performned.

OSQPFloat rho_estimate

Best rho estimate so far from residuals.

OSQPFloat setup_time

Setup phase time (seconds)

OSQPFloat solve_time

Solve phase time (seconds)

OSQPFloat update_time

Update phase time (seconds)

OSQPFloat polish_time

Polish phase time (seconds)

OSQPFloat run_time

Total solve time (seconds)

OSQPFloat primdual_int

Integral of duality gap over time (Primal-dual integral), requires profiling.

OSQPFloat rel_kkt_error

Relative KKT error.

struct OSQPSolution

Structure to hold the computed solution (if any), and any certificates of infeasibility (if any) found by the solver.

Public Members

OSQPFloat *x

Primal solution.

OSQPFloat *y

Lagrange multiplier associated with \(l \le Ax \le u\).

OSQPFloat *prim_inf_cert

Primal infeasibility certificate.

OSQPFloat *dual_inf_cert

Dual infeasibility certificate.

struct OSQPCodegenDefines

Structure to hold the settings for the generated code

Public Members

OSQPInt embedded_mode

Embedded mode (1 = vector update, 2 = vector + matrix update)

OSQPInt float_type

Use floats if 1, doubles if 0.

OSQPInt printing_enable

Enable printing if 1.

OSQPInt profiling_enable

Enable timing of code sections if 1.

OSQPInt interrupt_enable

Enable interrupt checking if 1.

OSQPInt derivatives_enable

Enable deriatives if 1.

Solver capabilities

enum osqp_capabilities_type

Flags used to represent what capabilities the compiled OSQP solver has.

These flags are used as bit flags, so capabilities can be tested using bit-wise operations.

Values:

enumerator OSQP_CAPABILITY_DIRECT_SOLVER

< A direct linear solver is present in the algebra.

enumerator OSQP_CAPABILITY_INDIRECT_SOLVER

< An indirect linear solver is present in the algebra.

enumerator OSQP_CAPABILITY_CODEGEN

< Code generation is present.

enumerator OSQP_CAPABILITY_UPDATE_MATRICES

< The problem matrices can be updated.

enumerator OSQP_CAPABILITY_DERIVATIVES

< Solution derivatives w.r.t P/q/A/l/u are available.

Solver status

enum osqp_status_type

Values:

enumerator OSQP_SOLVED
enumerator OSQP_SOLVED_INACCURATE
enumerator OSQP_PRIMAL_INFEASIBLE
enumerator OSQP_PRIMAL_INFEASIBLE_INACCURATE
enumerator OSQP_DUAL_INFEASIBLE
enumerator OSQP_DUAL_INFEASIBLE_INACCURATE
enumerator OSQP_MAX_ITER_REACHED
enumerator OSQP_TIME_LIMIT_REACHED
enumerator OSQP_NON_CVX
enumerator OSQP_SIGINT
enumerator OSQP_UNSOLVED

Polish status

enum osqp_polish_status_type

Values:

enumerator OSQP_POLISH_LINSYS_ERROR
enumerator OSQP_POLISH_FAILED
enumerator OSQP_POLISH_NOT_PERFORMED
enumerator OSQP_POLISH_SUCCESS
enumerator OSQP_POLISH_NO_ACTIVE_SET_FOUND

Linear system solvers

enum osqp_linsys_solver_type

Values:

enumerator OSQP_UNKNOWN_SOLVER
enumerator OSQP_DIRECT_SOLVER
enumerator OSQP_INDIRECT_SOLVER

Preconditioners for CG method

enum osqp_precond_type

Values:

enumerator OSQP_NO_PRECONDITIONER
enumerator OSQP_DIAGONAL_PRECONDITIONER

Solver errors

enum osqp_error_type

Values:

enumerator OSQP_NO_ERROR
enumerator OSQP_DATA_VALIDATION_ERROR
enumerator OSQP_SETTINGS_VALIDATION_ERROR
enumerator OSQP_LINSYS_SOLVER_INIT_ERROR
enumerator OSQP_NONCVX_ERROR
enumerator OSQP_MEM_ALLOC_ERROR
enumerator OSQP_WORKSPACE_NOT_INIT_ERROR
enumerator OSQP_ALGEBRA_LOAD_ERROR
enumerator OSQP_FOPEN_ERROR
enumerator OSQP_CODEGEN_DEFINES_ERROR
enumerator OSQP_DATA_NOT_INITIALIZED
enumerator OSQP_FUNC_NOT_IMPLEMENTED

Function not implemented in this library

enumerator OSQP_LAST_ERROR_PLACE

Solver parameters and settings

OSQP_VERBOSE
OSQP_WARM_STARTING
OSQP_SCALING
OSQP_POLISHING

ADMM parameters:

OSQP_RHO
OSQP_SIGMA
OSQP_ALPHA
OSQP_RHO_MIN
OSQP_RHO_MAX
OSQP_RHO_TOL

tolerance for detecting if an inequality is set to equality

OSQP_RHO_EQ_OVER_RHO_INEQ
OSQP_RHO_IS_VEC

CG parameters:

OSQP_CG_MAX_ITER
OSQP_CG_TOL_REDUCTION
OSQP_CG_TOL_FRACTION

Adaptive rho logic .. doxygendefine:: OSQP_ADAPTIVE_RHO

OSQP_ADAPTIVE_RHO_INTERVAL

Default interval for iteration-based rho update.

OSQP_ADAPTIVE_RHO_TOLERANCE

tolerance for adopting new rho; minimum ratio between new rho and the current one

OSQP_ADAPTIVE_RHO_FRACTION

fraction of setup time after which we update rho

OSQP_ADAPTIVE_RHO_MULTIPLE_TERMINATION

Multiple of check_termination after which we update rho if using interval-based rho adaptation and adaptive_rho_interval == 0.

OSQP_ADAPTIVE_RHO_FIXED

Number of iterations after which we update rho if using interval-based rho adaptation and adaptive_rho_interval == 0 and termination_check is disabled.

Termination parameters .. doxygendefine:: OSQP_MAX_ITER

OSQP_EPS_ABS
OSQP_EPS_REL
OSQP_EPS_PRIM_INF
OSQP_EPS_DUAL_INF
OSQP_SCALED_TERMINATION
OSQP_TIME_LIMIT

Disable time limit by default.

OSQP_CHECK_TERMINATION
OSQP_DELTA
OSQP_POLISH_REFINE_ITER

Hard-coded values and settings:

OSQP_NULL
OSQP_NAN
OSQP_INFTY
OSQP_DIVISION_TOL
OSQP_PRINT_INTERVAL
OSQP_MIN_SCALING

minimum scaling value

OSQP_MAX_SCALING

maximum scaling value

OSQP_CG_TOL_MIN
OSQP_CG_POLISH_TOL

Main solver API

OSQPInt osqp_capabilities(void)

Return the capabilities of the OSQP solver.

Returns:

const char *osqp_version(void)

Return OSQP version

Returns:

OSQP version string

const char *osqp_error_message(OSQPInt error_flag)

Return the error string for a given error code.

Parameters:

error_flag – Error code to get description of

Returns:

String describing the error code

void osqp_get_dimensions(OSQPSolver *solver, OSQPInt *m, OSQPInt *n)

Return the number of variables and constraints

Parameters:
  • solver[in] Solver

  • m[out] Pointer to integer that will hold m

  • n[out] Pointer to integer that will hold n

OSQPInt osqp_setup(OSQPSolver **solverp, const OSQPCscMatrix *P, const OSQPFloat *q, const OSQPCscMatrix *A, const OSQPFloat *l, const OSQPFloat *u, OSQPInt m, OSQPInt n, const OSQPSettings *settings)

Initialize OSQP solver allocating memory.

It performs:

  • data and settings validation

  • problem data scaling

  • setup linear system solver:

    • direct solver: KKT matrix factorization is performed here

    • indirect solver: reduced KKT matrix preconditioning is performed here

NB: This is the only function that allocates dynamic memory and is not used during code generation

Parameters:
  • solverp – Solver pointer

  • P – Problem data (upper triangular part of quadratic cost term, csc format)

  • q – Problem data (linear cost term)

  • A – Problem data (constraint matrix, csc format)

  • l – Problem data (constraint lower bound)

  • u – Problem data (constraint upper bound)

  • m – Problem data (number of constraints)

  • n – Problem data (number of variables)

  • settings – Solver settings

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_solve(OSQPSolver *solver)

Solve quadratic program

The final solver information is stored in the solver->info structure

The solution is stored in the solver->solution structure

If the problem is primal infeasible, the certificate is stored in solver->delta_y

If the problem is dual infeasible, the certificate is stored in solver->delta_x

Parameters:

solver – Solver

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_get_solution(OSQPSolver *solver, OSQPSolution *solution)

Store the optimization problem result from solver solver into the solution solution. Note that solution must already be allocated with the component vectors the correct lengths.

Parameters:
  • solver – Solver

  • solution – Solution object to store result in

Returns:

Error flag

OSQPInt osqp_cleanup(OSQPSolver *solver)

Cleanup workspace by deallocating memory

This function is not used in code generation

Parameters:

solver – Solver

Returns:

Exitflag for errors (0 if no errors)

Sublevel API

These functions can be called without performing setup again.

OSQPInt osqp_warm_start(OSQPSolver *solver, const OSQPFloat *x, const OSQPFloat *y)

Warm start primal and dual variables

Parameters:
  • solver – Solver

  • x – Primal variable, NULL if none

  • y – Dual variable, NULL if none

Returns:

Exitflag for errors (0 if no errors)

void osqp_cold_start(OSQPSolver *solver)

Cold start workspace variables xz and y

Parameters:

solver – Solver

void osqp_set_default_settings(OSQPSettings *settings)

Get the default settings from the osqp_api_constants.h file.

Note

the settings parameter must already be allocated in memory.

Parameters:

settings – Settings structure to populate

OSQPInt osqp_update_settings(OSQPSolver *solver, const OSQPSettings *new_settings)

Update settings in solver with the new settings from new_settings.

The following settings can only be set at problem setup time through osqp_setup and are ignored in this function:

  • scaling

  • rho_is_vec

  • sigma

  • adaptive_rho

  • adaptive_rho_interval

  • adaptive_rho_fraction

  • adaptive_rho_tolerance

The rho setting must be updated using osqp_update_rho, and is ignored by this function.

Note

Every setting from new_settings is copied to solver.

Parameters:
  • solver – Solver

  • new_settings – New solver settings

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_update_rho(OSQPSolver *solver, OSQPFloat rho_new)

Update the ADMM parameter rho.

Limit it between OSQP_RHO_MIN and OSQP_RHO_MAX.

Parameters:
  • solver – Solver

  • rho_new – New rho value

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_update_data_vec(OSQPSolver *solver, const OSQPFloat *q_new, const OSQPFloat *l_new, const OSQPFloat *u_new)

Update problem data vectors

Parameters:
  • solver – Solver

  • q_new – New linear cost, NULL if none

  • l_new – New lower bound, NULL if none

  • u_new – New upper bound, NULL if none

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_update_data_mat(OSQPSolver *solver, const OSQPFloat *Px_new, const OSQPInt *Px_new_idx, OSQPInt P_new_n, const OSQPFloat *Ax_new, const OSQPInt *Ax_new_idx, OSQPInt A_new_n)

Update elements of matrices P (upper triangular) and A by preserving their sparsity structures.

If Px_new_idx (Ax_new_idx) is OSQP_NULL, Px_new (Ax_new) is assumed to be as long as P->x (A->x) and the whole P->x (A->x) is replaced.

Parameters:
  • solver – Solver

  • Px_new – Vector of new elements in P->x (upper triangular), NULL if none

  • Px_new_idx – Index mapping new elements to positions in P->x

  • P_new_n – Number of new elements to be changed

  • Ax_new – Vector of new elements in A->x, NULL if none

  • Ax_new_idx – Index mapping new elements to positions in A->x

  • A_new_n – Number of new elements to be changed

Returns:

output flag: 0: OK 1: P_new_n > nnzP 2: A_new_n > nnzA <0: error in the update

Settings structure

The setting structure has the following fields.

struct OSQPSettings

User settings

Public Members

OSQPInt device

device identifier; currently used for CUDA devices

enum osqp_linsys_solver_type linsys_solver

linear system solver to use

OSQPInt allocate_solution

boolean; allocate solution in OSQPSolver during osqp_setup

OSQPInt verbose

boolean; write out progress

OSQPInt profiler_level

integer; level of detail for profiler annotations

OSQPInt warm_starting

boolean; warm start

OSQPInt scaling

data scaling iterations; if 0, then disabled

OSQPInt polishing

boolean; polish ADMM solution

OSQPFloat rho

ADMM penalty parameter.

OSQPInt rho_is_vec

boolean; is rho scalar or vector?

OSQPFloat sigma

ADMM penalty parameter.

OSQPFloat alpha

ADMM relaxation parameter.

OSQPInt cg_max_iter

maximum number of CG iterations per solve

OSQPInt cg_tol_reduction

number of consecutive zero CG iterations before the tolerance gets halved

OSQPFloat cg_tol_fraction

CG tolerance (fraction of ADMM residuals)

osqp_precond_type cg_precond

Preconditioner to use in the CG method.

OSQPInt adaptive_rho

rho stepsize adaption method

OSQPInt adaptive_rho_interval

Interval between rho adaptations

When adaptive_rho == OSQP_ADAPTIVE_RHO_UPDATE_ITERATIONS, this is the number of iterations between rho updates.

Not used when adaptive_rho is any other value.

OSQPFloat adaptive_rho_fraction

Adaptation parameter controlling when non-fixed rho adaptations occur.

  • When adaptive_rho == OSQP_ADAPTIVE_RHO_UPDATE_TIME, this is the fraction of the setup time to use as the rho adaptation interval.

  • When adaptive_rho == OSQP_ADAPTIVE_RHO_UPDATE_KKT_ERROR, this is the fraction of the previous KKT error to adapt rho at.

  • Not used for any other adaptive_rho value.

OSQPFloat adaptive_rho_tolerance

Tolerance applied when adapting rho.

New rho must be X times larger or smaller than the current one to change it

OSQPInt max_iter

maximum number of iterations

OSQPFloat eps_abs

absolute solution tolerance

OSQPFloat eps_rel

relative solution tolerance

OSQPFloat eps_prim_inf

primal infeasibility tolerance

OSQPFloat eps_dual_inf

dual infeasibility tolerance

OSQPInt scaled_termination

boolean; use scaled termination criteria

OSQPInt check_termination

integer, check termination interval; if 0, checking is disabled

OSQPInt check_dualgap

Boolean; use duality gap termination criteria.

OSQPFloat time_limit

maximum time to solve the problem (seconds)

OSQPFloat delta

regularization parameter for polishing

OSQPInt polish_refine_iter

number of iterative refinement steps in polishing

The Settings structure can be created and freed using the following functions.

OSQPSettings *OSQPSettings_new()

Allocate a new OSQPSettings object with the default settings.

Note

Not available in embedded mode (requires malloc)

Returns:

Pointer to new settings object, or null on error

void OSQPSettings_free(OSQPSettings *settings)

Free an OSQPSettings object.

Note

Not available in embedded mode (requires free)

Parameters:

settings – The settings object to free

Compute solution derivatives

Adjoint derivatives of the QP problem can be computed at the current solution.

OSQPInt osqp_adjoint_derivative_compute(OSQPSolver *solver, OSQPFloat *dx, OSQPFloat *dy)

Compute internal data structures needed for calculation of the adjoint derivatives of P/q/A/l/u.

Note

An optimal solution must be obtained before calling this function.

Parameters:
  • solver[in] Solver

  • dx[in] Vector of dx values (observed - true) of length n

  • dy[in] Vector of dy values (observed - true) of length m

Returns:

Exitflag for errors (0 if no errors)

OSQPInt osqp_adjoint_derivative_get_mat(OSQPSolver *solver, OSQPCscMatrix *dP, OSQPCscMatrix *dA)

Calculate adjoint derivatives of P/A.

Note

osqp_adjoint_derivative_compute must be called first.

Parameters:
  • solver[in] Solver

  • dP[out] Matrix of dP values (n x n)

  • dA[out] Matrix of dA values (m x n)

Returns:

Exitflag for errors (0 if no errors; dP, dA are filled in)

OSQPInt osqp_adjoint_derivative_get_vec(OSQPSolver *solver, OSQPFloat *dq, OSQPFloat *dl, OSQPFloat *du)

Calculate adjoint derivatives of q/l/u.

Note

osqp_adjoint_derivative_compute must be called first.

Parameters:
  • solver[in] Solver

  • dq[out] Vector of dq values of length n

  • dl[out] Matrix of dl values of length m

  • du[out] Matrix of du values of length m

Returns:

Exitflag for errors (0 if no errors; dq, dl, du are filled in)

OSQPInt osqp_codegen(OSQPSolver *solver, const char *output_dir, const char *prefix, OSQPCodegenDefines *defines)

Generate source files with a statically allocated OSQPSolver structure.

Note

osqp_setup must be called before a problem can be code generated.

Parameters:
  • solver – Solver

  • output_dir – Path to directory to output the files to. This string must include the trailing directory separator, and an empty string means output to the current directory.

  • prefix – String prefix for the variables and generated files.

  • defines – The defines to use in the generated code.

Returns:

Exitflag for errors (0 if no errors)

void osqp_set_default_codegen_defines(OSQPCodegenDefines *defines)

Set default codegen define values.

Note

The defines structure must already be allocated in memory.

Parameters:

defines – Structure to set to default values.

OSQPCodegenDefines *OSQPCodegenDefines_new()

Allocate a new OSQPCodegenDefines object with the default options.

Note

Not available in embedded mode (requires malloc)

Returns:

Pointer to new codegen defines object, or null on error

void OSQPCodegenDefines_free(OSQPCodegenDefines *defs)

Free an OSQPCodegenDefines object.

Note

Not available in embedded mode (requires free)

Parameters:

defs – The defines object to free

Data types

Primitive types

The most basic used datatypes are

typedef int OSQPInt
typedef double OSQPFloat
  • OSQPInt: can be long or int if the compiler flag OSQP_USE_LONG is set or not

  • OSQPFloat: can be a float or a double if the compiler flag OSQP_USE_FLOAT is set or not.

Matrices

The matrices are defined in Compressed Sparse Column (CSC) format using zero-based indexing, using the OSQPCscMatrix datatype.

struct OSQPCscMatrix

Matrix in compressed-column form. The structure is used internally to store matrices in the triplet form as well, but the API requires that the matrices are in the CSC format.

Public Members

OSQPInt m

number of rows

OSQPInt n

number of columns

OSQPInt *p

column pointers (size n+1); col indices (size nzmax) starting from 0 for triplet format

OSQPInt *i

row indices, size nzmax starting from 0

OSQPFloat *x

numerical values, size nzmax

OSQPInt nzmax

maximum number of entries

OSQPInt nz

number of entries in triplet matrix, -1 for csc

OSQPInt owned

1 if the p, i and x pointers were allocated automatically, 0 if they are owned by the user

Data

As a helper, OSQP provides the OSQPCscMatrix_set_data() function to assign existing data to an existing OSQPCscMatrix.

Note

When using the OSQPCscMatrix_set_data() function, the user is responsible for managing the memory used by the OSQPCscMatrix.x, OSQPCscMatrix.p OSQPCscMatrix.i arrays.

void OSQPCscMatrix_set_data(OSQPCscMatrix *M, OSQPInt m, OSQPInt n, OSQPInt nzmax, OSQPFloat *x, OSQPInt *i, OSQPInt *p)

Populates a Compressed-Column-Sparse matrix from existing data arrays. This will just assign the pointers - no malloc or copying is done.

Parameters:
  • M – Matrix pointer

  • m – First dimension

  • n – Second dimension

  • nzmax – Maximum number of nonzero elements

  • x – Vector of data

  • i – Vector of row indices

  • p – Vector of column pointers

Memory management

In non-embedded versions, the CSC matrix objects can be created with existing data, and free’d

OSQPCscMatrix *OSQPCscMatrix_new(OSQPInt m, OSQPInt n, OSQPInt nzmax, OSQPFloat *x, OSQPInt *i, OSQPInt *p)

Allocates a new Compressed-Column-Sparse (CSC) matrix from existing arrays.

This will malloc the new matrix structure, but use the arrays passed in as the backing data for the matrix (e.g. not copy the actual matrix data, just reference the existing data.)

Note

Not available in embedded mode (requires malloc)

Parameters:
  • m – Number of rows

  • n – number of columns

  • nzmax – Maximum number of nonzero elements

  • x – Vector of data

  • i – Vector of row indices

  • p – Vector of column pointers

Returns:

Pointer to new CSC matrix, or null on error

void OSQPCscMatrix_free(OSQPCscMatrix *mat)

Free a CSC matrix object allocated by OSQPCscMatrix_new.

Note

This function will only free the internal OSQPCscMatrix::x, OSQPCscMatrix::p, OSQPCscMatrix::i arrays if the OSQPCscMatrix::owned variable is set to 1, otherwise only the outer CSC wrapper is free’d.

Note

Not available in embedded mode (requires free)

Parameters:

mat – Matrix to free

Common matrices

In non-embedded versions, several helper functions are provided to create common matrices in the CSC format, including a matrix of all structural zeros, the square identity matrix, and diagonal matrices.

Note

The OSQPCscMatrix.x, OSQPCscMatrix.p and OSQPCscMatrix.i arrays are managed by OSQP when using these functions, meaning they will automatically be freed inside OSQPCscMatrix_free().

OSQPCscMatrix *OSQPCscMatrix_zeros(OSQPInt m, OSQPInt n)

Allocates a new Compressed-Column-Sparse (CSC) matrix with zero entries.

This matrix will contain all structural zeros, so no non-zero elements will be present.

Note

Not available in embedded mode (requires malloc)

Parameters:
  • m – Number of rows

  • n – Number of columns

Returns:

Pointer to new CSC matrix, or null on error

OSQPCscMatrix *OSQPCscMatrix_identity(OSQPInt m)

Allocates a new Compressed-Column-Sparse (CSC) identity with 1s on the diagonal.

Note

Not available in embedded mode (requires malloc)

Parameters:

m – Number of rows/columns

Returns:

Pointer to new CSC matrix, or null on error

The following diagonal matrix generation functions can generate tall, wide and square matrices. They will always start in the upper-left corner and fill the main diagonal until the smallest dimension runs out, leaving the remaining rows/columns as all 0.

OSQPCscMatrix *OSQPCscMatrix_diag_scalar(OSQPInt m, OSQPInt n, OSQPFloat scalar)

Allocates a new Compressed-Column-Sparse (CSC) diagonal matrix with a given value.

This will start populating at the upper-left element and continue down the main diagonal until either the last row or column (in the case of non-square matrices).

Note

Not available in embedded mode (requires malloc)

Parameters:
  • m – Number of rows

  • n – Number of columns

  • scalar – Scalar value to put on the diagonal

Returns:

Pointer to new CSC matrix, or null on error

OSQPCscMatrix *OSQPCscMatrix_diag_vec(OSQPInt m, OSQPInt n, OSQPFloat *vals)

Allocates a new Compressed-Column-Sparse (CSC) diagonal matrix with given values on the diagonal.

This will start populating at the upper-left element and continue down the main diagonal until either the last row or column (in the case of non-square matrices).

Note

Not available in embedded mode (requires malloc)

Parameters:
  • m – Number of rows

  • n – Number of columns

  • vals – Values to put on the diagonal - length min(n,m)

Returns:

Pointer to new CSC matrix, or null on error