SpM Handbook 1.2.4
Loading...
Searching...
No Matches
spm.h File Reference
#include <stdlib.h>
#include <assert.h>
#include <stdio.h>
#include "spm/config.h"
#include "spm/const.h"
#include "spm/datatypes.h"
#include "spm/mpi.h"
#include "spm/z_spm.h"
#include "spm/c_spm.h"
#include "spm/d_spm.h"
#include "spm/s_spm.h"
#include "spm/p_spm.h"

Go to the source code of this file.

Data Structures

struct  spmatrix_s
 The sparse matrix data structure. More...
 

Typedefs

typedef struct spmatrix_s spmatrix_t
 The sparse matrix data structure.
 

Functions

int spmBlasGetNumThreads (void)
 Return the current number of threads used for blas calls.
 
int spmBlasSetNumThreads (int nt)
 Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) and return the previous number of blas threads.
 
int spmBlasSetNumThreadsOne (void)
 Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) to 1 and return the previous number of blas threads.
 
SPM basic subroutines
void spmInit (spmatrix_t *spm)
 Init the spm structure.
 
void spmInitDist (spmatrix_t *spm, SPM_Comm comm)
 Init the spm structure with a specific communicator.
 
void spmAlloc (spmatrix_t *spm)
 Allocate the arrays of an spm structure.
 
void spmExit (spmatrix_t *spm)
 Cleanup the spm structure but do not free the spm pointer.
 
void spmCopy (const spmatrix_t *spm, spmatrix_t *newspm)
 Create a copy of the spm.
 
void spmBase (spmatrix_t *spm, int baseval)
 Rebase the arrays of the spm to the given value.
 
spm_int_t spmFindBase (const spmatrix_t *spm)
 Search the base used in the spm structure.
 
spm_int_t spmGetDegree (const spmatrix_t *spm)
 
int spmConvert (int ofmttype, spmatrix_t *spm)
 Convert the storage format of the spm.
 
void spmUpdateComputedFields (spmatrix_t *spm)
 Update all the computed fields based on the static values stored.
 
void spmGenFakeValues (spmatrix_t *spm)
 Generate the fake values array such that $ M =  \alpha * D - \beta * A $.
 
SPM distribution subroutines
int spmScatter (spmatrix_t *newspm, int root, const spmatrix_t *oldspm, spm_int_t n, const spm_int_t *loc2glob, int distByColumn, SPM_Comm comm)
 Scatter the SPM thanks to loc2glob.
 
int spmGather (const spmatrix_t *spm_scattered, int root, spmatrix_t *opt_spm_gathered)
 Gather a distributed Sparse Matrix on a single node.
 
int spmGatherInPlace (spmatrix_t *spm)
 This routine performs a allgather of a distributed spm in place.
 
int spmRedistribute (const spmatrix_t *spm, spm_int_t new_n, const spm_int_t *newl2g, spmatrix_t *newspm)
 Create a copy of a given SPM with a new distribution corresponding to the given loc2glob array.
 
SPM BLAS subroutines
double spmNorm (spm_normtype_t ntype, const spmatrix_t *spm)
 Compute the norm of the spm.
 
double spmNormVec (spm_normtype_t ntype, const spmatrix_t *spm, const void *x, spm_int_t inc)
 Compute the norm of the spm.
 
double spmNormMat (spm_normtype_t ntype, const spmatrix_t *spm, spm_int_t n, const void *A, spm_int_t lda)
 Compute the norm of the spm.
 
int spmMatVec (spm_trans_t trans, double alpha, const spmatrix_t *spm, const void *x, double beta, void *y)
 Compute a matrix-vector product.
 
int spmMatMat (spm_trans_t trans, spm_int_t n, double alpha, const spmatrix_t *A, const void *B, spm_int_t ldb, double beta, void *C, spm_int_t ldc)
 Compute a matrix-matrix product.
 
void spmScal (double alpha, spmatrix_t *spm)
 Scale the spm.
 
void spmScalVec (double alpha, const spmatrix_t *spm, void *x, spm_int_t incx)
 Scale a vector associated to a sparse matrix. Deprecated function replaced by spmScalVec() or spmScalMat().
 
void spmScalMat (double alpha, const spmatrix_t *spm, spm_int_t n, void *A, spm_int_t lda)
 Scale a matrix associated to a sparse matrix.
 
void spmScalMatrix (double alpha, spmatrix_t *spm)
 Scale the spm.
 
void spmScalVector (spm_coeftype_t flt, double alpha, spm_int_t n, void *x, spm_int_t incx)
 Scale a vector according to the spm type.
 
SPM subroutines to check format
int spmSort (spmatrix_t *spm)
 Sort the subarray of edges of each vertex.
 
spm_int_t spmMergeDuplicate (spmatrix_t *spm)
 Merge multiple entries in a spm by summing their values together.
 
spm_int_t spmSymmetrize (spmatrix_t *spm)
 Symmetrize the pattern of the spm.
 
int spmCheckAndCorrect (const spmatrix_t *spm_in, spmatrix_t *spm_out)
 Check the correctness of a spm.
 
SPM subroutines to check factorization/solve
int spmGenMat (spm_rhstype_t type, spm_int_t nrhs, const spmatrix_t *spm, void *alpha, unsigned long long int seed, void *A, spm_int_t lda)
 Generate a set of vectors associated to a given matrix.
 
int spmGenVec (spm_rhstype_t type, const spmatrix_t *spm, void *alpha, unsigned long long int seed, void *x, spm_int_t incx)
 Generate a vector associated to a given matrix.
 
int spmGenRHS (spm_rhstype_t type, spm_int_t nrhs, const spmatrix_t *spm, void *x, spm_int_t ldx, void *b, spm_int_t ldb)
 Generate right hand side vectors associated to a given matrix.
 
int spmCheckAxb (double eps, spm_int_t nrhs, const spmatrix_t *spm, void *x0, spm_int_t ldx0, void *b, spm_int_t ldb, const void *x, spm_int_t ldx)
 Check the backward error, and the forward error if x0 is provided.
 
SPM subroutines to manipulate RHS
int spmExtractLocalRHS (spm_int_t nrhs, const spmatrix_t *spm, const void *bglob, spm_int_t ldbg, void *bloc, spm_int_t ldbl)
 Stores the local values of a global RHS in the local one thanks to spm distribution.
 
int spmReduceRHS (spm_int_t nrhs, const spmatrix_t *spm, void *bglob, spm_int_t ldbg, void *bloc, spm_int_t ldbl)
 Reduce an RHS thanks to spm distribution.
 
int spmGatherRHS (spm_int_t nrhs, const spmatrix_t *spm, const void *bloc, spm_int_t ldbl, int root, void *bglob, spm_int_t ldbg)
 Gather an RHS thanks to spm distribution.
 
SPM subroutines to manipulate integers arrays
void spmIntConvert (spm_int_t n, const int *input, spm_int_t *output)
 Convert integer array to spm_int_t format.
 
void spmIntSort1Asc1 (void *const pbase, spm_int_t n)
 
void spmIntSort2Asc1 (void *const pbase, spm_int_t n)
 
void spmIntSort2Asc2 (void *const pbase, spm_int_t n)
 
void spmIntMSortIntAsc (void **const pbase, spm_int_t n)
 
SPM IO subroutines
int spmLoadDist (spmatrix_t *spm, const char *filename, SPM_Comm comm)
 Load the spm structure from a file (internal format).
 
int spmLoad (spmatrix_t *spm, const char *filename)
 Load the spm structure from a file (internal format).
 
int spmSave (const spmatrix_t *spm, const char *filename)
 Save the spm structure into a file (internal format).
 
SPM driver
int spmReadDriver (spm_driver_t driver, const char *filename, spmatrix_t *spm)
 Import a matrix file into a spm structure.
 
int spmReadDriverDist (spm_driver_t driver, const char *filename, spmatrix_t *spm, SPM_Comm comm)
 Import a matrix file into an spm structure for a specific communicator.
 
int spmParseLaplacianInfo (const char *filename, spm_coeftype_t *flttype, spm_int_t *dim1, spm_int_t *dim2, spm_int_t *dim3, double *alpha, double *beta, spm_int_t *dof)
 Parse information given through the filename string to configure the laplacian matrix to generate.
 
SPM debug subroutines
void spm2Dense (const spmatrix_t *spm, void *A)
 Convert the spm matrix into a dense matrix for test purpose.
 
void spmPrint (const spmatrix_t *spm, FILE *stream)
 Print an spm matrix into into a given file.
 
void spmPrintRHS (const spmatrix_t *spm, int nrhs, const void *x, spm_int_t ldx, FILE *stream)
 Print a set of vector associated to an spm matrix.
 
void spmPrintInfo (const spmatrix_t *spm, FILE *stream)
 Print basic informations about the spm matrix into a given stream.
 
void spmExpand (const spmatrix_t *spm_in, spmatrix_t *spm_out)
 Expand a multi-dof spm matrix into an spm with constant dof set to 1.
 
int spmDofExtend (const spmatrix_t *spm, const int type, const int dof, spmatrix_t *newspm)
 Generate a random multidof spm from a given spm (with dof=1).
 

SPM dev printing subroutines

#define p_spmPrintElt(f, i, j, A)
 Subroutines to print one element of an spm structure.
 
static void z_spmPrintElt (FILE *f, spm_int_t i, spm_int_t j, spm_complex64_t A)
 Subroutines to print one element of an spm structure.
 
static void c_spmPrintElt (FILE *f, spm_int_t i, spm_int_t j, spm_complex32_t A)
 Subroutines to print one element of an spm structure.
 
static void d_spmPrintElt (FILE *f, spm_int_t i, spm_int_t j, double A)
 Subroutines to print one element of an spm structure.
 
static void s_spmPrintElt (FILE *f, spm_int_t i, spm_int_t j, float A)
 Subroutines to print one element of an spm structure.
 

Detailed Description

SParse Matrix package header.

Version
1.2.4
Author
Mathieu Faverge
Pierre Ramet
Tony Delarue
Matias Hastaran
Florent Pruvost
Date
2024-06-25

Definition in file spm.h.

Macro Definition Documentation

◆ p_spmPrintElt

#define p_spmPrintElt (   f,
  i,
  j,
 
)
Value:
{ \
fprintf( f, "%ld %ld\n", (long)( i ), (long)( j ) ); \
}

Subroutines to print one element of an spm structure.

Parameters
[in]fPointer to the file
[in]iRow index of the element
[in]jColumn index of the element
[in]AValue of the element A|i,j]

Double complex case

Pattern case

Remarks
: uses a macro to avoid accessing A that would generate segfault.

Definition at line 344 of file spm.h.

Function Documentation

◆ c_spmPrintElt()

static void c_spmPrintElt ( FILE *  f,
spm_int_t  i,
spm_int_t  j,
spm_complex32_t  A 
)
static

Subroutines to print one element of an spm structure.

Parameters
[in]fPointer to the file
[in]iRow index of the element
[in]jColumn index of the element
[in]AValue of the element A|i,j]

Double complex case

Single complex case

Definition at line 316 of file spm.h.

Referenced by c_spm_print_elt_gen_col(), c_spm_print_elt_gen_row(), c_spm_print_elt_sym_diag(), c_spmDensePrint(), and c_spmPrintRHS().

◆ d_spmPrintElt()

static void d_spmPrintElt ( FILE *  f,
spm_int_t  i,
spm_int_t  j,
double  A 
)
static

Subroutines to print one element of an spm structure.

Parameters
[in]fPointer to the file
[in]iRow index of the element
[in]jColumn index of the element
[in]AValue of the element A|i,j]

Double complex case

Double real case

Definition at line 325 of file spm.h.

Referenced by d_spm_print_elt_gen_col(), d_spm_print_elt_gen_row(), d_spm_print_elt_sym_diag(), d_spmDensePrint(), and d_spmPrintRHS().

◆ s_spmPrintElt()

static void s_spmPrintElt ( FILE *  f,
spm_int_t  i,
spm_int_t  j,
float  A 
)
static

Subroutines to print one element of an spm structure.

Parameters
[in]fPointer to the file
[in]iRow index of the element
[in]jColumn index of the element
[in]AValue of the element A|i,j]

Double complex case

Single real case

Definition at line 334 of file spm.h.

Referenced by s_spm_print_elt_gen_col(), s_spm_print_elt_gen_row(), s_spm_print_elt_sym_diag(), s_spmDensePrint(), and s_spmPrintRHS().

◆ spmBlasGetNumThreads()

int spmBlasGetNumThreads ( void  )

Return the current number of threads used for blas calls.

Returns
The number of threads.

Definition at line 2280 of file spm.c.

Referenced by spmBlasSetNumThreads().

◆ spmBlasSetNumThreads()

int spmBlasSetNumThreads ( int  nt)

Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) and return the previous number of blas threads.

Parameters
[in]ntNumber of threads.
Returns
The previous number of blas threads.

Definition at line 2308 of file spm.c.

References spmBlasGetNumThreads().

Referenced by spmBlasSetNumThreadsOne().

◆ spmBlasSetNumThreadsOne()

int spmBlasSetNumThreadsOne ( void  )

Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) to 1 and return the previous number of blas threads.

Returns
The previous number of blas threads.

Definition at line 2334 of file spm.c.

References spmBlasSetNumThreads().