SpM Handbook 1.2.4
All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
SPM: SParse Matrix package

Functions to easily manipulate SPM data structure. More...

Modules

 SpM Constants
 This group describes all the constants of the SpM package.
 

Data Structures

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

Macros

#define SPM_MPI_INT   MPI_INT
 The MPI type associated to spm_int_t.
 
#define SPM_INT_MAX   INT_MAX
 The maximum spm_int_t value.
 

Typedefs

typedef struct spmatrix_s spmatrix_t
 The sparse matrix data structure.
 
typedef int spm_int_t
 The main integer datatype used in spm arrays.
 
typedef unsigned int spm_uint_t
 The main unsigned integer datatype used in spm arrays.
 

Functions

static int spm_load_local (spmatrix_t *spm, const char *filename)
 Load the spm structure from a file (internal format). One node only.
 
static int spm_save_local (const spmatrix_t *spm, const char *filename)
 Save the global spm structure into a file (internal format).
 
static int spm_read_driver (int scatter, spm_driver_t driver, const char *filename, spmatrix_t *spm, SPM_Comm comm)
 Import a matrix file into a spm structure for a specific communicator.
 

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).
 

Detailed Description

Functions to easily manipulate SPM data structure.

Describe all the internals routines of the SParse Matrix package.

This library provides a set of subroutines to manipulate sparse matrices in different format such as compressed sparse column (CSC), compressed sparse row (CSR), or coordinate (IJV) with single or multiple degrees of freedom per unknown. It provides basic BLAS 1 and BLAS 2 functions for those matrices, as well as norms computations and converter tools.


Data Type Documentation

◆ spmatrix_s

struct spmatrix_s

The sparse matrix data structure.

This structure describes matrices with different characteristics that can be useful to any solver:

  • the storage format (SpmCSC, SpmCSR or SpmIJV)
  • the properties (SpmGeneral, SpmHermitian, SpmSymmetric)
  • the base value (0 in C or 1 in Fortran)

All arrays indicated with [+baseval] are using indices in the baseval numbering (C, or Fortran). It is also possible to describe a matrix with constant or variable degrees of freedom.

Definition at line 64 of file spm.h.

Data Fields
spm_mtxtype_t mtxtype

Matrix structure: SpmGeneral, SpmSymmetric or SpmHermitian.

spm_coeftype_t flttype

values datatype: SpmPattern, SpmFloat, SpmDouble, SpmComplex32 or SpmComplex64

spm_fmttype_t fmttype

Matrix storage format: SpmCSC, SpmCSR, SpmIJV

spm_int_t baseval

The base value [C, Fortran] of the indices in the arrays

spm_int_t gN

Global number of vertices in the compressed graph (Computed)

spm_int_t n

Local number of vertices in the compressed graph

spm_int_t gnnz

Global number of non zeroes in the compressed graph (Computed)

spm_int_t nnz

Local number of non zeroes in the compressed graph

spm_int_t gNexp

Global number of vertices in the compressed graph (Computed)

spm_int_t nexp

Local number of vertices in the compressed graph (Computed)

spm_int_t gnnzexp

Global number of non zeroes in the compressed graph (Computed)

spm_int_t nnzexp

Local number of non zeroes in the compressed graph (Computed)

spm_int_t dof

Number of degrees of freedom per unknown, if > 0, constant degree of freedom otherwise, irregular degree of freedom (refer to dofs)

spm_int_t * dofs

Array of the first column of each element in the expanded matrix [+baseval]

spm_layout_t layout

SpmColMajor, or SpmRowMajor

spm_int_t * colptr

List of indirections to rows for each vertex [+baseval]

spm_int_t * rowptr

List of edges for each vertex [+baseval]

spm_int_t * loc2glob

Corresponding numbering from local to global [+baseval]

void * values

Values stored in the matrix

spm_int_t * glob2loc

Corresponding numbering from global to local [0-based], -(owner+1) if remote

int clustnum

Rank of the MPI node

int clustnbr

Number of MPI nodes in the communicator

SPM_Comm comm

Spm communicator to exhange datas

int replicated

0 if the matrix is distributed, 1 if it is replicated on all nodes

Macro Definition Documentation

◆ SPM_MPI_INT

#define SPM_MPI_INT   MPI_INT

The MPI type associated to spm_int_t.

Definition at line 72 of file datatypes.h.

◆ SPM_INT_MAX

#define SPM_INT_MAX   INT_MAX

The maximum spm_int_t value.

Definition at line 73 of file datatypes.h.

Typedef Documentation

◆ spmatrix_t

typedef struct spmatrix_s spmatrix_t

The sparse matrix data structure.

This structure describes matrices with different characteristics that can be useful to any solver:

  • the storage format (SpmCSC, SpmCSR or SpmIJV)
  • the properties (SpmGeneral, SpmHermitian, SpmSymmetric)
  • the base value (0 in C or 1 in Fortran)

All arrays indicated with [+baseval] are using indices in the baseval numbering (C, or Fortran). It is also possible to describe a matrix with constant or variable degrees of freedom.

◆ spm_int_t

The main integer datatype used in spm arrays.

Definition at line 70 of file datatypes.h.

◆ spm_uint_t

The main unsigned integer datatype used in spm arrays.

Definition at line 71 of file datatypes.h.

Function Documentation

◆ spmInit()

void spmInit ( spmatrix_t spm)

Init the spm structure.

Parameters
[in,out]spmThe sparse matrix to init.

Definition at line 152 of file spm.c.

References spmatrix_s::replicated, and spmInitDist().

Referenced by spm_scatter_init(), and spmGather().

◆ spmInitDist()

◆ spmAlloc()

void spmAlloc ( spmatrix_t spm)

Allocate the arrays of an spm structure.

This function must be called after initialization of the non-computed fields, and the call to spmUpdateComputedFields(). It allocates the colptr, rowptr, dof, and values arrays with C malloc function. This is highly recommended to use this function when using PaStiX from Fortran.

Parameters
[in,out]spmThe sparse matrix fr which the internal arrays needs to be allocated.

Definition at line 175 of file spm.c.

References spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, spmatrix_s::rowptr, spm_size_of(), SpmCSC, SpmCSR, SpmPattern, and spmatrix_s::values.

Referenced by spm_scatter_init(), and spmGather().

◆ spmExit()

◆ spmCopy()

void spmCopy ( const spmatrix_t spm,
spmatrix_t newspm 
)

Create a copy of the spm.

Duplicate the spm data structure given as parameter. All new arrays are allocated and copied from the original matrix. Both matrices need to be freed. Procedure version.

Parameters
[in]spmThe sparse matrix to copy.
[out]newspmThe pre-allocated spm to copy into.

Definition at line 969 of file spm.c.

References spmatrix_s::colptr, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::glob2loc, spmatrix_s::gN, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, spmatrix_s::rowptr, spm_size_of(), SpmCSC, SpmCSR, and spmatrix_s::values.

Referenced by c_spmExpand(), d_spmExpand(), p_spmExpand(), s_spmExpand(), spmCheckAndCorrect(), spmDofExtend(), spmGather(), spmScatter(), and z_spmExpand().

◆ spmBase()

void spmBase ( spmatrix_t spm,
int  baseval 
)

Rebase the arrays of the spm to the given value.

If the value is equal to the original base, then nothing is performed.

Parameters
[in,out]spmThe sparse matrix to rebase.
[in]basevalThe new base value to use in the graph (0 or 1).

Definition at line 271 of file spm.c.

References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dofs, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::replicated, spmatrix_s::rowptr, SpmCSC, and SpmCSR.

◆ spmFindBase()

spm_int_t spmFindBase ( const spmatrix_t spm)

Search the base used in the spm structure.

Parameters
[in]spmThe sparse matrix structure.
Return values
Thebaseval used in the given sparse matrix structure.

Definition at line 347 of file spm.c.

References spmatrix_s::clustnbr, spmatrix_s::colptr, spmatrix_s::comm, spmatrix_s::fmttype, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::replicated, spmatrix_s::rowptr, spm_imin(), SPM_MPI_INT, and SpmIJV.

◆ spmConvert()

int spmConvert ( int  ofmttype,
spmatrix_t spm 
)

Convert the storage format of the spm.

Parameters
[in]ofmttypeThe output format of the sparse matrix. It must be:
  • SpmCSC
  • SpmCSR
  • SpmIJV
[in,out]spmThe sparse matrix structure to convert.
Return values
SPM_SUCCESSif the conversion happened successfully.
SPM_ERR_BADPARAMETERif one the parameter is incorrect.
SPM_ERR_NOTIMPLEMENTEDif the case is not yet implemented.

Definition at line 418 of file spm.c.

References spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::replicated, and SPM_SUCCESS.

Referenced by spm_redist_finalize(), and spmCheckAndCorrect().

◆ spmUpdateComputedFields()

◆ spmGenFakeValues()

void spmGenFakeValues ( spmatrix_t spm)

Generate the fake values array such that $ M =  \alpha * D - \beta * A $.

D is the degree matrix, and A the adjacency matrix. The resulting matrix uses real double.

Parameters
[in,out]spmThe spm structure for which the values array must be generated.

Definition at line 510 of file spm_gen_fake_values.c.

References spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::flttype, spmatrix_s::gN, spmatrix_s::n, spmatrix_s::replicated, spm_add_diag(), spm_compute_degrees(), spm_generate_fake_values(), SPM_MPI_INT, SpmPattern, spmUpdateComputedFields(), and spmatrix_s::values.

◆ spmScatter()

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.

Parameters
[in,out]newspmThe pre-allocated spm filled by the scattered spm.
[in]rootThe root process of the scatter operation. -1 if everyone hold a copy of the oldspm.
[in]oldspmThe old sparse matrix to scatter. If spm is a CSC matrix, distribution by row is not possible. If spm is a CSR matrix, distribution by column is not possible. Referenced only on root node(s).
[in]nSize of the loc2glob array if provided. Unused otherwise.
[in]loc2globDistribution array of the matrix. Will be copied. If NULL, the columns/rows are evenly distributed among the processes.
[in]distByColumnBoolean to decide if the matrix is distributed by rows or columns. If false, distribution by rows. If true, distribution by columns.
[in]commMPI communicator.
Return values
SPM_SUCCESSon success,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 1457 of file spm_scatter.c.

References spmatrix_s::clustnbr, spmatrix_s::clustnum, spmatrix_s::comm, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, SPM_ERR_BADPARAMETER, SPM_MPI_INT, spm_scatter_csx(), spm_scatter_ijv(), spm_scatter_init(), SPM_SUCCESS, spmCopy(), SpmCSC, SpmCSR, SpmIJV, and spmUpdateComputedFields().

Referenced by spm_read_driver(), spmLoadDist(), and spmRedistribute().

◆ spmGather()

int spmGather ( const spmatrix_t oldspm,
int  root,
spmatrix_t newspm 
)

Gather a distributed Sparse Matrix on a single node.

The matrix must be distributed continously on each node to be gathered. If it's not the case, please apply first a permutation where the loc2glob array corresponds to the inverse permutation.

Parameters
[in]oldspmThe distributed sparse matrix to gather.
[in]rootThe root node that gather the final sparse matrix. If -1, all nodes gather a copy of the matrix.
[in,out]newspmOn entry, the allocated spm structure if we are on one root of the gather. Not referenced otherwise. On exit, contains the gathered spm on root node(s).
Return values
SPM_SUCCESSon success.
SPM_ERR_BADPARAMETERon incorrect parameter.

Definition at line 436 of file spm_gather.c.

References spmatrix_s::baseval, spmatrix_s::clustnbr, spmatrix_s::clustnum, spmatrix_s::colptr, spmatrix_s::comm, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::fmttype, spmatrix_s::glob2loc, spmatrix_s::gN, spmatrix_s::gNexp, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::nexp, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, spmatrix_s::rowptr, spm_create_loc2glob_continuous(), SPM_ERR_BADPARAMETER, spm_gather_check(), spm_gather_csx(), spm_gather_ijv(), spm_gather_init(), SPM_SUCCESS, spmAlloc(), spmCopy(), SpmCSC, SpmCSR, spmExit(), SpmIJV, spmInit(), spmRedistribute(), and spmatrix_s::values.

Referenced by spm_read_driver(), spmGatherInPlace(), spmRedistribute(), and spmSave().

◆ spmGatherInPlace()

int spmGatherInPlace ( spmatrix_t spm)

This routine performs a allgather of a distributed spm in place.

Parameters
[in]spmOn entry, the distributed spm. On exit, the gathered (replicated) spm.
Return values
1if the spm has been gathered, 0 if untouched Note that untouched is not necessarily a failure if the spm was not distributed.

Definition at line 565 of file spm_gather.c.

References spmatrix_s::clustnbr, spmatrix_s::replicated, SPM_SUCCESS, spmExit(), and spmGather().

◆ spmRedistribute()

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.

Parameters
[in]spmThe original sparse matrix to redistribute. If spm->loc2glob == NULL, the spm is just scattered based on loc2glob distribution to create the new one.
[in]new_nSize of the newl2g array. Must be > 0. Not referenced if newl2g == NULL.
[in]newl2gArray of size new_ with the same base value as the input spm. Contains the distribution array for the new distributed sparse matrix. If NULL, the spm will be gathered to create the new one.
[in,out]newspmOn entry, the allocated spm structure. On exit, the structure contains the redistributed matrix based on the newl2g distribution.
Return values
SPM_SUCCESSon success,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 930 of file spm_redistribute.c.

References spmatrix_s::baseval, spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::fmttype, spmatrix_s::loc2glob, spmatrix_s::replicated, spm_get_distribution(), spm_redist_extract_local(), spm_redist_finalize(), spm_redist_get_newg2l(), spm_redist_post_recv(), spm_redist_send(), SPM_SUCCESS, SpmCSR, spmGather(), and spmScatter().

Referenced by spmGather().

◆ spmNorm()

double spmNorm ( spm_normtype_t  ntype,
const spmatrix_t spm 
)

Compute the norm of the spm.

Return the ntype norm of the sparse matrix spm.

spmNorm = ( max(abs(spm(i,j))), NORM = SpmMaxNorm
          (
          ( norm1(spm),         NORM = SpmOneNorm
          (
          ( normI(spm),         NORM = SpmInfNorm
          (
          ( normF(spm),         NORM = SpmFrobeniusNorm

where norm1 denotes the one norm of a matrix (maximum column sum), normI denotes the infinity norm of a matrix (maximum row sum) and normF denotes the Frobenius norm of a matrix (square root of sum of squares). Note that max(abs(spm(i,j))) is not a consistent matrix norm.

Parameters
[in]ntype
  • SpmMaxNorm
  • SpmOneNorm
  • SpmInfNorm
  • SpmFrobeniusNorm
[in]spmThe sparse matrix structure.
Return values
normThe norm described above. Note that for simplicity, even if the norm of single real or single complex matrix is computed with single precision, the returned norm is stored in double precision.
-1If the floating point of the sparse matrix is undefined.

Definition at line 514 of file spm.c.

References c_spmNorm(), d_spmNorm(), spmatrix_s::flttype, spmatrix_s::replicated, s_spmNorm(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmNorm().

Referenced by c_spmCheckAxb(), d_spmCheckAxb(), s_spmCheckAxb(), and z_spmCheckAxb().

◆ spmNormVec()

double spmNormVec ( spm_normtype_t  ntype,
const spmatrix_t spm,
const void *  x,
spm_int_t  inc 
)

Compute the norm of the spm.

Return the ntype norm of the sparse matrix spm.

spmNormVec = ( max(abs(x(i))), NORM = SpmMaxNorm
             (
             ( norm1(x),       NORM = SpmOneNorm
             (
             ( normI(x),       NORM = SpmInfNorm
             (
             ( normF(x),       NORM = SpmFrobeniusNorm

where norm1 denotes the one norm of a matrix (maximum column sum), normI denotes the infinity norm of a matrix (maximum row sum) and normF denotes the Frobenius norm of a matrix (square root of sum of squares). Note that max(abs(x(i))) is not a consistent matrix norm.

Parameters
[in]ntype
  • SpmMaxNorm
  • SpmOneNorm
  • SpmInfNorm
  • SpmFrobeniusNorm
[in]spmThe sparse matrix structure.
[in]xThe vector for which to compute the norm that is distributed by row as described in the spm. The arithmetic type used is the one described by spm->flttype.
[in]incThe incremental step between each element of the vector.
Return values
normThe norm described above. Note that for simplicity, even if the norm of single real or single complex matrix is computed with single precision, the returned norm is stored in double precision.
-1If the floating point of the sparse matrix is undefined.

Definition at line 596 of file spm.c.

References c_spmNormMat(), d_spmNormMat(), spmatrix_s::flttype, spmatrix_s::nexp, spmatrix_s::replicated, s_spmNormMat(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmNormMat().

◆ spmNormMat()

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.

Return the ntype norm of the sparse matrix spm.

spmNormMat = ( max(abs(A(i,j))), NORM = SpmMaxNorm
             (
             ( norm1(A),         NORM = SpmOneNorm
             (
             ( normI(A),         NORM = SpmInfNorm
             (
             ( normF(A),         NORM = SpmFrobeniusNorm

where norm1 denotes the one norm of a matrix (maximum column sum), normI denotes the infinity norm of a matrix (maximum row sum) and normF denotes the Frobenius norm of a matrix (square root of sum of squares). Note that max(abs(A(i,j))) is not a consistent matrix norm.

Parameters
[in]ntype
  • SpmMaxNorm
  • SpmOneNorm
  • SpmInfNorm
  • SpmFrobeniusNorm
[in]spmThe sparse matrix structure.
[in]nThe number of columns of the matrix A.
[in]AThe matrix A of size lda-by-n.
[in]ldaThe leading dimension of the matrix A. Must be >= max(1, spm->nexp).
Return values
normThe norm described above. Note that for simplicity, even if the norm of single real or single complex matrix is computed with single precision, the returned norm is stored in double precision.
-1If the floating point of the sparse matrix is undefined.

Definition at line 692 of file spm.c.

References c_spmNormMat(), d_spmNormMat(), spmatrix_s::flttype, spmatrix_s::replicated, s_spmNormMat(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmNormMat().

◆ spmMatVec()

int spmMatVec ( spm_trans_t  trans,
double  alpha,
const spmatrix_t spm,
const void *  x,
double  beta,
void *  y 
)

Compute a matrix-vector product.

Performs $ y = alpha * op(A) * x + beta * y $, where $ op(A) $ is one of

$ op( A ) = A  $ or $ op( A ) = A' $ or $ op( A ) = conjg( A' ) $

alpha and beta are scalars, and x and y are vectors.

Parameters
[in]transSpecifies whether the matrix spm is transposed, not transposed or conjugate transposed:
  • SpmTrans
  • SpmNoTrans
  • SpmConjTrans
[in]alphaalpha specifies the scalar alpha.
[in]spmThe SpmGeneral spm.
[in]xThe vector x.
[in]betabeta specifies the scalar beta.
[in,out]yThe vector y.
Return values
SPM_SUCCESSif the y vector has been computed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 1234 of file spm.c.

References spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::replicated, spm_cspmv(), spm_dspmv(), SPM_ERR_BADPARAMETER, spm_sspmv(), SPM_SUCCESS, spm_zspmv(), SpmComplex32, SpmComplex64, SpmCSC, SpmCSR, SpmDouble, spmExit(), SpmFloat, SpmIJV, and SpmPattern.

◆ spmMatMat()

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.

y = alpha * op(A) * B + beta * C

where op(A) is one of:

op( A ) = A or op( A ) = A' or op( A ) = conjg( A' )

alpha and beta are scalars, and x and y are vectors.

Parameters
[in]transSpecifies whether the matrix spm is transposed, not transposed or conjugate transposed:
  • SpmTrans
  • SpmNoTrans
  • SpmConjTrans
[in]nThe number of columns of the matrices B and C.
[in]alphaalpha specifies the scalar alpha.
[in]AThe square sparse matrix A
[in]BThe matrix B of size ldb-by-n
[in]ldbThe leading dimension of the matrix B. ldb >= max( A->nexp, 1 )
[in]betabeta specifies the scalar beta.
[in,out]CThe matrix C of size ldc-by-n
[in]ldcThe leading dimension of the matrix C. ldc >= max( A->nexp, 1 )
Return values
SPM_SUCCESSif the y vector has been computed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 1327 of file spm.c.

References spmatrix_s::flttype, spmatrix_s::nexp, spmatrix_s::replicated, spm_cspmm(), spm_dspmm(), SPM_ERR_BADPARAMETER, spm_imax(), spm_sspmm(), SPM_SUCCESS, spm_zspmm(), SpmComplex32, SpmComplex64, SpmDouble, spmExit(), SpmFloat, SpmLeft, and SpmNoTrans.

◆ spmScal()

void spmScal ( double  alpha,
spmatrix_t spm 
)

Scale the spm.

A = alpha * A

Parameters
[in]alphaThe scaling parameter.
[in,out]spmThe sparse matrix to scal.

Definition at line 1481 of file spm.c.

References c_spmScal(), d_spmScal(), spmatrix_s::flttype, spmatrix_s::replicated, s_spmScal(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmScal().

Referenced by spmScalMatrix().

◆ spmScalVec()

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().

x = alpha * x

Parameters
[in]alphaThe scaling parameter.
[in]spmThe sparse matrix structure associated to the vector to describe the size, the arithmetic used and the distribution of x.
[in,out]xThe local vector to scal of size ( 1 + (spm->nexp-1) * abs(incx) ), and of type defined by spm->flttype.
[in]incxStorage spacing between elements of x. Usually 1, unless corner cases (see blas/lapack).

Definition at line 1545 of file spm.c.

References spmatrix_s::flttype, spmatrix_s::nexp, spmatrix_s::replicated, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and SpmPattern.

◆ 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.

A = alpha * A Rk: We do consider that the matrix A is stored in column major, and distributed by rows similarly to the distribution of the spm if in distributed.

Parameters
[in]alphaThe scaling parameter.
[in]spmThe sparse matrix structure associated to the matrix to describe the size, the arithmetic used and the distribution of A.
[in]nThe number of columns of the matrix A.
[in,out]AThe local matrix A of size LDA -by- n.
[in]ldaThe leading dimension of the matrix A. lda >= max( 1, spm->nexp ).

Definition at line 1602 of file spm.c.

References spmatrix_s::flttype, LAPACKE_clascl_work, LAPACKE_dlascl_work, LAPACKE_slascl_work, LAPACKE_zlascl_work, spmatrix_s::nexp, spmatrix_s::replicated, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and SpmPattern.

◆ spmScalMatrix()

void spmScalMatrix ( double  alpha,
spmatrix_t spm 
)

Scale the spm.

A = alpha * A

Parameters
[in]alphaThe scaling parameter.
[in,out]spmThe sparse matrix to scal.

Deprecated function replaced by spmScal().

Definition at line 1512 of file spm.c.

References spmScal().

◆ spmScalVector()

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.

x = alpha * x

Parameters
[in]fltDatatype of the elements in the vector that must be:
  • SpmFloat
  • SpmDouble
  • SpmComplex32
  • SpmComplex64
[in]nNumber of elements in the input vectors
[in]alphaThe scaling parameter.
[in,out]xThe vector to scal of size ( 1 + (n-1) * abs(incx) ), and of type defined by flt.
[in]incxStorage spacing between elements of x.

Definition at line 1665 of file spm.c.

References SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and SpmPattern.

◆ spmSort()

int spmSort ( spmatrix_t spm)

Sort the subarray of edges of each vertex.

Parameters
[in,out]spmOn entry, the pointer to the sparse matrix structure. On exit, the same sparse matrix with subarrays of edges sorted by ascending order.
Return values
SPM_SUCCESSif the sort was called
SPM_ERR_BADPARAMETERif the given spm was incorrect.

Definition at line 746 of file spm.c.

References c_spmSort(), d_spmSort(), spmatrix_s::flttype, p_spmSort(), spmatrix_s::replicated, s_spmSort(), SPM_ERR_BADPARAMETER, SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmSort().

Referenced by c_spmSortMultidof(), d_spmSortMultidof(), p_spmSortMultidof(), s_spmSortMultidof(), spm_redist_finalize(), spmCheckAndCorrect(), and z_spmSortMultidof().

◆ spmMergeDuplicate()

spm_int_t spmMergeDuplicate ( spmatrix_t spm)

Merge multiple entries in a spm by summing their values together.

The sparse matrix needs to be sorted first (see z_spmSort()). In distributed, only local entries are merged together.

Warning
Not implemented for IJV format.
Parameters
[in,out]spmOn entry, the pointer to the sparse matrix structure. On exit, the reduced sparse matrix of multiple entries were present in it. The multiple values for a same vertex are sum up together.
Return values
>=0the number of vertices that were merged,
SPM_ERR_BADPARAMETERif the given spm was incorrect.

Definition at line 796 of file spm.c.

References c_spmMergeDuplicate(), spmatrix_s::clustnbr, spmatrix_s::comm, d_spmMergeDuplicate(), spmatrix_s::flttype, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::nnz, spmatrix_s::nnzexp, p_spmMergeDuplicate(), spmatrix_s::replicated, s_spmMergeDuplicate(), SPM_ERR_BADPARAMETER, SPM_MPI_INT, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmMergeDuplicate().

Referenced by spmCheckAndCorrect().

◆ spmSymmetrize()

spm_int_t spmSymmetrize ( spmatrix_t spm)

Symmetrize the pattern of the spm.

This routine corrects the sparse matrix structure if it's pattern is not symmetric. It returns the new symmetric pattern with zeroes on the new entries.

First, we look for local structure and store missing entries,and/or entries that should be verified on other nodes if any. Second, if we are in a distributed case, we exhange the entries to check and verifies if they are available locally or not. The missing entries are added to the local buffer of missing entries. Third, we extend the matrix structure with the missing entries, as well as the values array in which 0 values are added.

Parameters
[in,out]spmOn entry, the pointer to the sparse matrix structure. On exit, the same sparse matrix with extra entries that makes it pattern symmetric.
Return values
>=0the number of entries added to the matrix,
SPM_ERR_BADPARAMETERif the given spm was incorrect.

Definition at line 879 of file spm_symmetrize.c.

References spmatrix_s::clustnbr, spmatrix_s::clustnum, spmatrix_s::comm, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, spm_imax(), SPM_MPI_INT, spm_symm_check_local_pattern(), and spm_symm_local_add().

Referenced by spmCheckAndCorrect().

◆ spmCheckAndCorrect()

int spmCheckAndCorrect ( const spmatrix_t spm_in,
spmatrix_t spm_out 
)

Check the correctness of a spm.

This routine initializes the sparse matrix to fit the PaStiX requirements. If needed, the format is changed to CSC, the duplicated vertices are merged together by summing their values; the graph is made symmetric for matrices with unsymmetric pattern, new values are set to 0. Only the lower part is kept for symmetric matrices.

Parameters
[in]spm_inThe pointer to the sparse matrix structure to check, and correct.
[in,out]spm_outOn entry, an allocated structure to hold the corrected spm. On exit, holds the pointer to spm corrected.
Return values
0if no changes have been made to the spm matrix.
1if corrections have been applied and a new spm is returned.

Definition at line 877 of file spm.c.

References spmatrix_s::clustnum, spmatrix_s::fmttype, spmatrix_s::mtxtype, spmatrix_s::replicated, SPM_SUCCESS, spmConvert(), spmCopy(), SpmCSC, spmExit(), SpmGeneral, spmMergeDuplicate(), spmSort(), and spmSymmetrize().

◆ spmGenMat()

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.

Parameters
[in]typeDefines how to compute the vector b.
  • SpmRhsOne: x = 1 [ + I ]
  • SpmRhsI: x = i [ + i * I ]
  • SpmRhsRndX: x is random
  • SpmRhsRndB: x is random
[in]nrhsNumber of columns of the generated vectors
[in]spmThe sparse matrix used to generate the right hand side, and the solution of the full problem.
[in]alphaScaling factor of x.
[in]seedThe seed for the random generator
[out]AThe generated matrix. It has to be preallocated with a size lda -by- nrhs.
[in]ldaDefines the leading dimension of A when multiple right hand sides are available. lda >= spm->nexp.
Return values
SPM_SUCCESSif the b vector has been computed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 1732 of file spm.c.

References c_spmGenMat(), d_spmGenMat(), spmatrix_s::flttype, spmatrix_s::nexp, spmatrix_s::replicated, s_spmGenMat(), SPM_ERR_BADPARAMETER, spm_imax(), SpmFloat, and z_spmGenMat().

Referenced by spmGenVec().

◆ spmGenVec()

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.

Parameters
[in]typeDefines how to compute the vector b.
  • SpmRhsOne: x = 1 [ + I ]
  • SpmRhsI: x = i [ + i * I ]
  • SpmRhsRndX: x is random
  • SpmRhsRndB: x is random
[in]spmThe sparse matrix used to generate the right hand side, and the solution of the full problem.
[in]alphaScaling factor of x.
[in]seedThe seed for the random generator
[out]xThe generated vector. Its size has to be preallocated.
[in]incxDefines the increment of x. Must be superior to 0.
Warning
For the moement, only incx = 1 is supported..
Return values
SPM_SUCCESSif the b vector has been computed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 1802 of file spm.c.

References spmatrix_s::nexp, SPM_ERR_BADPARAMETER, spm_imax(), and spmGenMat().

◆ spmGenRHS()

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.

The vectors can be initialized randomly or to get a specific solution.

Parameters
[in]typeDefines how to compute the vector b.
  • SpmRhsOne: b is computed such that x = 1 [ + I ]
  • SpmRhsI: b is computed such that x = i [ + i * I ]
  • SpmRhsRndX: b is computed by matrix-vector product, such that x is a random vector in the range [-0.5, 0.5]
  • SpmRhsRndB: b is computed randomly and x is not computed.
[in]nrhsDefines the number of right hand side that must be generated.
[in]spmThe sparse matrix used to generate the right hand side, and the solution of the full problem.
[out]xOn exit, if x != NULL, then the x vector(s) generated to compute b is returned. Must be of size at least ldx * nrhs.
[in]ldxDefines the leading dimension of x when multiple right hand sides are available. ldx >= max( 1, spm->nexp ).
[in,out]bb must be an allocated matrix of size at least ldb * nrhs. On exit, b is initialized as defined by the type parameter.
[in]ldbDefines the leading dimension of b when multiple right hand sides are available. ldb >= max( 1, spm->nexp ).
Return values
SPM_SUCCESSif the b vector has been computed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 124 of file spm_rhs.c.

References c_spmGenRHS(), d_spmGenRHS(), spmatrix_s::flttype, spmatrix_s::nexp, s_spmGenRHS(), SPM_ERR_BADPARAMETER, spm_imax(), SpmFloat, and z_spmGenRHS().

◆ spmCheckAxb()

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.

Parameters
[in]epsThe epsilon threshold used for the refinement step. -1. to use the machine precision.
[in]nrhsDefines the number of right hand side that must be generated.
[in]spmThe sparse matrix used to generate the right hand side, and the solution of the full problem.
[in,out]x0If x0 != NULL, the forward error is computed. On exit, x0 stores x0-x
[in]ldx0Defines the leading dimension of x0 when multiple right hand sides are available. ldx0 >= max( 1, spm->nexp ).
[in,out]bb is a matrix of size at least ldb * nrhs. On exit, b stores Ax-b.
[in]ldbDefines the leading dimension of b when multiple right hand sides are available. ldb >= max( 1, spm->nexp ).
[in]xContains the solution computed by the solver.
[in]ldxDefines the leading dimension of x when multiple right hand sides are available. ldx >= max( 1, spm->nexp ).
Return values
SPM_SUCCESSif the tests are succesfull
SPM_ERR_BADPARAMETERif the input matrix is incorrect
1,ifone of the test failed

Definition at line 1423 of file spm.c.

References c_spmCheckAxb(), d_spmCheckAxb(), spmatrix_s::flttype, spmatrix_s::nexp, spmatrix_s::replicated, s_spmCheckAxb(), SPM_ERR_BADPARAMETER, spm_imax(), SpmFloat, and z_spmCheckAxb().

◆ spmExtractLocalRHS()

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.

Parameters
[in]nrhsDefines the number of right hand side.
[in]spmThe sparse matrix used to generate the right hand side.
[in]bglobThe ldbg -by- nrhs global dense matrix
[in]ldbgThe leading dimension of bglob. ldbg >= max( 1, spm->gNexp ).
[in,out]blocThe ldbl -by- nrhs local dense matrix Must be allocated on entry, on exit contains the local values of bglob.
[in]ldblThe leading dimension of bloc. ldbl >= max( 1, spm->nexp ).
Return values
SPM_SUCCESSif the b vector has been Reduceed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 193 of file spm_rhs.c.

References c_spmExtractLocalRHS(), d_spmExtractLocalRHS(), spmatrix_s::flttype, spmatrix_s::gNexp, s_spmExtractLocalRHS(), SPM_ERR_BADPARAMETER, spm_imax(), SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, spmatrix_s::values, and z_spmExtractLocalRHS().

◆ spmReduceRHS()

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.

Parameters
[in]nrhsDefines the number of right hand side.
[in]spmThe sparse matrix used to generate the right hand side.
[in,out]bglobThe ldbg -by- nrhs global dense matrix On entry, the partial accumulation of a matrix vector product. On exit, each node conatins the reduction sum of the matrices.
[in]ldbgThe leading dimension of bglob. ldbg >= max( 1, spm->gNexp ).
[in,out]blocThe ldbl -by- nrhs local dense matrix On entry, must be allocated. On exit, contains the local extraction of the final bglob.
[in]ldblThe leading dimension of bloc. ldbl >= max( 1, spm->nexp ).
Return values
SPM_SUCCESSif the bglob vector has been Reduceed successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 273 of file spm_rhs.c.

References c_spmReduceRHS(), d_spmReduceRHS(), spmatrix_s::flttype, spmatrix_s::gNexp, s_spmReduceRHS(), SPM_ERR_BADPARAMETER, spm_imax(), SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, spmatrix_s::values, and z_spmReduceRHS().

◆ spmGatherRHS()

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.

Parameters
[in]nrhsDefines the number of right hand side.
[in]spmThe sparse matrix used to generate the right hand side.
[in]blocThe ldbl -by- nrhs local dense matrix On entry, the local dense matrix to gather
[in]ldblThe leading dimension of bloc. ldbl >= max( 1, spm->nexp ).
[in]rootClustnum where the complete vector will be gathered. Set it to -1 if you want to gather the global RHS on all nodes.
[in,out]bglobThe ldbg -by- nrhs global dense matrix On entry, the allocated matrix if on root. Non referenced otherwise. On exit, on the root node, contains the global dense matrix.
[in]ldbgThe leading dimension of bglob. ldbg >= max( 1, spm->gNexp ).
Return values
SPM_SUCCESSif the b vector has been gathered successfully,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 356 of file spm_rhs.c.

References c_spmGatherRHS(), spmatrix_s::clustnum, d_spmGatherRHS(), spmatrix_s::flttype, spmatrix_s::gNexp, spmatrix_s::nexp, s_spmGatherRHS(), SPM_ERR_BADPARAMETER, spm_imax(), SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, spmatrix_s::values, and z_spmGatherRHS().

◆ spmIntConvert()

void spmIntConvert ( spm_int_t  n,
const int *  input,
spm_int_t output 
)

Convert integer array to spm_int_t format.

Parameters
[in]nThe number of elements in the array.
[in]inputThe input int array.
[in,out]outputThe output spm_int_t array.

Definition at line 43 of file spm_integers.c.

◆ spmLoadDist()

int spmLoadDist ( spmatrix_t spm,
const char *  filename,
SPM_Comm  comm 
)

Load the spm structure from a file (internal format).

Parameters
[in,out]spmOn entry, an allocated spm data structure. On exit, the spm filled with the information read in the file.
[in]filenameThe filename where the spm is stored. If filename == NULL, "./matrix.spm" is used.
[in]commThe MPI communicator which share the file to load.
Return values
SPM_SUCCESSif the load happened successfully,
SPM_ERR_FILEif the input format is incorrect.

Definition at line 664 of file spm_io.c.

References spmatrix_s::fmttype, spm_load_local(), SPM_SUCCESS, SpmCSR, spmExit(), spmInitDist(), and spmScatter().

Referenced by spmLoad().

◆ spmLoad()

int spmLoad ( spmatrix_t spm,
const char *  filename 
)

Load the spm structure from a file (internal format).

Parameters
[in,out]spmOn entry, an allocated spm data structure. On exit, the spm filled with the information read in the file.
[in]filenameThe filename where the spm is stored. If filename == NULL, "./matrix.spm" is used.
Return values
SPM_SUCCESSif the load happened successfully,
SPM_ERR_FILEif the input format is incorrect.

Definition at line 749 of file spm_io.c.

References spmLoadDist().

Referenced by spm_read_driver().

◆ spmSave()

int spmSave ( const spmatrix_t spm,
const char *  filename 
)

Save the spm structure into a file (internal format).

Parameters
[in]spmThe sparse matrix to write into the file.
[in]filenameThe path where to store the spm. If filename == NULL, the spm is saved into the "./matrix.spm" file.
Return values
SPM_SUCCESSif the save happened successfully.
SPM_ERR_FILEif the input format is incorrect.

Definition at line 1065 of file spm_io.c.

References spmatrix_s::clustnbr, spmatrix_s::clustnum, spmatrix_s::comm, spmatrix_s::replicated, spm_save_local(), spmExit(), and spmGather().

◆ spmReadDriver()

int spmReadDriver ( spm_driver_t  driver,
const char *  filename,
spmatrix_t spm 
)

Import a matrix file into a spm structure.

This function read or generate a sparse matrix from a file to store it into an spm structure. The different formats accepted by this driver are described by the driver field.

Parameters
[in]driverThis defines the driver to use to create the spm structure:
  • SpmDriverRSA
  • SpmDriverHB
  • SpmDriverIJV
  • SpmDriverMM
  • SpmDriverLaplacian
  • SpmDriverXLaplacian
  • SpmDriverGraph
  • SpmDriverSPM
[in]filenameThe name of the file that stores the matrix (see driver).
[in,out]spmOn entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file.
Return values
SPM_SUCCESSif the file reading happened successfully,
SPM_ERR_BADPARAMETERif one the parameter is incorrect.

Definition at line 339 of file spm_read_driver.c.

References spm_read_driver().

◆ spmReadDriverDist()

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.

This function read or generate a sparse matrix from a file to store it into an spm structure. The different formats accepted by this driver are described by the driver field.

Parameters
[in]driverThis defines the driver to use to create the spm structure:
  • SpmDriverRSA
  • SpmDriverHB
  • SpmDriverIJV
  • SpmDriverMM
  • SpmDriverLaplacian
  • SpmDriverXLaplacian
  • SpmDriverGraph
  • SpmDriverSPM
[in]filenameThe name of the file that stores the matrix (see driver).
[in,out]spmOn entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file.
[in]commThe MPI communicator of the problem
Return values
SPM_SUCCESSif the file reading happened successfully,
SPM_ERR_BADPARAMETERif one the parameter is incorrect.

Definition at line 291 of file spm_read_driver.c.

References spm_read_driver().

◆ spmParseLaplacianInfo()

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.

The laplacian will be of size dim1 * dim2 * dim3, and will be equal to $ M = \alpha * D - \beta * A $

where D is the degree matrix, and A the adjacency matrix.

Parameters
[in]filenameConfiguration string of the Laplacian. See laplacian_usage() for more information.
[out]flttypeThe floating type of the elements in the matrix.
[out]dim1The first dimension of the laplacian
[out]dim2The second dimension of the laplacian
[out]dim3The third dimension of the laplacian
[out]alphaThe alpha coefficient for the degree matrix
[out]betaThe beta coefficient for the adjacency matrix
[out]dofThe degree of freedom of each unknown
Return values
SPM_SUCCESSif the matrix has been generated successfully
SPM_ERR_BADPARAMETERif the configuration string is incorrect

Definition at line 96 of file laplacian.c.

References laplacian_usage(), SPM_ERR_BADPARAMETER, spm_imax(), SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and SpmPattern.

Referenced by genExtendedLaplacian(), and genLaplacian().

◆ spm2Dense()

void spm2Dense ( const spmatrix_t spm,
void *  A 
)

Convert the spm matrix into a dense matrix for test purpose.

Remarks
DO NOT USE with large matrices. This is for test purpose only.
Parameters
[in]spmThe sparse matrix structure to convert.
[in,out]AOn entry, an allocated matrix of size spm->gNexp-by-spm->gNexp in the datatype of spm->flttype. On exit, the matrix A is set to the sparse matrix spm.

Definition at line 450 of file spm.c.

References c_spm2dense(), d_spm2dense(), spmatrix_s::flttype, spmatrix_s::replicated, s_spm2dense(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and z_spm2dense().

◆ spmPrint()

void spmPrint ( const spmatrix_t spm,
FILE *  stream 
)

Print an spm matrix into into a given file.

Parameters
[in]spmThe sparse matrix to print.
[in]streamFile to print the spm matrix. stdout, if stream == NULL.

Definition at line 1119 of file spm.c.

References c_spmPrint(), d_spmPrint(), spmatrix_s::flttype, spmatrix_s::replicated, s_spmPrint(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmPrint().

◆ spmPrintRHS()

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.

Parameters
[in]spmThe sparse matrix.
[in]nrhsThe number of columns of x.
[in]xThe set of vectors associated to the spm of size ldx-by-nrhs.
[in]ldxThe local leading dimension of the set of vectors (ldx >= spm->n).
[in]streamFile to print the spm matrix. stdout, if stream == NULL.

Definition at line 45 of file spm_rhs.c.

References c_spmPrintRHS(), d_spmPrintRHS(), spmatrix_s::flttype, s_spmPrintRHS(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmPrintRHS().

◆ spmPrintInfo()

void spmPrintInfo ( const spmatrix_t spm,
FILE *  stream 
)

Print basic informations about the spm matrix into a given stream.

Parameters
[in]spmThe sparse matrix to print.
[in,out]streamStream to print the spm matrix. stdout is used if stream == NULL.

Definition at line 1026 of file spm.c.

References spmatrix_s::clustnbr, spmatrix_s::clustnum, spmatrix_s::comm, spmatrix_s::dof, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::gNexp, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::nexp, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, SpmCSC, SpmGeneral, and SpmPattern.

◆ spmExpand()

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.

Duplicate the spm data structure given as parameter. All new arrays are allocated and copied from the original matrix. Both matrices need to be freed.

Parameters
[in]spm_inThe original non expanded matrix in CSC format used as the template for the muti-dof matrix.
[in,out]spm_outThe output expanded matrix generated from the template.

Definition at line 1168 of file spm.c.

References c_spmExpand(), d_spmExpand(), spmatrix_s::flttype, p_spmExpand(), spmatrix_s::replicated, s_spmExpand(), SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, SpmPattern, and z_spmExpand().

◆ spmDofExtend()

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).

Parameters
[in]spmThe sparse matrix used to generate the new multidof spm.
[in]typeDefines how to generate dofs.
  • 0: Generate a constant dof vector,
  • else: Generate a variable dof vector.
[in]dofThe maximum value for dofs.
[in,out]newspmOn entry, the allocated spm structure, on exit the structure filled with the extended multidof spm.
Return values
SPM_SUCCESSon success,
SPM_ERR_BADPARAMETERotherwise.

Definition at line 51 of file spm_dof_extend.c.

References spmatrix_s::baseval, c_spmDofExtend(), spmatrix_s::clustnum, spmatrix_s::comm, d_spmDofExtend(), spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::gN, s_spmDofExtend(), SPM_ERR_BADPARAMETER, SPM_MPI_INT, SPM_SUCCESS, SpmComplex32, SpmComplex64, spmCopy(), SpmDouble, SpmFloat, SpmPattern, spmUpdateComputedFields(), and z_spmDofExtend().

Referenced by genExtendedLaplacian(), and genLaplacian().

◆ spm_load_local()

static int spm_load_local ( spmatrix_t spm,
const char *  filename 
)
static

Load the spm structure from a file (internal format). One node only.

Parameters
[in,out]spmOn entry, an allocated spm data structure. On exit, the spm filled with the information read in the file.
[in]filenameThe opened file in which the spm is stored. If filename == NULL, matrix.spm is opened.
Return values
SPM_SUCCESSif the spm was load successfully.
SPM_ERR_FILEif the input format is incorrect.

Definition at line 464 of file spm_io.c.

References spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::layout, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::nnzexp, readArrayOfComplex32(), readArrayOfComplex64(), readArrayOfDouble(), readArrayOfFloat(), readArrayOfInteger(), spmatrix_s::replicated, spmatrix_s::rowptr, SPM_ERR_FILE, spm_size_of(), SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmCSC, SpmCSR, SpmDouble, spmExit(), SpmFloat, SpmIJV, SpmPattern, spmUpdateComputedFields(), and spmatrix_s::values.

Referenced by spmLoadDist().

◆ spm_save_local()

static int spm_save_local ( const spmatrix_t spm,
const char *  filename 
)
static

Save the global spm structure into a file (internal format).

Parameters
[in]spmThe sparse matrix to write into the file.
[in]filenameThe filename in which to store the spm. If filename == NULL, data is saved into matrix.spm file.
Return values
SPM_SUCCESSif the save happened successfully.
SPM_ERR_FILEif the input format is incorrect.

Definition at line 942 of file spm_io.c.

References spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::gN, spmatrix_s::layout, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::rowptr, SPM_ERR_FILE, SPM_SUCCESS, SpmComplex32, SpmComplex64, SpmCSC, SpmCSR, SpmDouble, SpmFloat, SpmIJV, SpmPattern, spmatrix_s::values, writeArrayOfComplex32(), writeArrayOfComplex64(), writeArrayOfDouble(), and writeArrayOfFloat().

Referenced by spmSave().

◆ spm_read_driver()

static int spm_read_driver ( int  scatter,
spm_driver_t  driver,
const char *  filename,
spmatrix_t spm,
SPM_Comm  comm 
)
static

Import a matrix file into a spm structure for a specific communicator.

This function read or generate a sparse matrix from a file to store it into a spm structure. The different formats accepted by this driver are described by the driver field.

Parameters
[in]scatterBoolean to specify if the final spm must be scattered or not.
[in]driverThis defines the driver to use to create the spm structure:
  • SpmDriverRSA
  • SpmDriverHB
  • SpmDriverIJV
  • SpmDriverMM
  • SpmDriverLaplacian
  • SpmDriverXLaplacian
  • SpmDriverGraph
  • SpmDriverSPM
[in]filenameThe name of the file that stores the matrix (see driver).
[in,out]spmOn entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file.
[in]commThe MPI communicator of the problem
Return values
SPM_SUCCESSif the file reading happened successfully,
SPM_ERR_BADPARAMETERif one the parameter is incorrect.

Definition at line 139 of file spm_read_driver.c.

References spmatrix_s::clustnbr, spmatrix_s::comm, genExtendedLaplacian(), genLaplacian(), spmatrix_s::replicated, SPM_ERR_BADPARAMETER, SPM_ERR_UNKNOWN, spm_read_scotch(), SPM_SUCCESS, SpmDriverGraph, SpmDriverHB, SpmDriverIJV, SpmDriverLaplacian, SpmDriverMM, SpmDriverRSA, SpmDriverSPM, SpmDriverXLaplacian, spmExit(), spmGather(), spmInitDist(), spmLoad(), and spmScatter().

Referenced by spmReadDriver(), and spmReadDriverDist().