SpM Handbook 1.2.4
|
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 ![]() | |
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). | |
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.
struct spmatrix_s |
The sparse matrix data structure.
This structure describes matrices with different characteristics that can be useful to any solver:
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.
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 |
#define SPM_MPI_INT MPI_INT |
The MPI type associated to spm_int_t.
Definition at line 72 of file datatypes.h.
#define SPM_INT_MAX INT_MAX |
The maximum spm_int_t value.
Definition at line 73 of file datatypes.h.
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:
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.
The main integer datatype used in spm arrays.
Definition at line 70 of file datatypes.h.
The main unsigned integer datatype used in spm arrays.
Definition at line 71 of file datatypes.h.
void spmInit | ( | spmatrix_t * | spm | ) |
Init the spm structure.
[in,out] | spm | The 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().
void spmInitDist | ( | spmatrix_t * | spm, |
SPM_Comm | comm | ||
) |
Init the spm structure with a specific communicator.
[in,out] | spm | The sparse matrix to init. |
[in] | comm | The MPI communicator used for the sparse matrix. MPI_COMM_WORLD by default. |
Definition at line 101 of file spm.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::flttype, spmatrix_s::fmttype, spmatrix_s::glob2loc, spmatrix_s::gN, spmatrix_s::gNexp, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::layout, spmatrix_s::loc2glob, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::nexp, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, spmatrix_s::rowptr, SpmColMajor, SpmCSC, SpmDouble, SpmGeneral, and spmatrix_s::values.
Referenced by spm_read_driver(), spm_redist_extract_local(), spm_redist_get_newg2l(), spmInit(), and spmLoadDist().
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.
[in,out] | spm | The 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().
void spmExit | ( | spmatrix_t * | spm | ) |
Cleanup the spm structure but do not free the spm pointer.
[in,out] | spm | The sparse matrix to free. |
Definition at line 224 of file spm.c.
References spmatrix_s::colptr, spmatrix_s::dofs, spmatrix_s::glob2loc, spmatrix_s::loc2glob, spmatrix_s::replicated, spmatrix_s::rowptr, and spmatrix_s::values.
Referenced by c_spmConvertCSR2CSC_gen(), d_spmConvertCSR2CSC_gen(), genExtendedLaplacian(), genLaplacian(), p_spmConvertCSR2CSC_gen(), s_spmConvertCSR2CSC_gen(), spm_load_local(), spm_read_driver(), spmCheckAndCorrect(), spmGather(), spmGatherInPlace(), spmLoadDist(), spmMatMat(), spmMatVec(), spmSave(), and z_spmConvertCSR2CSC_gen().
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.
[in] | spm | The sparse matrix to copy. |
[out] | newspm | The 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().
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.
[in,out] | spm | The sparse matrix to rebase. |
[in] | baseval | The 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.
spm_int_t spmFindBase | ( | const spmatrix_t * | spm | ) |
Search the base used in the spm structure.
[in] | spm | The sparse matrix structure. |
The | baseval 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.
int spmConvert | ( | int | ofmttype, |
spmatrix_t * | spm | ||
) |
Convert the storage format of the spm.
[in] | ofmttype | The output format of the sparse matrix. It must be:
|
[in,out] | spm | The sparse matrix structure to convert. |
SPM_SUCCESS | if the conversion happened successfully. |
SPM_ERR_BADPARAMETER | if one the parameter is incorrect. |
SPM_ERR_NOTIMPLEMENTED | if 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().
void spmUpdateComputedFields | ( | spmatrix_t * | spm | ) |
Update all the computed fields based on the static values stored.
[in,out] | spm | The sparse matrix to update. |
Definition at line 228 of file spm_update_compute_fields.c.
References spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::dof, spmatrix_s::gN, spmatrix_s::gNexp, spmatrix_s::gnnz, spmatrix_s::gnnzexp, spmatrix_s::n, spmatrix_s::nexp, spmatrix_s::nnz, spmatrix_s::nnzexp, spmatrix_s::replicated, SPM_MPI_INT, spm_ucf_variadic_mpi(), and spm_ucf_variadic_shm().
Referenced by c_spmExpand(), d_spmExpand(), genExtendedLaplacian(), genLaplacian(), p_spmExpand(), s_spmExpand(), spm_getandset_glob2loc(), spm_load_local(), spm_read_scotch(), spm_redist_finalize(), spmDofExtend(), spmGenFakeValues(), spmScatter(), and z_spmExpand().
void spmGenFakeValues | ( | spmatrix_t * | spm | ) |
Generate the fake values array such that
D is the degree matrix, and A the adjacency matrix. The resulting matrix uses real double.
[in,out] | spm | The 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.
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.
[in,out] | newspm | The pre-allocated spm filled by the scattered spm. |
[in] | root | The root process of the scatter operation. -1 if everyone hold a copy of the oldspm. |
[in] | oldspm | The 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] | n | Size of the loc2glob array if provided. Unused otherwise. |
[in] | loc2glob | Distribution array of the matrix. Will be copied. If NULL, the columns/rows are evenly distributed among the processes. |
[in] | distByColumn | Boolean to decide if the matrix is distributed by rows or columns. If false, distribution by rows. If true, distribution by columns. |
[in] | comm | MPI communicator. |
SPM_SUCCESS | on success, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | oldspm | The distributed sparse matrix to gather. |
[in] | root | The root node that gather the final sparse matrix. If -1, all nodes gather a copy of the matrix. |
[in,out] | newspm | On 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). |
SPM_SUCCESS | on success. |
SPM_ERR_BADPARAMETER | on 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().
int spmGatherInPlace | ( | spmatrix_t * | spm | ) |
This routine performs a allgather of a distributed spm in place.
[in] | spm | On entry, the distributed spm. On exit, the gathered (replicated) spm. |
1 | if 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().
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.
[in] | spm | The 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_n | Size of the newl2g array. Must be > 0. Not referenced if newl2g == NULL. |
[in] | newl2g | Array 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] | newspm | On entry, the allocated spm structure. On exit, the structure contains the redistributed matrix based on the newl2g distribution. |
SPM_SUCCESS | on success, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | ntype |
|
[in] | spm | The sparse matrix structure. |
norm | The 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. |
-1 | If 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().
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.
[in] | ntype |
|
[in] | spm | The sparse matrix structure. |
[in] | x | The 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] | inc | The incremental step between each element of the vector. |
norm | The 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. |
-1 | If 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().
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.
[in] | ntype |
|
[in] | spm | The sparse matrix structure. |
[in] | n | The number of columns of the matrix A. |
[in] | A | The matrix A of size lda-by-n. |
[in] | lda | The leading dimension of the matrix A. Must be >= max(1, spm->nexp). |
norm | The 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. |
-1 | If 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().
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
alpha and beta are scalars, and x and y are vectors.
[in] | trans | Specifies whether the matrix spm is transposed, not transposed or conjugate transposed:
|
[in] | alpha | alpha specifies the scalar alpha. |
[in] | spm | The SpmGeneral spm. |
[in] | x | The vector x. |
[in] | beta | beta specifies the scalar beta. |
[in,out] | y | The vector y. |
SPM_SUCCESS | if the y vector has been computed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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.
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.
[in] | trans | Specifies whether the matrix spm is transposed, not transposed or conjugate transposed:
|
[in] | n | The number of columns of the matrices B and C. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | The square sparse matrix A |
[in] | B | The matrix B of size ldb-by-n |
[in] | ldb | The leading dimension of the matrix B. ldb >= max( A->nexp, 1 ) |
[in] | beta | beta specifies the scalar beta. |
[in,out] | C | The matrix C of size ldc-by-n |
[in] | ldc | The leading dimension of the matrix C. ldc >= max( A->nexp, 1 ) |
SPM_SUCCESS | if the y vector has been computed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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.
void spmScal | ( | double | alpha, |
spmatrix_t * | spm | ||
) |
Scale the spm.
A = alpha * A
[in] | alpha | The scaling parameter. |
[in,out] | spm | The 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().
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
[in] | alpha | The scaling parameter. |
[in] | spm | The sparse matrix structure associated to the vector to describe the size, the arithmetic used and the distribution of x. |
[in,out] | x | The local vector to scal of size ( 1 + (spm->nexp-1) * abs(incx) ), and of type defined by spm->flttype. |
[in] | incx | Storage 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.
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.
[in] | alpha | The scaling parameter. |
[in] | spm | The sparse matrix structure associated to the matrix to describe the size, the arithmetic used and the distribution of A. |
[in] | n | The number of columns of the matrix A. |
[in,out] | A | The local matrix A of size LDA -by- n. |
[in] | lda | The 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.
void spmScalMatrix | ( | double | alpha, |
spmatrix_t * | 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.
x = alpha * x
[in] | flt | Datatype of the elements in the vector that must be:
|
[in] | n | Number of elements in the input vectors |
[in] | alpha | The scaling parameter. |
[in,out] | x | The vector to scal of size ( 1 + (n-1) * abs(incx) ), and of type defined by flt. |
[in] | incx | Storage spacing between elements of x. |
Definition at line 1665 of file spm.c.
References SpmComplex32, SpmComplex64, SpmDouble, SpmFloat, and SpmPattern.
int spmSort | ( | spmatrix_t * | spm | ) |
Sort the subarray of edges of each vertex.
[in,out] | spm | On entry, the pointer to the sparse matrix structure. On exit, the same sparse matrix with subarrays of edges sorted by ascending order. |
SPM_SUCCESS | if the sort was called |
SPM_ERR_BADPARAMETER | if 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().
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.
[in,out] | spm | On 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. |
>=0 | the number of vertices that were merged, |
SPM_ERR_BADPARAMETER | if 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().
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.
[in,out] | spm | On entry, the pointer to the sparse matrix structure. On exit, the same sparse matrix with extra entries that makes it pattern symmetric. |
>=0 | the number of entries added to the matrix, |
SPM_ERR_BADPARAMETER | if 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().
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.
[in] | spm_in | The pointer to the sparse matrix structure to check, and correct. |
[in,out] | spm_out | On entry, an allocated structure to hold the corrected spm. On exit, holds the pointer to spm corrected. |
0 | if no changes have been made to the spm matrix. |
1 | if 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().
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.
[in] | type | Defines how to compute the vector b.
|
[in] | nrhs | Number of columns of the generated vectors |
[in] | spm | The sparse matrix used to generate the right hand side, and the solution of the full problem. |
[in] | alpha | Scaling factor of x. |
[in] | seed | The seed for the random generator |
[out] | A | The generated matrix. It has to be preallocated with a size lda -by- nrhs. |
[in] | lda | Defines the leading dimension of A when multiple right hand sides are available. lda >= spm->nexp. |
SPM_SUCCESS | if the b vector has been computed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | type | Defines how to compute the vector b.
|
[in] | spm | The sparse matrix used to generate the right hand side, and the solution of the full problem. |
[in] | alpha | Scaling factor of x. |
[in] | seed | The seed for the random generator |
[out] | x | The generated vector. Its size has to be preallocated. |
[in] | incx | Defines the increment of x. Must be superior to 0. |
SPM_SUCCESS | if the b vector has been computed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
Definition at line 1802 of file spm.c.
References spmatrix_s::nexp, SPM_ERR_BADPARAMETER, spm_imax(), and spmGenMat().
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.
[in] | type | Defines how to compute the vector b.
|
[in] | nrhs | Defines the number of right hand side that must be generated. |
[in] | spm | The sparse matrix used to generate the right hand side, and the solution of the full problem. |
[out] | x | On exit, if x != NULL, then the x vector(s) generated to compute b is returned. Must be of size at least ldx * nrhs. |
[in] | ldx | Defines the leading dimension of x when multiple right hand sides are available. ldx >= max( 1, spm->nexp ). |
[in,out] | b | b must be an allocated matrix of size at least ldb * nrhs. On exit, b is initialized as defined by the type parameter. |
[in] | ldb | Defines the leading dimension of b when multiple right hand sides are available. ldb >= max( 1, spm->nexp ). |
SPM_SUCCESS | if the b vector has been computed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | eps | The epsilon threshold used for the refinement step. -1. to use the machine precision. |
[in] | nrhs | Defines the number of right hand side that must be generated. |
[in] | spm | The sparse matrix used to generate the right hand side, and the solution of the full problem. |
[in,out] | x0 | If x0 != NULL, the forward error is computed. On exit, x0 stores x0-x |
[in] | ldx0 | Defines the leading dimension of x0 when multiple right hand sides are available. ldx0 >= max( 1, spm->nexp ). |
[in,out] | b | b is a matrix of size at least ldb * nrhs. On exit, b stores Ax-b. |
[in] | ldb | Defines the leading dimension of b when multiple right hand sides are available. ldb >= max( 1, spm->nexp ). |
[in] | x | Contains the solution computed by the solver. |
[in] | ldx | Defines the leading dimension of x when multiple right hand sides are available. ldx >= max( 1, spm->nexp ). |
SPM_SUCCESS | if the tests are succesfull |
SPM_ERR_BADPARAMETER | if the input matrix is incorrect |
1,if | one 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().
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.
[in] | nrhs | Defines the number of right hand side. |
[in] | spm | The sparse matrix used to generate the right hand side. |
[in] | bglob | The ldbg -by- nrhs global dense matrix |
[in] | ldbg | The leading dimension of bglob. ldbg >= max( 1, spm->gNexp ). |
[in,out] | bloc | The ldbl -by- nrhs local dense matrix Must be allocated on entry, on exit contains the local values of bglob. |
[in] | ldbl | The leading dimension of bloc. ldbl >= max( 1, spm->nexp ). |
SPM_SUCCESS | if the b vector has been Reduceed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | nrhs | Defines the number of right hand side. |
[in] | spm | The sparse matrix used to generate the right hand side. |
[in,out] | bglob | The 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] | ldbg | The leading dimension of bglob. ldbg >= max( 1, spm->gNexp ). |
[in,out] | bloc | The ldbl -by- nrhs local dense matrix On entry, must be allocated. On exit, contains the local extraction of the final bglob. |
[in] | ldbl | The leading dimension of bloc. ldbl >= max( 1, spm->nexp ). |
SPM_SUCCESS | if the bglob vector has been Reduceed successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
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.
[in] | nrhs | Defines the number of right hand side. |
[in] | spm | The sparse matrix used to generate the right hand side. |
[in] | bloc | The ldbl -by- nrhs local dense matrix On entry, the local dense matrix to gather |
[in] | ldbl | The leading dimension of bloc. ldbl >= max( 1, spm->nexp ). |
[in] | root | Clustnum where the complete vector will be gathered. Set it to -1 if you want to gather the global RHS on all nodes. |
[in,out] | bglob | The 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] | ldbg | The leading dimension of bglob. ldbg >= max( 1, spm->gNexp ). |
SPM_SUCCESS | if the b vector has been gathered successfully, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
Convert integer array to spm_int_t format.
[in] | n | The number of elements in the array. |
[in] | input | The input int array. |
[in,out] | output | The output spm_int_t array. |
Definition at line 43 of file spm_integers.c.
int spmLoadDist | ( | spmatrix_t * | spm, |
const char * | filename, | ||
SPM_Comm | comm | ||
) |
Load the spm structure from a file (internal format).
[in,out] | spm | On entry, an allocated spm data structure. On exit, the spm filled with the information read in the file. |
[in] | filename | The filename where the spm is stored. If filename == NULL, "./matrix.spm" is used. |
[in] | comm | The MPI communicator which share the file to load. |
SPM_SUCCESS | if the load happened successfully, |
SPM_ERR_FILE | if 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().
int spmLoad | ( | spmatrix_t * | spm, |
const char * | filename | ||
) |
Load the spm structure from a file (internal format).
[in,out] | spm | On entry, an allocated spm data structure. On exit, the spm filled with the information read in the file. |
[in] | filename | The filename where the spm is stored. If filename == NULL, "./matrix.spm" is used. |
SPM_SUCCESS | if the load happened successfully, |
SPM_ERR_FILE | if the input format is incorrect. |
Definition at line 749 of file spm_io.c.
References spmLoadDist().
Referenced by spm_read_driver().
int spmSave | ( | const spmatrix_t * | spm, |
const char * | filename | ||
) |
Save the spm structure into a file (internal format).
[in] | spm | The sparse matrix to write into the file. |
[in] | filename | The path where to store the spm. If filename == NULL, the spm is saved into the "./matrix.spm" file. |
SPM_SUCCESS | if the save happened successfully. |
SPM_ERR_FILE | if 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().
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.
[in] | driver | This defines the driver to use to create the spm structure:
|
[in] | filename | The name of the file that stores the matrix (see driver). |
[in,out] | spm | On entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file. |
SPM_SUCCESS | if the file reading happened successfully, |
SPM_ERR_BADPARAMETER | if one the parameter is incorrect. |
Definition at line 339 of file spm_read_driver.c.
References spm_read_driver().
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.
[in] | driver | This defines the driver to use to create the spm structure:
|
[in] | filename | The name of the file that stores the matrix (see driver). |
[in,out] | spm | On entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file. |
[in] | comm | The MPI communicator of the problem |
SPM_SUCCESS | if the file reading happened successfully, |
SPM_ERR_BADPARAMETER | if one the parameter is incorrect. |
Definition at line 291 of file spm_read_driver.c.
References spm_read_driver().
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
where D is the degree matrix, and A the adjacency matrix.
[in] | filename | Configuration string of the Laplacian. See laplacian_usage() for more information. |
[out] | flttype | The floating type of the elements in the matrix. |
[out] | dim1 | The first dimension of the laplacian |
[out] | dim2 | The second dimension of the laplacian |
[out] | dim3 | The third dimension of the laplacian |
[out] | alpha | The alpha coefficient for the degree matrix |
[out] | beta | The beta coefficient for the adjacency matrix |
[out] | dof | The degree of freedom of each unknown |
SPM_SUCCESS | if the matrix has been generated successfully |
SPM_ERR_BADPARAMETER | if 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().
void spm2Dense | ( | const spmatrix_t * | spm, |
void * | A | ||
) |
Convert the spm matrix into a dense matrix for test purpose.
[in] | spm | The sparse matrix structure to convert. |
[in,out] | A | On 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().
void spmPrint | ( | const spmatrix_t * | spm, |
FILE * | stream | ||
) |
Print an spm matrix into into a given file.
[in] | spm | The sparse matrix to print. |
[in] | stream | File 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().
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.
[in] | spm | The sparse matrix. |
[in] | nrhs | The number of columns of x. |
[in] | x | The set of vectors associated to the spm of size ldx-by-nrhs. |
[in] | ldx | The local leading dimension of the set of vectors (ldx >= spm->n). |
[in] | stream | File 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().
void spmPrintInfo | ( | const spmatrix_t * | spm, |
FILE * | stream | ||
) |
Print basic informations about the spm matrix into a given stream.
[in] | spm | The sparse matrix to print. |
[in,out] | stream | Stream 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.
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.
[in] | spm_in | The original non expanded matrix in CSC format used as the template for the muti-dof matrix. |
[in,out] | spm_out | The 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().
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).
[in] | spm | The sparse matrix used to generate the new multidof spm. |
[in] | type | Defines how to generate dofs.
|
[in] | dof | The maximum value for dofs. |
[in,out] | newspm | On entry, the allocated spm structure, on exit the structure filled with the extended multidof spm. |
SPM_SUCCESS | on success, |
SPM_ERR_BADPARAMETER | otherwise. |
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().
|
static |
Load the spm structure from a file (internal format). One node only.
[in,out] | spm | On entry, an allocated spm data structure. On exit, the spm filled with the information read in the file. |
[in] | filename | The opened file in which the spm is stored. If filename == NULL, matrix.spm is opened. |
SPM_SUCCESS | if the spm was load successfully. |
SPM_ERR_FILE | if 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().
|
static |
Save the global spm structure into a file (internal format).
[in] | spm | The sparse matrix to write into the file. |
[in] | filename | The filename in which to store the spm. If filename == NULL, data is saved into matrix.spm file. |
SPM_SUCCESS | if the save happened successfully. |
SPM_ERR_FILE | if 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().
|
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.
[in] | scatter | Boolean to specify if the final spm must be scattered or not. |
[in] | driver | This defines the driver to use to create the spm structure:
|
[in] | filename | The name of the file that stores the matrix (see driver). |
[in,out] | spm | On entry, an allocated sparse matrix structure. On exit, the filled sparse matrix structure with the matrix from the file. |
[in] | comm | The MPI communicator of the problem |
SPM_SUCCESS | if the file reading happened successfully, |
SPM_ERR_BADPARAMETER | if 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().