My Project
programmer's documentation
Data Structures | Macros | Typedefs | Functions
cs_sdm.h File Reference
#include <stdio.h>
#include <string.h>
#include <math.h>
Include dependency graph for cs_sdm.h:

Go to the source code of this file.

Data Structures

struct  cs_sdm_block_t
 
struct  cs_sdm_t
 

Macros

#define CS_SDM_BY_BLOCK   (1 << 0) /* Matrix is defined by block */
 
#define CS_SDM_SYMMETRIC   (1 << 1) /* Matrix is symmetric by construction */
 
#define CS_SDM_SHARED_VAL   (1 << 2) /* Matrix is not owner of its values */
 

Typedefs

typedef void() cs_sdm_product_t(const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Generic prototype for computing a local dense matrix-product c = a*b where c has been previously allocated. More...
 
typedef void() cs_sdm_matvec_t(const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Generic prototype for computing a dense matrix-vector product mv has been previously allocated. More...
 

Functions

static cs_real_t cs_sdm_dot (int n, const cs_real_t x[], const cs_real_t y[])
 Basic dot product for a small vector For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas. More...
 
static void cs_sdm_scalvect (int n, const cs_real_t s, const cs_real_t x[], cs_real_t y[])
 Multiply a small vector by a scalar coefficient: y = s x For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas. More...
 
static void cs_sdm_add_scalvect (int n, const cs_real_t s, const cs_real_t x[], cs_real_t y[])
 Multiply a small vector by a scalar coefficient: y += s x For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas. More...
 
cs_sdm_t * cs_sdm_create (cs_flag_t flag, int n_max_rows, int n_max_cols)
 Allocate and initialize a cs_sdm_t structure Most generic function to create a cs_sdm_t structure. More...
 
cs_sdm_t * cs_sdm_square_create (int n_max_rows)
 Allocate and initialize a cs_sdm_t structure Case of a square matrix. More...
 
cs_sdm_t * cs_sdm_create_copy (const cs_sdm_t *m)
 Allocate a cs_sdm_t structure and initialized it with the copy of the matrix m in input. More...
 
cs_sdm_t * cs_sdm_create_transpose (cs_sdm_t *mat)
 Define a new matrix which is its transpose. More...
 
cs_sdm_t * cs_sdm_block_create (int n_max_blocks_by_row, int n_max_blocks_by_col, const int max_row_block_sizes[], const int max_col_block_sizes[])
 Allocate and initialize a cs_sdm_t structure. More...
 
cs_sdm_t * cs_sdm_block33_create (int n_max_blocks_by_row, int n_max_blocks_by_col)
 Allocate and initialize a cs_sdm_t structure by block when the block size is constant and equal to 3. More...
 
static void cs_sdm_map_array (int n_max_rows, int n_max_cols, cs_sdm_t *m, cs_real_t *array)
 Map an array into a predefined cs_sdm_t structure. This array is shared and the lifecycle of this array is not managed by the cs_sdm_t structure. More...
 
cs_sdm_t * cs_sdm_free (cs_sdm_t *mat)
 Free a cs_sdm_t structure. More...
 
static void cs_sdm_init (int n_rows, int n_cols, cs_sdm_t *mat)
 Initialize a cs_sdm_t structure Case of a square matrix. More...
 
static void cs_sdm_square_init (int n_rows, cs_sdm_t *mat)
 Initialize a cs_sdm_t structure Case of a square matrix. More...
 
void cs_sdm_block_init (cs_sdm_t *m, int n_row_blocks, int n_col_blocks, const int row_block_sizes[], const int col_block_sizes[])
 Initialize the pattern of cs_sdm_t structure defined by block The matrix should have been allocated before calling this function. More...
 
void cs_sdm_block33_init (cs_sdm_t *m, int n_row_blocks, int n_col_blocks)
 Initialize the pattern of cs_sdm_t structure defined by 3x3 block The matrix should have been allocated before calling this function. More...
 
static void cs_sdm_copy (cs_sdm_t *recv, const cs_sdm_t *send)
 Copy a cs_sdm_t structure into another cs_sdm_t structure which has been already allocated. More...
 
cs_sdm_t * cs_sdm_block_create_copy (const cs_sdm_t *mref)
 Allocate and initialize a cs_sdm_t structure w.r.t. to a given matrix. More...
 
static cs_sdm_t * cs_sdm_get_block (const cs_sdm_t *const m, int row_block_id, int col_block_id)
 Get a specific block in a cs_sdm_t structure defined by block. More...
 
static void cs_sdm_get_col (const cs_sdm_t *m, int col_id, cs_real_t *col_vals)
 Get a copy of a column in a preallocated vector. More...
 
static void cs_sdm_copy_block (const cs_sdm_t *m, const short int r_id, const short int c_id, const short int nr, const short int nc, cs_sdm_t *b)
 copy a block of a matrix into a sub-matrix starting from (r_id, c_id) with a size of nr rows and nc cols More...
 
static void cs_sdm_transpose_and_update (const cs_sdm_t *m, cs_sdm_t *mt)
 transpose and copy a matrix into another one already shaped sub-matrix starting from (r_id, c_id) More...
 
void cs_sdm_multiply (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a local dense matrix-product c = a*b c has been previously allocated. More...
 
static void cs_sdm_multiply_r1c3_rowrow (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Specific version: a (1x3) and b (3x1) More...
 
void cs_sdm_multiply_rowrow (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Generic version: all compatible sizes. More...
 
void cs_sdm_multiply_rowrow_sym (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Generic version: all compatible sizes Result is known to be symmetric. More...
 
void cs_sdm_block_multiply_rowrow (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Case of a matrix defined by block. More...
 
void cs_sdm_block_multiply_rowrow_sym (const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)
 Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Case of matrices defined by block. Results is known to be symmetric. More...
 
void cs_sdm_square_matvec (const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Compute a dense matrix-vector product for a small square matrix mv has been previously allocated. More...
 
void cs_sdm_matvec (const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Compute a dense matrix-vector product for a rectangular matrix mv has been previously allocated. More...
 
void cs_sdm_update_matvec (const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Compute a dense matrix-vector product for a rectangular matrix mv has been previously allocated and initialized Thus mv is updated inside this function. More...
 
void cs_sdm_matvec_transposed (const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Compute a dense matrix-vector product for a rectangular matrix which is transposed. mv has been previously allocated. mv is updated inside this function. Don't forget to initialize mv if needed. More...
 
void cs_sdm_block_add (cs_sdm_t *mat, const cs_sdm_t *add)
 Add two matrices defined by block: loc += add. More...
 
void cs_sdm_block_add_mult (cs_sdm_t *mat, cs_real_t mult_coef, const cs_sdm_t *add)
 Add two matrices defined by block: loc += mult_coef * add. More...
 
void cs_sdm_block_matvec (const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)
 Compute a dense matrix-vector product for a rectangular matrix defined by block mv has been previously allocated. More...
 
void cs_sdm_add (cs_sdm_t *mat, const cs_sdm_t *add)
 Add two small dense matrices: loc += add. More...
 
void cs_sdm_add_mult (cs_sdm_t *mat, cs_real_t alpha, const cs_sdm_t *add)
 Add two small dense matrices: loc += alpha*add. More...
 
void cs_sdm_square_add_transpose (cs_sdm_t *mat, cs_sdm_t *tr)
 Define a new matrix by adding the given matrix with its transpose. Keep the transposed matrix for a future use. More...
 
void cs_sdm_square_2symm (cs_sdm_t *mat)
 Set the given matrix to two times its symmetric part mat --> mat + mat_tr = 2*symm(mat) More...
 
void cs_sdm_square_asymm (cs_sdm_t *mat)
 Set the given matrix into its anti-symmetric part. More...
 
void cs_sdm_block_square_asymm (cs_sdm_t *mat)
 Set the given block matrix into its anti-symmetric part. More...
 
void cs_sdm_33_sym_qr_compute (const cs_real_t m[9], cs_real_t Qt[9], cs_real_t R[6])
 Decompose a matrix into the matrix product Q.R Case of a 3x3 symmetric matrix. More...
 
void cs_sdm_33_lu_compute (const cs_sdm_t *m, cs_real_t facto[9])
 LU factorization of a small dense 3x3 matrix. More...
 
void cs_sdm_lu_compute (const cs_sdm_t *m, cs_real_t facto[])
 LU factorization of a small dense matrix. Small means that the number m->n_rows is less than 100 for instance. More...
 
void cs_sdm_33_lu_solve (const cs_real_t facto[9], const cs_real_t rhs[3], cs_real_t sol[3])
 Solve a system A.sol = rhs using a LU factorization of A (a small 3x3 dense matrix). More...
 
void cs_sdm_lu_solve (cs_lnum_t n_rows, const cs_real_t facto[], const cs_real_t *rhs, cs_real_t *sol)
 Solve a system A.sol = rhs using a LU factorization of A (a small dense matrix). More...
 
void cs_sdm_33_ldlt_compute (const cs_sdm_t *m, cs_real_t facto[6])
 LDL^T: Modified Cholesky decomposition of a 3x3 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_44_ldlt_compute (const cs_sdm_t *m, cs_real_t facto[10])
 LDL^T: Modified Cholesky decomposition of a 4x4 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_66_ldlt_compute (const cs_sdm_t *m, cs_real_t facto[21])
 LDL^T: Modified Cholesky decomposition of a 6x6 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_ldlt_compute (const cs_sdm_t *m, cs_real_t *facto, cs_real_t *dkk)
 LDL^T: Modified Cholesky decomposition of a SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_33_ldlt_solve (const cs_real_t facto[6], const cs_real_t rhs[3], cs_real_t sol[3])
 Solve a 3x3 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_44_ldlt_solve (const cs_real_t facto[10], const cs_real_t rhs[4], cs_real_t x[4])
 Solve a 4x4 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_66_ldlt_solve (const cs_real_t f[21], const cs_real_t b[6], cs_real_t x[6])
 Solve a 6x6 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
void cs_sdm_ldlt_solve (int n_rows, const cs_real_t *facto, const cs_real_t *rhs, cs_real_t *sol)
 Solve a SPD matrix with a L.D.L^T (Modified Cholesky decomposition) The solution should be already allocated Reference: http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf. More...
 
double cs_sdm_test_symmetry (const cs_sdm_t *mat)
 Test if a matrix is symmetric. Return 0. if the extradiagonal differences are lower thann the machine precision. More...
 
void cs_sdm_simple_dump (const cs_sdm_t *mat)
 Dump a small dense matrix. More...
 
void cs_sdm_dump (cs_lnum_t parent_id, const cs_lnum_t *row_ids, const cs_lnum_t *col_ids, const cs_sdm_t *mat)
 Dump a small dense matrix. More...
 
void cs_sdm_fprintf (FILE *fp, const char *fname, cs_real_t thd, const cs_sdm_t *m)
 Print a cs_sdm_t structure not defined by block Print into the file f if given otherwise open a new file named fname if given otherwise print into the standard output The usage of threshold allows one to compare more easier matrices without taking into account numerical roundoff. More...
 
void cs_sdm_block_dump (cs_lnum_t parent_id, const cs_sdm_t *mat)
 Dump a small dense matrix defined by blocks. More...
 
void cs_sdm_block_fprintf (FILE *fp, const char *fname, cs_real_t thd, const cs_sdm_t *m)
 Print a cs_sdm_t structure which is defined by block Print into the file f if given otherwise open a new file named fname if given otherwise print into the standard output The usage of threshold allows one to compare more easier matrices without taking into account numerical roundoff. More...
 

Macro Definition Documentation

◆ CS_SDM_BY_BLOCK

#define CS_SDM_BY_BLOCK   (1 << 0) /* Matrix is defined by block */

◆ CS_SDM_SHARED_VAL

#define CS_SDM_SHARED_VAL   (1 << 2) /* Matrix is not owner of its values */

◆ CS_SDM_SYMMETRIC

#define CS_SDM_SYMMETRIC   (1 << 1) /* Matrix is symmetric by construction */

Typedef Documentation

◆ cs_sdm_matvec_t

typedef void() cs_sdm_matvec_t(const cs_sdm_t *mat, const cs_real_t *vec, cs_real_t *mv)

Generic prototype for computing a dense matrix-vector product mv has been previously allocated.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use
[in,out]mvresult of the local matrix-vector product

◆ cs_sdm_product_t

typedef void() cs_sdm_product_t(const cs_sdm_t *a, const cs_sdm_t *b, cs_sdm_t *c)

Generic prototype for computing a local dense matrix-product c = a*b where c has been previously allocated.

Parameters
[in]alocal dense matrix to use
[in]blocal dense matrix to use
[in,out]cresult of the local matrix-product

Function Documentation

◆ cs_sdm_33_ldlt_compute()

void cs_sdm_33_ldlt_compute ( const cs_sdm_t *  m,
cs_real_t  facto[6] 
)

LDL^T: Modified Cholesky decomposition of a 3x3 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factovector of the coefficient of the decomposition
Note
: facto is a lower triangular matrix. The first value of the j-th (zero-based) row is easily accessed: its index is j*(j+1)/2 (cf sum of the first j natural numbers). Instead of 1 on the diagonal we store the inverse of D mat in the L.D.L^T decomposition

◆ cs_sdm_33_ldlt_solve()

void cs_sdm_33_ldlt_solve ( const cs_real_t  facto[6],
const cs_real_t  rhs[3],
cs_real_t  sol[3] 
)

Solve a 3x3 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]factovector of the coefficients of the decomposition
[in]rhsright-hand side
[in,out]solsolution

◆ cs_sdm_33_lu_compute()

void cs_sdm_33_lu_compute ( const cs_sdm_t *  m,
cs_real_t  facto[9] 
)

LU factorization of a small dense 3x3 matrix.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factocompact storage of coefficients for the LU factorization
Note
: facto stores L the lower triangular matrix (without its diagonal entries assumed to be equal to 1) and U the upper triangular matrix.

◆ cs_sdm_33_lu_solve()

void cs_sdm_33_lu_solve ( const cs_real_t  facto[9],
const cs_real_t  rhs[3],
cs_real_t  sol[3] 
)

Solve a system A.sol = rhs using a LU factorization of A (a small 3x3 dense matrix).

Parameters
[in]factocompact storage of coefficients for the LU factorization (should be allocated to the right size, i.e. n_rows*n_rows)
[in]rhsright-hand side
[in,out]solsolution
Note
: facto stores L the lower triangular matrix (without its diagonal entries assumed to be equal to 1) and U the upper triangular matrix.

◆ cs_sdm_33_sym_qr_compute()

void cs_sdm_33_sym_qr_compute ( const cs_real_t  m[9],
cs_real_t  Qt[9],
cs_real_t  R[6] 
)

Decompose a matrix into the matrix product Q.R Case of a 3x3 symmetric matrix.

Parameters
[in]mmatrix values
[in,out]Qttransposed of matrix Q
[in,out]Rvector of the coefficient of the decomposition
Note
: R is an upper triangular matrix. Stored in a compact way.

j= 0, 1, 2 i=0| 0| 1| 2| i=1 | 4| 5| i=2 6|

◆ cs_sdm_44_ldlt_compute()

void cs_sdm_44_ldlt_compute ( const cs_sdm_t *  m,
cs_real_t  facto[10] 
)

LDL^T: Modified Cholesky decomposition of a 4x4 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factovector of the coefficient of the decomposition
Note
: facto is a lower triangular matrix. The first value of the j-th (zero-based) row is easily accessed: its index is j*(j+1)/2 (cf sum of the first j natural numbers). Instead of 1 on the diagonal we store the inverse of D mat in the L.D.L^T decomposition

◆ cs_sdm_44_ldlt_solve()

void cs_sdm_44_ldlt_solve ( const cs_real_t  facto[10],
const cs_real_t  rhs[4],
cs_real_t  x[4] 
)

Solve a 4x4 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]factovector of the coefficients of the decomposition
[in]rhsright-hand side
[in,out]xsolution

◆ cs_sdm_66_ldlt_compute()

void cs_sdm_66_ldlt_compute ( const cs_sdm_t *  m,
cs_real_t  facto[21] 
)

LDL^T: Modified Cholesky decomposition of a 6x6 SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factovector of the coefficient of the decomposition
Note
: facto is a lower triangular matrix. The first value of the j-th (zero-based) row is easily accessed: its index is j*(j+1)/2 (cf sum of the first j natural numbers). Instead of 1 on the diagonal we store the inverse of D mat in the L.D.L^T decomposition

◆ cs_sdm_66_ldlt_solve()

void cs_sdm_66_ldlt_solve ( const cs_real_t  f[21],
const cs_real_t  b[6],
cs_real_t  x[6] 
)

Solve a 6x6 matrix with a modified Cholesky decomposition (L.D.L^T) The solution should be already allocated Ref. http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]fvector of the coefficients of the decomposition
[in]bright-hand side
[in,out]xsolution

◆ cs_sdm_add()

void cs_sdm_add ( cs_sdm_t *  mat,
const cs_sdm_t *  add 
)

Add two small dense matrices: loc += add.

Parameters
[in,out]matlocal matrix storing the result
[in]addvalues to add to mat

◆ cs_sdm_add_mult()

void cs_sdm_add_mult ( cs_sdm_t *  mat,
cs_real_t  alpha,
const cs_sdm_t *  add 
)

Add two small dense matrices: loc += alpha*add.

Parameters
[in,out]matlocal matrix storing the result
[in]alphamultiplicative coefficient
[in]addvalues to add to mat

◆ cs_sdm_add_scalvect()

static void cs_sdm_add_scalvect ( int  n,
const cs_real_t  s,
const cs_real_t  x[],
cs_real_t  y[] 
)
inlinestatic

Multiply a small vector by a scalar coefficient: y += s x For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas.

Parameters
[in]nsize of arrays x and y (small)
[in]sscalar coefficient
[in]xin array
[in,out]yout array

◆ cs_sdm_block33_create()

cs_sdm_t* cs_sdm_block33_create ( int  n_max_blocks_by_row,
int  n_max_blocks_by_col 
)

Allocate and initialize a cs_sdm_t structure by block when the block size is constant and equal to 3.

Parameters
[in]n_max_blocks_by_rowmax number of blocks in a row
[in]n_max_blocks_by_colmax number of blocks in a column
Returns
a new allocated cs_sdm_t structure

◆ cs_sdm_block33_init()

void cs_sdm_block33_init ( cs_sdm_t *  m,
int  n_row_blocks,
int  n_col_blocks 
)

Initialize the pattern of cs_sdm_t structure defined by 3x3 block The matrix should have been allocated before calling this function.

Parameters
[in,out]m
[in]n_row_blocksnumber of blocks in a row
[in]n_col_blocksnumber of blocks in a column

◆ cs_sdm_block_add()

void cs_sdm_block_add ( cs_sdm_t *  mat,
const cs_sdm_t *  add 
)

Add two matrices defined by block: loc += add.

Parameters
[in,out]matlocal matrix storing the result
[in]addvalues to add to mat

◆ cs_sdm_block_add_mult()

void cs_sdm_block_add_mult ( cs_sdm_t *  mat,
cs_real_t  mult_coef,
const cs_sdm_t *  add 
)

Add two matrices defined by block: loc += mult_coef * add.

Parameters
[in,out]matlocal matrix storing the result
[in]mult_coefmultiplicative coefficient
[in]addvalues to add to mat

◆ cs_sdm_block_create()

cs_sdm_t* cs_sdm_block_create ( int  n_max_blocks_by_row,
int  n_max_blocks_by_col,
const int  max_row_block_sizes[],
const int  max_col_block_sizes[] 
)

Allocate and initialize a cs_sdm_t structure.

Parameters
[in]n_max_blocks_by_rowmax number of blocks in a row
[in]n_max_blocks_by_colmax number of blocks in a column
[in]max_row_block_sizesmax number of rows by block in a column
[in]max_col_block_sizesmax number of columns by block in a row
Returns
a new allocated cs_sdm_t structure

◆ cs_sdm_block_create_copy()

cs_sdm_t* cs_sdm_block_create_copy ( const cs_sdm_t *  mref)

Allocate and initialize a cs_sdm_t structure w.r.t. to a given matrix.

Parameters
[in]mrefpointer to a matrix to copy
Returns
a new allocated cs_sdm_t structure which is a copy of mref

◆ cs_sdm_block_dump()

void cs_sdm_block_dump ( cs_lnum_t  parent_id,
const cs_sdm_t *  mat 
)

Dump a small dense matrix defined by blocks.

Parameters
[in]parent_idid of the related parent entity
[in]matpointer to the cs_sdm_t structure

◆ cs_sdm_block_fprintf()

void cs_sdm_block_fprintf ( FILE *  fp,
const char *  fname,
cs_real_t  thd,
const cs_sdm_t *  m 
)

Print a cs_sdm_t structure which is defined by block Print into the file f if given otherwise open a new file named fname if given otherwise print into the standard output The usage of threshold allows one to compare more easier matrices without taking into account numerical roundoff.

Parameters
[in]fppointer to a file structure or NULL
[in]fnamefilename or NULL
[in]thdthreshold (below this value --> set 0)
[in]mpointer to the cs_sdm_t structure

◆ cs_sdm_block_init()

void cs_sdm_block_init ( cs_sdm_t *  m,
int  n_row_blocks,
int  n_col_blocks,
const int  row_block_sizes[],
const int  col_block_sizes[] 
)

Initialize the pattern of cs_sdm_t structure defined by block The matrix should have been allocated before calling this function.

Parameters
[in,out]m
[in]n_row_blocksnumber of blocks in a row
[in]n_col_blocksnumber of blocks in a column
[in]row_block_sizesnumber of rows by block in a column
[in]col_block_sizesnumber of columns by block in a row

◆ cs_sdm_block_matvec()

void cs_sdm_block_matvec ( const cs_sdm_t *  mat,
const cs_real_t vec,
cs_real_t mv 
)

Compute a dense matrix-vector product for a rectangular matrix defined by block mv has been previously allocated.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use (size = mat->n_cols)
[in,out]mvresult of the operation (size = mat->n_rows)

◆ cs_sdm_block_multiply_rowrow()

void cs_sdm_block_multiply_rowrow ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Case of a matrix defined by block.

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Case of a matrix defined by block.

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated

◆ cs_sdm_block_multiply_rowrow_sym()

void cs_sdm_block_multiply_rowrow_sym ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Case of matrices defined by block. Results is known to be symmetric.

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated

◆ cs_sdm_block_square_asymm()

void cs_sdm_block_square_asymm ( cs_sdm_t *  mat)

Set the given block matrix into its anti-symmetric part.

Parameters
[in,out]matsmall dense matrix defined by block to transform

◆ cs_sdm_copy()

static void cs_sdm_copy ( cs_sdm_t *  recv,
const cs_sdm_t *  send 
)
inlinestatic

Copy a cs_sdm_t structure into another cs_sdm_t structure which has been already allocated.

Parameters
[in,out]recvpointer to a cs_sdm_t struct.
[in]sendpointer to a cs_sdm_t struct.

◆ cs_sdm_copy_block()

static void cs_sdm_copy_block ( const cs_sdm_t *  m,
const short int  r_id,
const short int  c_id,
const short int  nr,
const short int  nc,
cs_sdm_t *  b 
)
inlinestatic

copy a block of a matrix into a sub-matrix starting from (r_id, c_id) with a size of nr rows and nc cols

Parameters
[in]mpointer to cs_sdm_t structure
[in]r_idrow index
[in]c_idcolumn index
[in]nrnumber of rows to extract
[in]ncnumber of column to extract
[in,out]bsubmatrix

◆ cs_sdm_create()

cs_sdm_t* cs_sdm_create ( cs_flag_t  flag,
int  n_max_rows,
int  n_max_cols 
)

Allocate and initialize a cs_sdm_t structure Most generic function to create a cs_sdm_t structure.

Parameters
[in]flagmetadata related to a cs_sdm_t structure
[in]n_max_rowsmax number of rows
[in]n_max_colsmax number of columns
Returns
a new allocated cs_sdm_t structure

Allocate and initialize a cs_sdm_t structure Most generic function to create a cs_sdm_t structure.

Parameters
[in]flagmetadata related to a cs_sdm_t structure
[in]n_max_rowsmax number of rows
[in]n_max_colsmax number of columns
Returns
a new allocated cs_sdm_t structure

◆ cs_sdm_create_copy()

cs_sdm_t* cs_sdm_create_copy ( const cs_sdm_t *  m)

Allocate a cs_sdm_t structure and initialized it with the copy of the matrix m in input.

Parameters
[in]mpointer to a cs_sdm_t structure to copy
Returns
a pointer to a new allocated cs_sdm_t structure

◆ cs_sdm_create_transpose()

cs_sdm_t* cs_sdm_create_transpose ( cs_sdm_t *  mat)

Define a new matrix which is its transpose.

Parameters
[in]matlocal matrix to transpose
Returns
a pointer to the new allocated transposed matrix

◆ cs_sdm_dot()

static cs_real_t cs_sdm_dot ( int  n,
const cs_real_t  x[],
const cs_real_t  y[] 
)
inlinestatic

Basic dot product for a small vector For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas.

Parameters
[in]nsize of arrays x and y (small)
[in]xfirst array
[in]ysecond array
Returns
the dot product

◆ cs_sdm_dump()

void cs_sdm_dump ( cs_lnum_t  parent_id,
const cs_lnum_t row_ids,
const cs_lnum_t col_ids,
const cs_sdm_t *  mat 
)

Dump a small dense matrix.

Parameters
[in]parent_idid of the related parent entity
[in]row_idslist of ids related to associated entities (or NULL)
[in]col_idslist of ids related to associated entities (or NULL)
[in]matpointer to the cs_sdm_t structure

◆ cs_sdm_fprintf()

void cs_sdm_fprintf ( FILE *  fp,
const char *  fname,
cs_real_t  thd,
const cs_sdm_t *  m 
)

Print a cs_sdm_t structure not defined by block Print into the file f if given otherwise open a new file named fname if given otherwise print into the standard output The usage of threshold allows one to compare more easier matrices without taking into account numerical roundoff.

Parameters
[in]fppointer to a file structure or NULL
[in]fnamefilename or NULL
[in]thdthreshold (below this value --> set 0)
[in]mpointer to the cs_sdm_t structure

◆ cs_sdm_free()

cs_sdm_t* cs_sdm_free ( cs_sdm_t *  mat)

Free a cs_sdm_t structure.

Parameters
[in]matpointer to a cs_sdm_t struct. to free
Returns
a NULL pointer

◆ cs_sdm_get_block()

static cs_sdm_t* cs_sdm_get_block ( const cs_sdm_t *const  m,
int  row_block_id,
int  col_block_id 
)
inlinestatic

Get a specific block in a cs_sdm_t structure defined by block.

Parameters
[in]mpointer to a cs_sdm_t struct.
[in]row_block_idid of the block row, zero-based.
[in]col_block_idid of the block column, zero-based.
Returns
a pointer to a cs_sdm_t structure corresponfing to a block

◆ cs_sdm_get_col()

static void cs_sdm_get_col ( const cs_sdm_t *  m,
int  col_id,
cs_real_t col_vals 
)
inlinestatic

Get a copy of a column in a preallocated vector.

Parameters
[in]mpointer to a cs_sdm_t struct.
[in]col_idid of the column, zero-based.
[in,out]col_valsextracted values

◆ cs_sdm_init()

static void cs_sdm_init ( int  n_rows,
int  n_cols,
cs_sdm_t *  mat 
)
inlinestatic

Initialize a cs_sdm_t structure Case of a square matrix.

Parameters
[in]n_rowscurrent number of rows
[in]n_colscurrent number of columns
[in,out]matpointer to the cs_sdm_t structure to initialize

◆ cs_sdm_ldlt_compute()

void cs_sdm_ldlt_compute ( const cs_sdm_t *  m,
cs_real_t facto,
cs_real_t dkk 
)

LDL^T: Modified Cholesky decomposition of a SPD matrix. For more reference, see for instance http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factovector of the coefficient of the decomposition
[in,out]dkkstore temporary the diagonal (size = n_rows)
Note
: facto is a lower triangular matrix. The first value of the j-th (zero-based) row is easily accessed: its index is j*(j+1)/2 (cf sum of the first j natural numbers). Instead of 1 on the diagonal we store the inverse of D mat in the L.D.L^T decomposition

◆ cs_sdm_ldlt_solve()

void cs_sdm_ldlt_solve ( int  n_rows,
const cs_real_t facto,
const cs_real_t rhs,
cs_real_t sol 
)

Solve a SPD matrix with a L.D.L^T (Modified Cholesky decomposition) The solution should be already allocated Reference: http://mathforcollege.com/nm/mws/gen/04sle/mws_gen_sle_txt_cholesky.pdf.

Parameters
[in]n_rowsdimension of the system to solve
[in]factovector of the coefficients of the decomposition
[in]rhsright-hand side
[in,out]solsolution

◆ cs_sdm_lu_compute()

void cs_sdm_lu_compute ( const cs_sdm_t *  m,
cs_real_t  facto[] 
)

LU factorization of a small dense matrix. Small means that the number m->n_rows is less than 100 for instance.

Parameters
[in]mpointer to a cs_sdm_t structure
[in,out]factocompact storage of coefficients for the LU factorization (should be allocated to the right size, i.e. m->n_rows*m->n_rows)
Note
: facto stores L the lower triangular matrix (without its diagonal entries assumed to be equal to 1) and U the upper triangular matrix.

◆ cs_sdm_lu_solve()

void cs_sdm_lu_solve ( cs_lnum_t  n_rows,
const cs_real_t  facto[],
const cs_real_t rhs,
cs_real_t sol 
)

Solve a system A.sol = rhs using a LU factorization of A (a small dense matrix).

Parameters
[in]n_rowsdimension of the system to solve
[in]factocompact storage of coefficients for the LU factorization (should be allocated to the right size, i.e. n_rows*n_rows)
[in]rhsright-hand side
[in,out]solsolution
Note
: facto stores L the lower triangular matrix (without its diagonal entries assumed to be equal to 1) and U the upper triangular matrix.

◆ cs_sdm_map_array()

static void cs_sdm_map_array ( int  n_max_rows,
int  n_max_cols,
cs_sdm_t *  m,
cs_real_t array 
)
inlinestatic

Map an array into a predefined cs_sdm_t structure. This array is shared and the lifecycle of this array is not managed by the cs_sdm_t structure.

Parameters
[in]n_max_rowsmax. number of rows
[in]n_max_colsmax. number of columns
[in,out]mpointer to a cs_sdm_t structure to set
[in,out]arraypointer to an array of values of size equal to n_max_rows x n_max_cols

◆ cs_sdm_matvec()

void cs_sdm_matvec ( const cs_sdm_t *  mat,
const cs_real_t vec,
cs_real_t mv 
)

Compute a dense matrix-vector product for a rectangular matrix mv has been previously allocated.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use (size = mat->n_cols)
[in,out]mvresult of the operation (size = mat->n_rows)

◆ cs_sdm_matvec_transposed()

void cs_sdm_matvec_transposed ( const cs_sdm_t *  mat,
const cs_real_t vec,
cs_real_t mv 
)

Compute a dense matrix-vector product for a rectangular matrix which is transposed. mv has been previously allocated. mv is updated inside this function. Don't forget to initialize mv if needed.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use (size = mat->n_cols)
[in,out]mvresult of the operation (size = mat->n_rows)

◆ cs_sdm_multiply()

void cs_sdm_multiply ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)

Compute a local dense matrix-product c = a*b c has been previously allocated.

Parameters
[in]alocal dense matrix to use
[in]blocal dense matrix to use
[in,out]cresult of the local matrix-product is updated

◆ cs_sdm_multiply_r1c3_rowrow()

static void cs_sdm_multiply_r1c3_rowrow ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)
inlinestatic

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Specific version: a (1x3) and b (3x1)

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated

◆ cs_sdm_multiply_rowrow()

void cs_sdm_multiply_rowrow ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Generic version: all compatible sizes.

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated.

◆ cs_sdm_multiply_rowrow_sym()

void cs_sdm_multiply_rowrow_sym ( const cs_sdm_t *  a,
const cs_sdm_t *  b,
cs_sdm_t *  c 
)

Compute a row-row matrix product of a and b. It is basically equal to the classical a*b^T. It is a fast (matrices are row-major) way of computing a*b if b is symmetric or if b^T is given. Generic version: all compatible sizes Result is known to be symmetric.

Parameters
[in]alocal matrix to use
[in]blocal matrix to use
[in,out]cresult of the local matrix-product (already allocated) is updated.

◆ cs_sdm_scalvect()

static void cs_sdm_scalvect ( int  n,
const cs_real_t  s,
const cs_real_t  x[],
cs_real_t  y[] 
)
inlinestatic

Multiply a small vector by a scalar coefficient: y = s x For very small array sizes (3, 4, 6) prefer functions in cs_math For large array sizes ( from 10^3 to ..) prefer functions in cs_blas.

Parameters
[in]nsize of arrays x and y (small)
[in]sscalar coefficient
[in]xin array
[in,out]yout array

◆ cs_sdm_simple_dump()

void cs_sdm_simple_dump ( const cs_sdm_t *  mat)

Dump a small dense matrix.

Parameters
[in]matpointer to the cs_sdm_t structure

◆ cs_sdm_square_2symm()

void cs_sdm_square_2symm ( cs_sdm_t *  mat)

Set the given matrix to two times its symmetric part mat --> mat + mat_tr = 2*symm(mat)

Parameters
[in,out]matsmall dense matrix to transform

◆ cs_sdm_square_add_transpose()

void cs_sdm_square_add_transpose ( cs_sdm_t *  mat,
cs_sdm_t *  tr 
)

Define a new matrix by adding the given matrix with its transpose. Keep the transposed matrix for a future use.

Parameters
[in,out]matlocal matrix to transpose and add
[in,out]trtransposed of the local matrix

◆ cs_sdm_square_asymm()

void cs_sdm_square_asymm ( cs_sdm_t *  mat)

Set the given matrix into its anti-symmetric part.

Parameters
[in,out]matsmall dense matrix to transform

◆ cs_sdm_square_create()

cs_sdm_t* cs_sdm_square_create ( int  n_max_rows)

Allocate and initialize a cs_sdm_t structure Case of a square matrix.

Parameters
[in]n_max_rowsmax number of rows
Returns
a new allocated cs_sdm_t structure

◆ cs_sdm_square_init()

static void cs_sdm_square_init ( int  n_rows,
cs_sdm_t *  mat 
)
inlinestatic

Initialize a cs_sdm_t structure Case of a square matrix.

Parameters
[in]n_rowscurrent number of rows
[in,out]matpointer to the cs_sdm_t structure to initialize

◆ cs_sdm_square_matvec()

void cs_sdm_square_matvec ( const cs_sdm_t *  mat,
const cs_real_t vec,
cs_real_t mv 
)

Compute a dense matrix-vector product for a small square matrix mv has been previously allocated.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use
[in,out]mvresult of the local matrix-vector product

◆ cs_sdm_test_symmetry()

double cs_sdm_test_symmetry ( const cs_sdm_t *  mat)

Test if a matrix is symmetric. Return 0. if the extradiagonal differences are lower thann the machine precision.

Parameters
[in]matpointer to the cs_sdm_t structure to test
Returns
0 if the matrix is symmetric at the machine tolerance otherwise the absolute max. value between two transposed terms

◆ cs_sdm_transpose_and_update()

static void cs_sdm_transpose_and_update ( const cs_sdm_t *  m,
cs_sdm_t *  mt 
)
inlinestatic

transpose and copy a matrix into another one already shaped sub-matrix starting from (r_id, c_id)

Parameters
[in]mpointer to cs_sdm_t structure
[in,out]mtmatrix to update with the transposed of m

◆ cs_sdm_update_matvec()

void cs_sdm_update_matvec ( const cs_sdm_t *  mat,
const cs_real_t vec,
cs_real_t mv 
)

Compute a dense matrix-vector product for a rectangular matrix mv has been previously allocated and initialized Thus mv is updated inside this function.

Parameters
[in]matlocal matrix to use
[in]veclocal vector to use (size = mat->n_cols)
[in,out]mvresult of the operation (size = mat->n_rows)