File numerics/src/tools/CSparseMatrix.h#

Go to the source code of this file

Structure definition and functions related to sparse matrix storage in Numerics

Defines

CS_INT#
css#
csn#
NSM_UNKNOWN_ERR(func, orig)#
NSM_NROW_CSR(mat)#
NSM_NCOL_CSR(mat)#

Typedefs

typedef struct cs_dl_symbolic cs_dls#
typedef struct cs_di_symbolic cs_dis#
typedef struct cs_dl_numeric cs_dln#
typedef struct cs_di_numeric cs_din#
typedef struct cs_di_sparse CSparseMatrix#

Functions

int CSparseMatrix_lu_factorization(CS_INT order, const CSparseMatrix *A, double tol, CSparseMatrix_factors *cs_lu_A)#

compute a LU factorization of A and store it in a workspace

Parameters:
  • order – control if ordering is used

  • A – the sparse matrix

  • tol – the tolerance

  • cs_lu_A – the parameter structure that eventually holds the factors

Returns:

1 if the factorization was successful, 1 otherwise

int CSparseMatrix_chol_factorization(CS_INT order, const CSparseMatrix *A, CSparseMatrix_factors *cs_chol_A)#

compute a Cholesky factorization of A and store it in a workspace

Parameters:
  • order – control if ordering is used

  • A – the sparse matrix

  • cs_chol_A – the parameter structure that eventually holds the factors

Returns:

1 if the factorization was successful, 1 otherwise

int CSparseMatrix_ldlt_factorization(CS_INT order, const CSparseMatrix *A, CSparseMatrix_factors *cs_ldlt_A)#

compute a LDLT factorization of A and store it in a workspace

Parameters:
  • order – control if ordering is used

  • A – the sparse matrix

  • cs_ldlt_A – the parameter structure that eventually holds the factors

Returns:

1 if the factorization was successful, 1 otherwise

CS_INT CSparseMatrix_solve(CSparseMatrix_factors *cs_lu_A, double *x, double *b)#

reuse a LU factorization (stored in the cs_lu_A) to solve a linear system Ax = b

Parameters:
  • cs_lu_A – contains the LU factors of A, permutation information

  • x – workspace

  • b[inout] on input RHS of the linear system; on output the solution

Returns:

0 if failed, 1 otherwise

CS_INT CSparseMatrix_spsolve(CSparseMatrix_factors *cs_lu_A, CSparseMatrix *X, CSparseMatrix *B)#

reuse a LU factorization (stored in the cs_lu_A) to solve a linear system Ax = B with a sparse r.h.s

Parameters:
  • cs_lu_A – contains the LU factors of A, permutation information

  • X – a csc sparse matrix workspace

  • B[inout] on input sparse RHS of the linear system; on output the solution

Returns:

0 if failed, 1 otherwise

CS_INT CSparseMatrix_chol_solve(CSparseMatrix_factors *cs_chol_A, double *x, double *b)#

reuse a Cholesky factorization (stored in the cs_chol_A) to solve a linear system Ax = b

Parameters:
  • cs_chol_A – contains the Cholesky factors of A, permutation information

  • x – workspace

  • b[inout] on input RHS of the linear system; on output the solution

Returns:

0 if failed, 1 otherwise

CS_INT CSparseMatrix_chol_spsolve(CSparseMatrix_factors *cs_chol_A, CSparseMatrix *X, CSparseMatrix *B)#

reuse a Cholesky factorization (stored in the cs_chol_A) to solve a linear system Ax = B with a sparse r.h.s

Parameters:
  • cs_chol_A – contains the Cholesky factors of A, permutation information

  • X – a csc sparse matrix workspace

  • b[inout] on input sparse RHS of the linear system; on output the solution

Returns:

0 if failed, 1 otherwise

CS_INT CSparseMatrix_ldlt_solve(CSparseMatrix_factors *cs_ldlt_A, double *x, double *b)#

reuse a LDLT factorization (stored in the cs_ldlt_A) to solve a linear system Ax = b

Parameters:
  • cs_ldlt_A – contains the LDLT factors of A, permutation information

  • x – workspace

  • b[inout] on input RHS of the linear system; on output the solution

Returns:

0 if failed, 1 otherwise

void CSparseMatrix_free_lu_factors(CSparseMatrix_factors *cs_lu_A)#

Free a workspace related to a LU factorization.

Parameters:

cs_lu_A – the structure to free

int CSparseMatrix_aaxpby(const double alpha, const CSparseMatrix *A, const double *x, const double beta, double *y)#

Matrix vector multiplication : y = alpha*A*x+beta*y.

Parameters:
  • alpha[in] matrix coefficient

  • A[in] the sparse matrix

  • x[in] pointer on a dense vector of size A->n

  • beta[in] vector coefficient

  • y[inout] pointer on a dense vector of size A->n

Returns:

0 if A x or y is NULL else 1

CSparseMatrix *CSparseMatrix_alloc_for_copy(const CSparseMatrix *const m)#

Allocate a CSparse matrix for future copy (as in NSM_copy)

Parameters:

m – the matrix used as model

Returns:

an newly allocated matrix

CS_INT CSparseMatrix_to_dense(const CSparseMatrix *const A, double *B)#
int CSparseMatrix_print(const CSparseMatrix *A, int brief)#

print a matrix to std output

Parameters:
  • A – matrix to print

  • brief – if positive, print only a portion of the matrix

int CSparseMatrix_print_in_file(const CSparseMatrix *A, int brief, FILE *file)#

print a matrix to a text file

Parameters:
  • A – matrix to print

  • brief – if positive, print only a portion of the matrix

  • file – file descriptor

int CSparseMatrix_print_in_Matlab_file(const CSparseMatrix *A, int brief, FILE *file)#
CSparseMatrix *CSparseMatrix_new_from_file(FILE *file)#
CS_INT CSparseMatrix_zentry(CSparseMatrix *T, CS_INT i, CS_INT j, double x, double threshold)#

Add an entry to a triplet matrix only if the absolute value is greater than threshold.

Parameters:
  • T – the sparse matrix

  • i – row index

  • j – column index

  • x – the value

Returns:

integer value : 1 if the absolute value is less than DBL_EPSILON, otherwise the return value of cs_entry.

CS_INT CSparseMatrix_block_dense_zentry(CSparseMatrix *T, CS_INT row_off, CS_INT col_off, double *x, CS_INT row_size, CS_INT col_size, double threshold)#
CS_INT CSparseMatrix_symmetric_zentry(CSparseMatrix *T, CS_INT i, CS_INT j, double x, double threshold)#

Add an entry to a symmetric triplet matrix only if the absolute value is greater than threshold.

Parameters:
  • T – the sparse matrix

  • i – row index

  • j – column index

  • x – the value

Returns:

integer value : 1 if the absolute value is less than DBL_EPSILON, otherwise the return value of cs_entry.

CS_INT CSparseMatrix_entry(CSparseMatrix *T, CS_INT i, CS_INT j, double x)#

Add an entry to a triplet matrix.

Parameters:
  • T – the sparse matrix

  • i – row index

  • j – column index

  • x – the value

CS_INT CSparseMatrix_symmetric_entry(CSparseMatrix *T, CS_INT i, CS_INT j, double x)#

Add an entry to a symmetric triplet matrix.

Parameters:
  • T – the sparse matrix

  • i – row index

  • j – column index

  • x – the value

int CSparseMatrix_check_triplet(CSparseMatrix *T)#

Check if the given triplet matrix is properly constructed (col and row indices are correct)

Parameters:

T – the sparse matrix to check

Returns:

0 if the matrix is fine, 1 otherwise

int CSparseMatrix_check_csc(CSparseMatrix *T)#

Check if the given triplet matrix is properly constructed (col and row indices are correct)

Parameters:

T – the sparse matrix to check

Returns:

0 if the matrix is fine, 1 otherwise

CSparseMatrix *CSparseMatrix_spfree_on_stack(CSparseMatrix *A)#

Free space allocated for a SparseMatrix.

note : cs_spfree also free the cs_struct this fails when the struct is allocated on the stack.

Parameters:

A – the sparse matrix

Returns:

NULL on success

void CSparseMatrix_copy(const CSparseMatrix *const A, CSparseMatrix *B)#

Copy a CSparseMatrix inside another CSparseMatrix.

Reallocations are performed if B cannot hold a copy of A

Parameters:
  • A[in] a CSparseMatrix

  • B[inout] a CSparseMatrix

int CSparseMatrix_scal(const double alpha, const CSparseMatrix *A)#

Multiply a matrix with a double alpha*A –> A.

Parameters:
  • alpha – the coefficient

  • A – the matrix

double CSparseMatrix_get_value(const CSparseMatrix *A, CS_INT i, CS_INT j)#

Return the element A(i,j)

Parameters:
  • A – the sparse matrix

  • i – the row index

  • j – the column index

void CSparseMatrix_write_in_file_python(const CSparseMatrix *const m, FILE *file)#

print a matrix to a text file in pyhton format

Parameters:
  • m – matrix to print

  • file – file descriptor

int CSparseMatrix_max_by_columns(const CSparseMatrix *A, double *max)#

Compute the max by columns of a sparse matrix.

Parameters:
  • A – the sparse matrix

  • max – the vector of maximum by columns

  • j – the column index

int CSparseMatrix_max_abs_by_columns(const CSparseMatrix *A, double *max)#

Compute the max in absolute value by columns of a sparse matrix.

Parameters:
  • A – the sparse matrix

  • max – the vector of maximum by columns

  • j – the column index

struct CSparseMatrix_factors#
#include <CSparseMatrix.h>

Information used and produced by CSparse for an LU factorization.

Public Members

CS_INT n#

size of linear system

css *S#

symbolic analysis

csn *N#

numerics factorization