GHOST  1.1.2
General, Hybrid, and Optimized Sparse Toolkit
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
Classes | Macros | Typedefs | Enumerations | Functions | Variables
sparsemat.h File Reference

Types and functions related to sparse matrices. More...

#include "config.h"
#include "types.h"
#include "spmv.h"
#include "context.h"
#include "densemat.h"
#include "sparsemat_src.h"
#include <stdarg.h>
Include dependency graph for sparsemat.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  ghost_sorting_helper
 Helper for sparse matrix row sorting. More...
 
struct  ghost_sparsemat_rowfunc_crs_arg
 
struct  ghost_spmv_opts
 
struct  ghost_kacz_opts
 
struct  ghost_carp_opts
 
struct  ghost_sellspmv_parameters
 The parameters to identify a SELL SpMV kernel. More...
 
struct  ghost_cusellspmv_parameters
 The parameters to identify a CUDA SELL SpMV kernel. More...
 
struct  ghost_kacz_parameters
 The parameters to identify a SELL Kaczmarz kernel. More...
 
struct  ghost_sparsemat_traits
 Sparse matrix traits. More...
 
struct  ghost_sparsemat
 A sparse matrix. More...
 
struct  ghost_sparsemat_rowfunc_file_initargs
 

Macros

#define GHOST_SPARSEMAT_SORT_GLOBAL   -1
 
#define GHOST_SPARSEMAT_SORT_LOCAL   -2
 
#define GHOST_SELL_CHUNKHEIGHT_ELLPACK   0
 Create only a single chunk, i.e., use the ELLPACK storage format. More...
 
#define GHOST_SELL_CHUNKHEIGHT_AUTO   -1
 A chunkheight should automatically be determined. More...
 
#define GHOST_SPARSEMAT_PERM_ANY   (GHOST_SPARSEMAT_PERM_ANY_LOCAL|GHOST_SPARSEMAT_PERM_ANY_GLOBAL)
 Combination of flags which apply any permutation to a ghost_sparsemat. More...
 
#define GHOST_SPARSEMAT_PERM_ANY_LOCAL   (GHOST_SPARSEMAT_COLOR|GHOST_SPARSEMAT_RCM|GHOST_SPARSEMAT_BLOCKCOLOR|GHOST_SPARSEMAT_SORT_ROWS)
 
#define GHOST_SPARSEMAT_PERM_ANY_GLOBAL   (GHOST_SPARSEMAT_SCOTCHIFY|GHOST_SPARSEMAT_ZOLTAN)
 
#define GHOST_SCOTCH_STRAT_DEFAULT   "g"
 
#define SPM_NROWS(mat)   mat->context->row_map->dim
 
#define SPM_NNZ(mat)   mat->context->nnz
 
#define SPM_NCOLS(mat)   mat->context->col_map->dim
 
#define SPM_GNCOLS(mat)   mat->context->col_map->gdim
 
#define SPM_NROWSPAD(mat)   mat->context->row_map->dimpad
 
#define SPM_NCHUNKS(mat)   (mat->nchunks)
 

Typedefs

typedef struct
ghost_sparsemat_traits 
ghost_sparsemat_traits
 
typedef struct ghost_sparsemat ghost_sparsemat
 
typedef ghost_error(* ghost_spmv_kernel )(ghost_densemat *, ghost_sparsemat *, ghost_densemat *, ghost_spmv_opts)
 
typedef ghost_error(* ghost_kacz_kernel )(ghost_densemat *, ghost_sparsemat *, ghost_densemat *, ghost_kacz_opts)
 
typedef ghost_error(* ghost_kacz_shift_kernel )(ghost_densemat *, ghost_densemat *, ghost_sparsemat *, ghost_densemat *, double, double, ghost_kacz_opts)
 

Enumerations

enum  ghost_sparsemat_symmetry { GHOST_SPARSEMAT_SYMM_GENERAL = 1, GHOST_SPARSEMAT_SYMM_SYMMETRIC = 2, GHOST_SPARSEMAT_SYMM_SKEW_SYMMETRIC = 4, GHOST_SPARSEMAT_SYMM_HERMITIAN = 8 }
 Possible sparse matrix symmetries. More...
 
enum  ghost_kacz_direction { GHOST_KACZ_DIRECTION_UNDEFINED, GHOST_KACZ_DIRECTION_FORWARD, GHOST_KACZ_DIRECTION_BACKWARD }
 
enum  ghost_kacz_normalize { GHOST_KACZ_NORMALIZE_YES = 0, GHOST_KACZ_NORMALIZE_NO = 1<<0, GHOST_KACZ_NORMALIZE_SAVE = 1<<1 }
 
enum  ghost_kacz_mode { GHOST_KACZ_MODE_NORMAL, GHOST_KACZ_MODE_ECO, GHOST_KACZ_MODE_PERFORMANCE }
 
enum  ghost_sparsemat_flags {
  GHOST_SPARSEMAT_DEFAULT = 0, GHOST_SPARSEMAT_HOST = 1<<0, GHOST_SPARSEMAT_DEVICE = 1<<1, GHOST_SPARSEMAT_NOT_PERMUTE_COLS = 1<<2,
  GHOST_SPARSEMAT_PERMUTE = 1<<3, GHOST_SPARSEMAT_NOT_SORT_COLS = 1<<4, GHOST_SPARSEMAT_NOT_STORE_SPLIT = 1<<5, GHOST_SPARSEMAT_NOT_STORE_FULL = 1<<6,
  GHOST_SPARSEMAT_SCOTCHIFY = 1<<7, GHOST_SPARSEMAT_SAVE_ORIG_COLS = 1<<8, GHOST_SPARSEMAT_COLOR = 1<<9, GHOST_SPARSEMAT_TRANSPOSE_MM = 1<<10,
  GHOST_SPARSEMAT_ZOLTAN = 1<<11, GHOST_SPARSEMAT_RCM = 1<<12, GHOST_SPARSEMAT_BLOCKCOLOR = 1<<13, GHOST_SOLVER_KACZ = 1<<14,
  GHOST_SPARSEMAT_SORT_ROWS = 1<<15, GHOST_SPARSEMAT_PERM_NO_DISTINCTION =1<<16, GHOST_SPARSEMAT_DIAG_FIRST =1<<17
}
 Flags to a sparse matrix. More...
 

Functions

ghost_error ghost_sparsemat_create (ghost_sparsemat **mat, ghost_context *ctx, ghost_sparsemat_traits *traits, int nTraits)
 Create a sparse matrix. More...
 
ghost_error ghost_sparsemat_info_string (char **str, ghost_sparsemat *matrix)
 Create a string holding information about the sparsemat. More...
 
ghost_error ghost_sparsemat_nnz (ghost_gidx *nnz, ghost_sparsemat *mat)
 Obtain the global number of nonzero elements of a sparse matrix. More...
 
ghost_error ghost_sparsemat_nrows (ghost_gidx *nrows, ghost_sparsemat *mat)
 Obtain the global number of rows of a sparse matrix. More...
 
const char * ghost_sparsemat_symmetry_string (ghost_sparsemat_symmetry symmetry)
 Convert the matrix' symmetry information to a string. More...
 
bool ghost_sparsemat_symmetry_valid (ghost_sparsemat_symmetry symmetry)
 Check if the symmetry information of a sparse matrix is valid. More...
 
ghost_error ghost_sparsemat_perm_scotch (ghost_context *ctx, ghost_sparsemat *mat)
 Create a matrix permutation based on (PT-)SCOTCH. More...
 
ghost_error ghost_sparsemat_perm_sort (ghost_context *ctx, ghost_sparsemat *mat, ghost_lidx scope)
 Create a matrix permutation based on row length sorting within a given scope. More...
 
ghost_error ghost_sparsemat_perm_spmp (ghost_context *ctx, ghost_sparsemat *mat)
 
ghost_error ghost_sparsemat_perm_color (ghost_context *ctx, ghost_sparsemat *mat)
 Create a matrix permutation based on 2-way coloring using ColPack. More...
 
ghost_error ghost_sparsemat_blockColor (ghost_context *ctx, ghost_sparsemat *mat)
 
ghost_error ghost_sparsemat_perm_zoltan (ghost_context *ctx, ghost_sparsemat *mat)
 
ghost_error ghost_sparsemat_fromfile_common (ghost_sparsemat *mat, char *matrixPath, ghost_lidx **rpt)
 Common function for matrix creation from a file. More...
 
ghost_error ghost_sparsematofile_header (ghost_sparsemat *mat, char *path)
 Write the sparse matrix header to a binary CRS file. More...
 
ghost_error ghost_sparsemat_registerrow (ghost_sparsemat *mat, ghost_gidx row, ghost_gidx *col, ghost_lidx rowlen, ghost_lidx stride)
 Store matrix information like bandwidth and nonzero distribution for a given matrix row. More...
 
ghost_error ghost_sparsemat_registerrow_finalize (ghost_sparsemat *mat)
 Finalize the storing of matrix information like bandwidth and nonzero distribution. More...
 
void ghost_sparsemat_destroy (ghost_sparsemat *mat)
 Destroy a sparsemat and free all memory. More...
 
ghost_error ghost_sell_spmv_selector (ghost_densemat *lhs, ghost_sparsemat *mat, ghost_densemat *rhs, ghost_spmv_opts traits)
 Select and call the right SELL SpMV kernel. More...
 
ghost_error ghost_cu_sell_spmv_selector (ghost_densemat *lhs, ghost_sparsemat *mat, ghost_densemat *rhs, ghost_spmv_opts traits)
 Select and call the right CUDA SELL SpMV kernel. More...
 
ghost_error ghost_sell_stringify_selector (ghost_sparsemat *mat, char **str, int dense)
 Select and call the right SELL stringification function. More...
 
ghost_error ghost_sparsemat_string (char **str, ghost_sparsemat *mat, int dense)
 Creates a string of the sparsemat's contents. More...
 
ghost_error ghost_sparsemat_init_rowfunc (ghost_sparsemat *mat, ghost_sparsemat_src_rowfunc *src, ghost_mpi_comm mpicomm, double weight)
 Initializes a sparsemat from a row-based callback function. More...
 
ghost_error ghost_sparsemat_init_bin (ghost_sparsemat *mat, char *path, ghost_mpi_comm mpicomm, double weight)
 Initializes a sparsemat from a binary CRS file. More...
 
ghost_error ghost_sparsemat_init_mm (ghost_sparsemat *mat, char *path, ghost_mpi_comm mpicomm, double weight)
 Initializes a sparsemat from a Matrix Market file. More...
 
ghost_error ghost_sparsemat_init_crs (ghost_sparsemat *mat, ghost_gidx offs, ghost_lidx n, ghost_gidx *col, void *val, ghost_lidx *rpt, ghost_mpi_comm mpicomm, double weight)
 Initializes a sparsemat from local CRS data. More...
 
ghost_error ghost_sparsemat_to_bin (ghost_sparsemat *mat, char *path)
 Write a matrix to a binary CRS file. More...
 
size_t ghost_sparsemat_bytesize (ghost_sparsemat *mat)
 Get the entire memory footprint of the matrix. More...
 
ghost_error ghost_cu_sell1_spmv_selector (ghost_densemat *lhs_in, ghost_sparsemat *mat, ghost_densemat *rhs_in, ghost_spmv_opts traits)
 
ghost_error ghost_sell_kacz_selector (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, ghost_kacz_opts opts)
 Select and call the right SELL KACZ kernel. More...
 
ghost_error ghost_sell_kacz_shift_selector (ghost_densemat *x_real, ghost_densemat *x_imag, ghost_sparsemat *mat, ghost_densemat *b, double sigma_r, double sigma_i, ghost_kacz_opts opts)
 Select and call the right SELL KACZ kernel with complex shifts. More...
 
ghost_error ghost_kacz (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, ghost_kacz_opts opts)
 Perform a Kaczmarz sweep with the SELL matrix. More...
 
ghost_error ghost_kacz_mc (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, ghost_kacz_opts opts)
 
ghost_error ghost_kacz_rb (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, ghost_kacz_opts opts)
 
ghost_error ghost_kacz_bmc (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, ghost_kacz_opts opts)
 
ghost_error ghost_kacz_rb_with_shift (ghost_densemat *x, ghost_sparsemat *mat, ghost_densemat *b, double *shift_r, ghost_kacz_opts opts)
 
ghost_error ghost_carp (ghost_sparsemat *mat, ghost_densemat *x, ghost_densemat *b, ghost_carp_opts opts)
 
ghost_error checker (ghost_sparsemat *mat)
 
ghost_error split_transition (ghost_sparsemat *mat)
 
ghost_error ghost_carp_init (ghost_sparsemat *mat, ghost_densemat *b, ghost_carp_opts *opts)
 Initialize CARP. More...
 
ghost_error ghost_carp_perf_init (ghost_sparsemat *mat, ghost_carp_opts *opts)
 Finds optimum parameters for CARP. More...
 
void ghost_carp_destroy (ghost_carp_opts *opts)
 Destroy(Finalize) CARP. More...
 
ghost_error kacz_analyze_print (ghost_sparsemat *mat)
 Prints the row distribution details of KACZ. More...
 
ghost_error ghost_sparsemat_to_mm (char *path, ghost_sparsemat *mat)
 
ghost_error ghost_context_comm_init (ghost_context *ctx, ghost_gidx *col_orig, ghost_sparsemat *mat, ghost_lidx *col, ghost_lidx *nhalo)
 Assemble communication information in the given context. More...
 
ghost_error ghost_sparsemat_perm_global_cols (ghost_gidx *cols, ghost_lidx ncols, ghost_context *context)
 
static int ghost_sparsemat_rowfunc_crs (ghost_gidx row, ghost_lidx *rowlen, ghost_gidx *col, void *val, void *crsdata)
 
ghost_error calculate_bw (ghost_sparsemat *mat, void *matrixSource, ghost_sparsemat_src srcType)
 
ghost_error set_kacz_ratio (ghost_sparsemat *mat, void *matrixSource, ghost_sparsemat_src srcType)
 

Variables

const ghost_spmv_opts GHOST_SPMV_OPTS_INITIALIZER
 
const ghost_kacz_opts GHOST_KACZ_OPTS_INITIALIZER
 
const ghost_carp_opts GHOST_CARP_OPTS_INITIALIZER
 
const ghost_sparsemat_traits GHOST_SPARSEMAT_TRAITS_INITIALIZER
 Initialize sparse matrix traits with default values. More...
 

Detailed Description

Types and functions related to sparse matrices.

Author
Moritz Kreutzer morit.nosp@m.z.kr.nosp@m.eutze.nosp@m.r@fa.nosp@m.u.de

Macro Definition Documentation

#define GHOST_SCOTCH_STRAT_DEFAULT   "g"
#define GHOST_SELL_CHUNKHEIGHT_AUTO   -1

A chunkheight should automatically be determined.

#define GHOST_SELL_CHUNKHEIGHT_ELLPACK   0

Create only a single chunk, i.e., use the ELLPACK storage format.

#define GHOST_SPARSEMAT_PERM_ANY   (GHOST_SPARSEMAT_PERM_ANY_LOCAL|GHOST_SPARSEMAT_PERM_ANY_GLOBAL)

Combination of flags which apply any permutation to a ghost_sparsemat.

#define GHOST_SPARSEMAT_PERM_ANY_GLOBAL   (GHOST_SPARSEMAT_SCOTCHIFY|GHOST_SPARSEMAT_ZOLTAN)
#define GHOST_SPARSEMAT_SORT_GLOBAL   -1
#define GHOST_SPARSEMAT_SORT_LOCAL   -2
#define SPM_GNCOLS (   mat)    mat->context->col_map->gdim
#define SPM_NCHUNKS (   mat)    (mat->nchunks)
#define SPM_NCOLS (   mat)    mat->context->col_map->dim
#define SPM_NNZ (   mat)    mat->context->nnz
#define SPM_NROWS (   mat)    mat->context->row_map->dim
#define SPM_NROWSPAD (   mat)    mat->context->row_map->dimpad

Typedef Documentation

typedef ghost_error(* ghost_kacz_shift_kernel)(ghost_densemat *, ghost_densemat *, ghost_sparsemat *, ghost_densemat *, double, double, ghost_kacz_opts)

Enumeration Type Documentation

Enumerator
GHOST_KACZ_DIRECTION_UNDEFINED 
GHOST_KACZ_DIRECTION_FORWARD 
GHOST_KACZ_DIRECTION_BACKWARD 
Enumerator
GHOST_KACZ_MODE_NORMAL 
GHOST_KACZ_MODE_ECO 
GHOST_KACZ_MODE_PERFORMANCE 
Enumerator
GHOST_KACZ_NORMALIZE_YES 
GHOST_KACZ_NORMALIZE_NO 
GHOST_KACZ_NORMALIZE_SAVE 

Flags to a sparse matrix.

Enumerator
GHOST_SPARSEMAT_DEFAULT 

A default sparse matrix.

GHOST_SPARSEMAT_HOST 

Matrix is stored on host.

GHOST_SPARSEMAT_DEVICE 

Matrix is store on device.

GHOST_SPARSEMAT_NOT_PERMUTE_COLS 

If the matrix rows have been re-ordered, do NOT permute the column indices accordingly.

GHOST_SPARSEMAT_PERMUTE 

The matrix rows should be re-ordered in a certain way (defined in the traits).

GHOST_SPARSEMAT_NOT_SORT_COLS 

Do not sort the matrix cols wrt. memory location.

GHOST_SPARSEMAT_NOT_STORE_SPLIT 

Do NOT store the local and remote part of the matrix.

GHOST_SPARSEMAT_NOT_STORE_FULL 

Do NOT store the full matrix (local and remote combined).

GHOST_SPARSEMAT_SCOTCHIFY 

Reduce the matrix bandwidth with PT-Scotch.

GHOST_SPARSEMAT_SAVE_ORIG_COLS 

Save the un-compressed original columns of a distributed matrix.

GHOST_SPARSEMAT_COLOR 

Create a matrix permutation reflecting a distance-2-coloring.

GHOST_SPARSEMAT_TRANSPOSE_MM 

If the matrix comes from a matrix market file, transpose it on read-in. If this is implemented for other rowfuncs, the _MM may get removed in the future.

GHOST_SPARSEMAT_ZOLTAN 

Re-order the matrix globally using Zoltan hypergraph partitioning.

GHOST_SPARSEMAT_RCM 

Re-order the local part of the matrix using parallel RCM re-ordering.

GHOST_SPARSEMAT_BLOCKCOLOR 

Re-order the local part of the matrix using a block coloring.

GHOST_SOLVER_KACZ 

SETS the sparsematrix permutation as needed by the KACZ solver depending on the bandwidth of the matrix.

GHOST_SPARSEMAT_SORT_ROWS 

Sort matrix rows according to their length (SELL-C-Sigma sorting)

GHOST_SPARSEMAT_PERM_NO_DISTINCTION 

Does not make a distinction between local and remote entries if set; this might lead to higher communication time.

GHOST_SPARSEMAT_DIAG_FIRST 

Store the diagonal entry first in each row. Store an explicit zero if the diagonal is zero.

Possible sparse matrix symmetries.

Enumerator
GHOST_SPARSEMAT_SYMM_GENERAL 

Non-symmetric (general) matrix.

GHOST_SPARSEMAT_SYMM_SYMMETRIC 

Symmetric matrix.

GHOST_SPARSEMAT_SYMM_SKEW_SYMMETRIC 

Skew-symmetric matrix.

GHOST_SPARSEMAT_SYMM_HERMITIAN 

Hermitian matrix.

Function Documentation

ghost_error calculate_bw ( ghost_sparsemat mat,
void *  matrixSource,
ghost_sparsemat_src  srcType 
)
ghost_error checker ( ghost_sparsemat mat)

Here is the call graph for this function:

ghost_error ghost_carp ( ghost_sparsemat mat,
ghost_densemat x,
ghost_densemat b,
ghost_carp_opts  opts 
)

Here is the call graph for this function:

void ghost_carp_destroy ( ghost_carp_opts opts)

Destroy(Finalize) CARP.

Parameters
optsCARP options
ghost_error ghost_carp_init ( ghost_sparsemat mat,
ghost_densemat b,
ghost_carp_opts opts 
)

Initialize CARP.

Parameters
[in]matThe sparsematrix
[in]bThe rhs
optsCARP Options.

Here is the call graph for this function:

ghost_error ghost_carp_perf_init ( ghost_sparsemat mat,
ghost_carp_opts opts 
)

Finds optimum parameters for CARP.

Parameters
[in]matThe sparsematrix
optsOptions.

Here is the call graph for this function:

ghost_error ghost_context_comm_init ( ghost_context ctx,
ghost_gidx col_orig,
ghost_sparsemat mat,
ghost_lidx col,
ghost_lidx nhalo 
)

Assemble communication information in the given context.

Parameters
[in,out]ctxThe context.
[in]col_origThe original column indices of the sparse matrix which is bound to the context.
[out]colThe compressed column indices of the sparse matrix which is bound to the context.
[out]nhaloThe number of halo elements.
Returns
GHOST_SUCCESS on success or an error indicator.

The following fields of ghost_context are being filled in this function: wishes, wishlist, dues, duelist, hput_pos, wishpartners, nwishpartners, duepartners, nduepartners. Additionally, the columns in col_orig are being compressed and stored in col.

Here is the call graph for this function:

ghost_error ghost_cu_sell1_spmv_selector ( ghost_densemat lhs_in,
ghost_sparsemat mat,
ghost_densemat rhs_in,
ghost_spmv_opts  traits 
)
ghost_error ghost_cu_sell_spmv_selector ( ghost_densemat lhs,
ghost_sparsemat mat,
ghost_densemat rhs,
ghost_spmv_opts  traits 
)

Select and call the right CUDA SELL SpMV kernel.

Parameters
matThe matrix.
lhsThe result densemat.
rhsThe input densemat.
traitsThe SpMV traits.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_kacz ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
ghost_kacz_opts  opts 
)

Perform a Kaczmarz sweep with the SELL matrix.

Parameters
xOutput densemat.
matThe matrix.
bInput densemat.
optsOptions.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_kacz_bmc ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
ghost_kacz_opts  opts 
)

Here is the call graph for this function:

ghost_error ghost_kacz_mc ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
ghost_kacz_opts  opts 
)

Here is the call graph for this function:

ghost_error ghost_kacz_rb ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
ghost_kacz_opts  opts 
)

Here is the call graph for this function:

ghost_error ghost_kacz_rb_with_shift ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
double *  shift_r,
ghost_kacz_opts  opts 
)

Here is the call graph for this function:

ghost_error ghost_sell_kacz_selector ( ghost_densemat x,
ghost_sparsemat mat,
ghost_densemat b,
ghost_kacz_opts  opts 
)

Select and call the right SELL KACZ kernel.

Parameters
matThe matrix.
xThe result densemat.
bThe input densemat.
optsThe kacz Options.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sell_kacz_shift_selector ( ghost_densemat x_real,
ghost_densemat x_imag,
ghost_sparsemat mat,
ghost_densemat b,
double  sigma_r,
double  sigma_i,
ghost_kacz_opts  opts 
)

Select and call the right SELL KACZ kernel with complex shifts.

Parameters
matThe matrix.
x_real(densemat) The real part of result.
x_imag(densemat) The imaginary part of result.
bThe input densemat.
sigma_rThe real part of shift
sigma_iThe imaginary part of shift
optsThe kacz Options.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sell_spmv_selector ( ghost_densemat lhs,
ghost_sparsemat mat,
ghost_densemat rhs,
ghost_spmv_opts  traits 
)

Select and call the right SELL SpMV kernel.

Parameters
matThe matrix.
lhsThe result densemat.
rhsThe input densemat.
traitsThe SpMV traits.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_sell_stringify_selector ( ghost_sparsemat mat,
char **  str,
int  dense 
)

Select and call the right SELL stringification function.

Parameters
matThe matrix.
strWhere to store the string.
densePrint in a dense or sparse manner.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sparsemat_blockColor ( ghost_context ctx,
ghost_sparsemat mat 
)

Here is the call graph for this function:

size_t ghost_sparsemat_bytesize ( ghost_sparsemat mat)

Get the entire memory footprint of the matrix.

Parameters
matThe matrix.
Returns
The memory footprint of the matrix in bytes or zero if the matrix is not valid.
ghost_error ghost_sparsemat_fromfile_common ( ghost_sparsemat mat,
char *  matrixPath,
ghost_lidx **  rpt 
)

Common function for matrix creation from a file.

This function does work which is independent of the storage format and it should be called at the beginning of a sparse matrix' fromFile function. This function also reads the row pointer from the file and stores it into the parameter rpt.

Parameters
[in]matThe sparse matrix.
[in]matrixPathThe matrix file path.
[out]rptWhere to store the row pointer information which may be needed afterwards.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sparsemat_nnz ( ghost_gidx nnz,
ghost_sparsemat mat 
)

Obtain the global number of nonzero elements of a sparse matrix.

Parameters
[out]nnzWhere to store the result.
[in]matThe sparse matrix.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sparsemat_nrows ( ghost_gidx nrows,
ghost_sparsemat mat 
)

Obtain the global number of rows of a sparse matrix.

Parameters
[out]nrowsWhere to store the result.
[in]matThe sparse matrix.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sparsemat_perm_color ( ghost_context ctx,
ghost_sparsemat mat 
)

Create a matrix permutation based on 2-way coloring using ColPack.

Parameters
[out]ctxThe context in which to store the permutations and color information.
[in]matThe unpermuted SELL-1-1 source sparse matrix.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_sparsemat_perm_global_cols ( ghost_gidx cols,
ghost_lidx  ncols,
ghost_context context 
)

Here is the call graph for this function:

ghost_error ghost_sparsemat_perm_scotch ( ghost_context ctx,
ghost_sparsemat mat 
)

Create a matrix permutation based on (PT-)SCOTCH.

Parameters
[out]ctxThe context in which to store the permutations.
[in]matThe unpermuted SELL-1-1 source sparse matrix.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_sparsemat_perm_sort ( ghost_context ctx,
ghost_sparsemat mat,
ghost_lidx  scope 
)

Create a matrix permutation based on row length sorting within a given scope.

Parameters
[out]ctxThe context in which to store the permutations.
[in]matThe unpermuted SELL-1-1 source sparse matrix.
[in]scopeThe sorting scope.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

ghost_error ghost_sparsemat_perm_spmp ( ghost_context ctx,
ghost_sparsemat mat 
)

Here is the call graph for this function:

ghost_error ghost_sparsemat_perm_zoltan ( ghost_context ctx,
ghost_sparsemat mat 
)

Here is the call graph for this function:

ghost_error ghost_sparsemat_registerrow ( ghost_sparsemat mat,
ghost_gidx  row,
ghost_gidx col,
ghost_lidx  rowlen,
ghost_lidx  stride 
)

Store matrix information like bandwidth and nonzero distribution for a given matrix row.

This function should be called at matrix creation for each row.

Parameters
[out]matThe matrix.
[in]rowThe row index.
[in]colThe column indices of the row.
[in]rowlenThe length of the row.
[in]strideThe stride for the parameter col.
Returns
GHOST_SUCCESS on success or an error indicator.
ghost_error ghost_sparsemat_registerrow_finalize ( ghost_sparsemat mat)

Finalize the storing of matrix information like bandwidth and nonzero distribution.

This function should be after matrix creation. It finalizes the processing of information obtained in ghost_sparsemat_registerrow.

Parameters
[out]matThe matrix.
Returns
GHOST_SUCCESS on success or an error indicator.

Here is the call graph for this function:

static int ghost_sparsemat_rowfunc_crs ( ghost_gidx  row,
ghost_lidx rowlen,
ghost_gidx col,
void *  val,
void *  crsdata 
)
inlinestatic
bool ghost_sparsemat_symmetry_valid ( ghost_sparsemat_symmetry  symmetry)

Check if the symmetry information of a sparse matrix is valid.

Parameters
[in]symmetryThe symmetry information.
Returns
True if it is valid, false otherwise.
ghost_error ghost_sparsemat_to_bin ( ghost_sparsemat mat,
char *  path 
)

Write a matrix to a binary CRS file.

Parameters
matThe matrix.
pathPath of the file.
ghost_error ghost_sparsemat_to_mm ( char *  path,
ghost_sparsemat mat 
)

Here is the call graph for this function:

ghost_error kacz_analyze_print ( ghost_sparsemat mat)

Prints the row distribution details of KACZ.

Parameters
[in]matThe sparsematrix

Here is the call graph for this function:

ghost_error set_kacz_ratio ( ghost_sparsemat mat,
void *  matrixSource,
ghost_sparsemat_src  srcType 
)
ghost_error split_transition ( ghost_sparsemat mat)

Here is the call graph for this function:

Variable Documentation

const ghost_carp_opts GHOST_CARP_OPTS_INITIALIZER
const ghost_kacz_opts GHOST_KACZ_OPTS_INITIALIZER
const ghost_sparsemat_traits GHOST_SPARSEMAT_TRAITS_INITIALIZER

Initialize sparse matrix traits with default values.

const ghost_spmv_opts GHOST_SPMV_OPTS_INITIALIZER