single precision
[Sparse auxiliary]

Functions

magma_int_t magma_s_precond (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is chosen.
magma_int_t magma_s_precondsetup (magma_s_matrix A, magma_s_matrix b, magma_s_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is preprocessed.
magma_int_t magma_s_applyprecond (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is applied.
magma_int_t magma_s_applyprecond_left (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective left preconditioner is applied.
magma_int_t magma_s_applyprecond_right (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective right-preconditioner is applied.
magma_int_t magma_s_solver (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_sopts *zopts, magma_queue_t queue)
 ALlows the user to choose a solver.
magma_int_t magma_sapplycustomprecond_l (magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 This is an interface to the left solve for any custom preconditioner.
magma_int_t magma_sapplycustomprecond_r (magma_s_matrix b, magma_s_matrix *x, magma_s_preconditioner *precond, magma_queue_t queue)
 This is an interface to the right solve for any custom preconditioner.
magma_int_t magma_sresidual (magma_s_matrix A, magma_s_matrix b, magma_s_matrix x, float *res, magma_queue_t queue)
 Computes the residual ||b-Ax|| for a solution approximation x.
magma_int_t magma_sresidualvec (magma_s_matrix A, magma_s_matrix b, magma_s_matrix x, magma_s_matrix *r, float *res, magma_queue_t queue)
 Computes the residual r = b-Ax for a solution approximation x.
magma_int_t magma_scsrsplit (magma_int_t bsize, magma_s_matrix A, magma_s_matrix *D, magma_s_matrix *R, magma_queue_t queue)
 Splits a CSR matrix into two matrices, one containing the diagonal blocks with the diagonal element stored first, one containing the rest of the original matrix.
magma_int_t magma_sdomainoverlap (magma_index_t num_rows, magma_index_t *num_indices, magma_index_t *rowptr, magma_index_t *colidx, magma_index_t *x, magma_queue_t queue)
 Generates the update list.
magma_int_t magma_smfree (magma_s_matrix *A, magma_queue_t queue)
 Free the memory of a magma_s_matrix.
magma_int_t magma_sfrobenius (magma_s_matrix A, magma_s_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B.
magma_int_t magma_snonlinres (magma_s_matrix A, magma_s_matrix L, magma_s_matrix U, magma_s_matrix *LU, real_Double_t *res, magma_queue_t queue)
 Computes the nonlinear residual A - LU and returns the difference as well es the Frobenius norm of the difference.
magma_int_t magma_silures (magma_s_matrix A, magma_s_matrix L, magma_s_matrix U, magma_s_matrix *LU, real_Double_t *res, real_Double_t *nonlinres, magma_queue_t queue)
 Computes the ILU residual A - LU and returns the difference as well es the Frobenius norm of the difference.
magma_int_t magma_sicres (magma_s_matrix A, magma_s_matrix C, magma_s_matrix CT, magma_s_matrix *LU, real_Double_t *res, real_Double_t *nonlinres, magma_queue_t queue)
 Computes the IC residual A - CC^T and returns the difference as well es the Frobenius norm of the difference.
magma_int_t magma_sinitguess (magma_s_matrix A, magma_s_matrix *L, magma_s_matrix *U, magma_queue_t queue)
 Computes an initial guess for the iterative ILU/IC.
magma_int_t magma_sinitrecursiveLU (magma_s_matrix A, magma_s_matrix *B, magma_queue_t queue)
 Using the iterative approach of computing ILU factorizations with increasing fill-in, it takes the input matrix A, containing the approximate factors, ( L and U as well ) computes a matrix with one higher level of fill-in, inserts the original approximation as initial guess, and provides the factors L and U also filled with the scaled initial guess.
magma_int_t magma_smLdiagadd (magma_s_matrix *L, magma_queue_t queue)
 Checks for a lower triangular matrix whether it is strictly lower triangular and in the negative case adds a unit diagonal.
magma_int_t magma_srowentries (magma_s_matrix *A, magma_queue_t queue)
 Checks the maximal number of nonzeros in a row of matrix A.
magma_int_t magma_sdiameter (magma_s_matrix *A, magma_queue_t queue)
 Computes the diameter of a sparse matrix and stores the value in diameter.
magma_int_t magma_s_csr_compressor (float **val, magma_index_t **row, magma_index_t **col, float **valn, magma_index_t **rown, magma_index_t **coln, magma_int_t *n, magma_queue_t queue)
 Helper function to compress CSR containing zero-entries.
magma_int_t magma_smconvert (magma_s_matrix A, magma_s_matrix *B, magma_storage_t old_format, magma_storage_t new_format, magma_queue_t queue)
 Converter between different sparse storage formats.
magma_int_t magma_smcsrcompressor (magma_s_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix.
magma_int_t magma_scsrset (magma_int_t m, magma_int_t n, magma_index_t *row, magma_index_t *col, float *val, magma_s_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA.
magma_int_t magma_scsrget (magma_s_matrix A, magma_int_t *m, magma_int_t *n, magma_index_t **row, magma_index_t **col, float **val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure.
magma_int_t magma_scsrset_gpu (magma_int_t m, magma_int_t n, magmaIndex_ptr row, magmaIndex_ptr col, magmaFloat_ptr val, magma_s_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA (located on DEV).
magma_int_t magma_scsrget_gpu (magma_s_matrix A, magma_int_t *m, magma_int_t *n, magmaIndex_ptr *row, magmaIndex_ptr *col, magmaFloat_ptr *val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure (located on DEV).
magma_int_t magma_smdiff (magma_s_matrix A, magma_s_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B.
magma_int_t magma_smgenerator (magma_int_t n, magma_int_t offdiags, magma_index_t *diag_offset, float *diag_vals, magma_s_matrix *A, magma_queue_t queue)
 Generate a symmetric n x n CSR matrix for a stencil.
magma_int_t magma_sm_27stencil (magma_int_t n, magma_s_matrix *A, magma_queue_t queue)
 Generate a 27-point stencil for a 3D FD discretization.
magma_int_t magma_sm_5stencil (magma_int_t n, magma_s_matrix *A, magma_queue_t queue)
 Generate a 5-point stencil for a 2D FD discretization.
magma_int_t magma_ssymbilu (magma_s_matrix *A, magma_int_t levels, magma_s_matrix *L, magma_s_matrix *U, magma_queue_t queue)
 This routine performs a symbolic ILU factorization.
magma_int_t read_s_csr_from_mtx (magma_storage_t *type, magma_location_t *location, magma_int_t *n_row, magma_int_t *n_col, magma_int_t *nnz, float **val, magma_index_t **row, magma_index_t **col, const char *filename, magma_queue_t queue)
 Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
magma_int_t magma_swrite_csr_mtx (magma_s_matrix A, magma_order_t MajorType, const char *filename, magma_queue_t queue)
 Writes a CSR matrix to a file using Matrix Market format.
magma_int_t magma_sprint_csr_mtx (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, float **val, magma_index_t **row, magma_index_t **col, magma_order_t MajorType, magma_queue_t queue)
 Prints a CSR matrix in Matrix Market format.
magma_int_t magma_sprint_csr (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, float **val, magma_index_t **row, magma_index_t **col, magma_queue_t queue)
 Prints a CSR matrix in CSR format.
magma_int_t magma_sprint_matrix (magma_s_matrix A, magma_queue_t queue)
 Prints a sparse matrix in CSR format.
magma_int_t magma_s_csr_mtx (magma_s_matrix *A, const char *filename, magma_queue_t queue)
 Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
magma_int_t magma_s_csr_mtxsymm (magma_s_matrix *A, const char *filename, magma_queue_t queue)
 Reads in a SYMMETRIC matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
magma_int_t magma_smlumerge (magma_s_matrix L, magma_s_matrix U, magma_s_matrix *A, magma_queue_t queue)
 Takes an strictly lower triangular matrix L and an upper triangular matrix U and merges them into a matrix A containing the upper and lower triangular parts.
magma_int_t magma_smscale (magma_s_matrix *A, magma_scale_t scaling, magma_queue_t queue)
 Scales a matrix.
magma_int_t magma_smdiagadd (magma_s_matrix *A, float add, magma_queue_t queue)
 Adds a multiple of the Identity matrix to a matrix: A = A+add * I.
magma_int_t magma_smtransfer (magma_s_matrix A, magma_s_matrix *B, magma_location_t src, magma_location_t dst, magma_queue_t queue)
 Copies a matrix from memory location src to memory location dst.
magma_int_t s_transpose_csr (magma_int_t n_rows, magma_int_t n_cols, magma_int_t nnz, float *values, magma_index_t *rowptr, magma_index_t *colind, magma_int_t *new_n_rows, magma_int_t *new_n_cols, magma_int_t *new_nnz, float **new_values, magma_index_t **new_rowptr, magma_index_t **new_colind, magma_queue_t queue)
 Transposes a matrix stored in CSR format on the CPU host.
magma_int_t magma_smtranspose (magma_s_matrix A, magma_s_matrix *B, magma_queue_t queue)
 Interface to cuSPARSE transpose.
magma_int_t magma_s_cucsrtranspose (magma_s_matrix A, magma_s_matrix *B, magma_queue_t queue)
 Helper function to transpose CSR matrix.
magma_int_t magma_ssolverinfo (magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Prints information about a previously called solver.
magma_int_t magma_ssolverinfo_free (magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Frees any memory assocoiated with the verbose mode of solver_par.
magma_int_t magma_ssolverinfo_init (magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Initializes all solver and preconditioner parameters.
magma_int_t magma_seigensolverinfo_init (magma_s_solver_par *solver_par, magma_queue_t queue)
 Initializes space for eigensolvers.
magma_int_t magma_sindexsort (magma_index_t *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers.
magma_int_t magma_sparse_opts (int argc, char **argv, magma_sopts *opts, int *matrices, magma_queue_t queue)
 Parses input options for a solver.
magma_int_t magma_svinit (magma_s_matrix *x, magma_location_t mem_loc, magma_int_t num_rows, magma_int_t num_cols, float values, magma_queue_t queue)
 Allocates memory for magma_s_matrix and initializes it with the passed value.
magma_int_t magma_sprint_vector (magma_s_matrix x, magma_int_t offset, magma_int_t visulen, magma_queue_t queue)
 Visualizes part of a vector of type magma_s_matrix.
magma_int_t magma_svread (magma_s_matrix *x, magma_int_t length, char *filename, magma_queue_t queue)
 Reads in a float vector of length "length".
magma_int_t magma_svspread (magma_s_matrix *x, const char *filename, magma_queue_t queue)
 Reads in a sparse vector-block stored in COO format.
magma_int_t magma_svset (magma_int_t m, magma_int_t n, float *val, magma_s_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA.
magma_int_t magma_svget (magma_s_matrix v, magma_int_t *m, magma_int_t *n, float **val, magma_queue_t queue)
 Passes a MAGMA vector back.
magma_int_t magma_svset_dev (magma_int_t m, magma_int_t n, magmaFloat_ptr val, magma_s_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA (located on DEV).
magma_int_t magma_svget_dev (magma_s_matrix v, magma_int_t *m, magma_int_t *n, magmaFloat_ptr *val, magma_queue_t queue)
 Passes a MAGMA vector back (located on DEV).
magma_int_t magma_svtranspose (magma_s_matrix x, magma_s_matrix *y, magma_queue_t queue)
 Transposes a vector from col to row major and vice versa.
magma_int_t magma_s_spmv_shift (float alpha, magma_s_matrix A, float lambda, magma_s_matrix x, float beta, magma_int_t offset, magma_int_t blocksize, magma_index_t *add_rows, magma_s_matrix y, magma_queue_t queue)
 For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * ( A - lambda I ) * x + beta * y.
magma_int_t magma_smcsrcompressor_gpu (magma_s_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix.
magma_int_t magma_slobpcg_res (magma_int_t num_rows, magma_int_t num_vecs, magmaFloat_ptr evalues, magmaFloat_ptr X, magmaFloat_ptr R, magmaFloat_ptr res, magma_queue_t queue)
 This routine computes for Block-LOBPCG, the set of residuals.
magma_int_t magma_slobpcg_shift (magma_int_t num_rows, magma_int_t num_vecs, magma_int_t shift, magmaFloat_ptr x, magma_queue_t queue)
 For a Block-LOBPCG, the set of residuals (entries consecutive in memory) shrinks and the vectors are shifted in case shift residuals drop below threshold.
magma_int_t magma_scopyscale (int n, int k, magmaFloat_ptr r, magmaFloat_ptr v, magmaFloat_ptr skp, magma_queue_t queue)
 Computes the correction term of the pipelined GMRES according to P.

Function Documentation

magma_int_t magma_s_applyprecond ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is applied.

E.g. for Jacobi: the scaling-vetor, for ILU the triangular solves.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector b
[in,out] x magma_s_matrix* output vector x
[in] precond magma_s_preconditioner preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_applyprecond_left ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and the preconditioner parameters, the respective left preconditioner is applied.

E.g. for Jacobi: the scaling-vetor, for ILU the left triangular solve.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector b
[in,out] x magma_s_matrix* output vector x
[in] precond magma_s_preconditioner preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_applyprecond_right ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and the preconditioner parameters, the respective right-preconditioner is applied.

E.g. for Jacobi: the scaling-vetor, for ILU the right triangular solve.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector b
[in,out] x magma_s_matrix* output vector x
[in] precond magma_s_preconditioner preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_csr_compressor ( float **  val,
magma_index_t **  row,
magma_index_t **  col,
float **  valn,
magma_index_t **  rown,
magma_index_t **  coln,
magma_int_t *  n,
magma_queue_t  queue 
)

Helper function to compress CSR containing zero-entries.

Parameters:
[in] val float** input val pointer to compress
[in] row magma_int_t** input row pointer to modify
[in] col magma_int_t** input col pointer to compress
[in] valn float** output val pointer
[out] rown magma_int_t** output row pointer
[out] coln magma_int_t** output col pointer
[out] n magma_int_t* number of rows in matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_csr_mtx ( magma_s_matrix *  A,
const char *  filename,
magma_queue_t  queue 
)

Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.

It duplicates the off-diagonal entries in the symmetric case.

Parameters:
[out] A magma_s_matrix* matrix in magma sparse matrix format
[in] filename const char* filname of the mtx matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_csr_mtxsymm ( magma_s_matrix *  A,
const char *  filename,
magma_queue_t  queue 
)

Reads in a SYMMETRIC matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.

It does not duplicate the off-diagonal entries!

Parameters:
[out] A magma_s_matrix* matrix in magma sparse matrix format
[in] filename const char* filname of the mtx matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_cucsrtranspose ( magma_s_matrix  A,
magma_s_matrix *  B,
magma_queue_t  queue 
)

Helper function to transpose CSR matrix.

Using the CUSPARSE CSR2CSC function.

Parameters:
[in] A magma_s_matrix input matrix (CSR)
[out] B magma_s_matrix* output matrix (CSR)
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_precond ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is chosen.

It approximates x for A x = y.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector b
[in] x magma_s_matrix* output vector x
[in,out] precond magma_s_preconditioner preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_precondsetup ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is preprocessed.

E.g. for Jacobi: the scaling-vetor, for ILU the factorization.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector y
[in,out] precond magma_s_preconditioner preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_solver ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_sopts *  zopts,
magma_queue_t  queue 
)

ALlows the user to choose a solver.

Parameters:
[in] A magma_s_matrix sparse matrix A
[in] b magma_s_matrix input vector b
[in] x magma_s_matrix* output vector x
[in] zopts magma_sopts options for solver and preconditioner
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_s_spmv_shift ( float  alpha,
magma_s_matrix  A,
float  lambda,
magma_s_matrix  x,
float  beta,
magma_int_t  offset,
magma_int_t  blocksize,
magma_index_t *  add_rows,
magma_s_matrix  y,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * ( A - lambda I ) * x + beta * y.

Parameters:
alpha float scalar alpha
A magma_s_matrix sparse matrix A
lambda float scalar lambda
x magma_s_matrix input vector x
beta float scalar beta
offset magma_int_t in case not the main diagonal is scaled
blocksize magma_int_t in case of processing multiple vectors
add_rows magma_int_t* in case the matrixpowerskernel is used
y magma_s_matrix output vector y
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sapplycustomprecond_l ( magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

This is an interface to the left solve for any custom preconditioner.

It should compute x = FUNCTION(b) The vectors are located on the device.

Parameters:
[in] b magma_s_matrix RHS
[in,out] x magma_s_matrix* vector to precondition
[in,out] precond magma_s_preconditioner* preconditioner parameters
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sapplycustomprecond_r ( magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_preconditioner *  precond,
magma_queue_t  queue 
)

This is an interface to the right solve for any custom preconditioner.

It should compute x = FUNCTION(b) The vectors are located on the device.

Parameters:
[in] b magma_s_matrix RHS
[in,out] x magma_s_matrix* vector to precondition
[in,out] precond magma_s_preconditioner* preconditioner parameters
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scopyscale ( int  n,
int  k,
magmaFloat_ptr  r,
magmaFloat_ptr  v,
magmaFloat_ptr  skp,
magma_queue_t  queue 
)

Computes the correction term of the pipelined GMRES according to P.

Ghysels and scales and copies the new search direction

Returns the vector v = r/ ( skp[k] - (sum_i=1^k skp[i]^2) ) .

Parameters:
[in] n int length of v_i
[in] k int # skp entries v_i^T * r ( without r )
[in] r magmaFloat_ptr vector of length n
[in] v magmaFloat_ptr vector of length n
[in] skp magmaFloat_ptr array of parameters
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scsrget ( magma_s_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magma_index_t **  row,
magma_index_t **  col,
float **  val,
magma_queue_t  queue 
)

Passes a MAGMA matrix to CSR structure.

Parameters:
[in] A magma_s_matrix magma sparse matrix in CSR format
[out] m magma_int_t number of rows
[out] n magma_int_t number of columns
[out] row magma_index_t* row pointer
[out] col magma_index_t* column indices
[out] val float* array containing matrix entries
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scsrget_gpu ( magma_s_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magmaIndex_ptr *  row,
magmaIndex_ptr *  col,
magmaFloat_ptr *  val,
magma_queue_t  queue 
)

Passes a MAGMA matrix to CSR structure (located on DEV).

Parameters:
[in] A magma_s_matrix magma sparse matrix in CSR format
[out] m magma_int_t number of rows
[out] n magma_int_t number of columns
[out] row magmaIndex_ptr row pointer
[out] col magmaIndex_ptr column indices
[out] val magmaFloat_ptr array containing matrix entries
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scsrset ( magma_int_t  m,
magma_int_t  n,
magma_index_t *  row,
magma_index_t *  col,
float *  val,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Passes a CSR matrix to MAGMA.

Parameters:
[in] m magma_int_t number of rows
[in] n magma_int_t number of columns
[in] row magma_index_t* row pointer
[in] col magma_index_t* column indices
[in] val float* array containing matrix entries
[out] A magma_s_matrix* matrix in magma sparse matrix format
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scsrset_gpu ( magma_int_t  m,
magma_int_t  n,
magmaIndex_ptr  row,
magmaIndex_ptr  col,
magmaFloat_ptr  val,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Passes a CSR matrix to MAGMA (located on DEV).

Parameters:
[in] m magma_int_t number of rows
[in] n magma_int_t number of columns
[in] row magmaIndex_ptr row pointer
[in] col magmaIndex_ptr column indices
[in] val magmaFloat_ptr array containing matrix entries
[out] A magma_s_matrix* matrix in magma sparse matrix format
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_scsrsplit ( magma_int_t  bsize,
magma_s_matrix  A,
magma_s_matrix *  D,
magma_s_matrix *  R,
magma_queue_t  queue 
)

Splits a CSR matrix into two matrices, one containing the diagonal blocks with the diagonal element stored first, one containing the rest of the original matrix.

Parameters:
[in] bsize magma_int_t size of the diagonal blocks
[in] A magma_s_matrix CSR input matrix
[out] D magma_s_matrix* CSR matrix containing diagonal blocks
[out] R magma_s_matrix* CSR matrix containing rest
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sdiameter ( magma_s_matrix *  A,
magma_queue_t  queue 
)

Computes the diameter of a sparse matrix and stores the value in diameter.

Parameters:
[in,out] A magma_s_matrix* sparse matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sdomainoverlap ( magma_index_t  num_rows,
magma_index_t *  num_indices,
magma_index_t *  rowptr,
magma_index_t *  colidx,
magma_index_t *  x,
magma_queue_t  queue 
)

Generates the update list.

Parameters:
[in] x magma_index_t* array to sort
[in] num_rows magma_int_t number of rows in matrix
[out] num_indices magma_int_t* number of indices in array
[in] rowptr magma_index_t* rowpointer of matrix
[in] colidx magma_index_t* colindices of matrix
[in] x magma_index_t* array containing indices for domain overlap
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_seigensolverinfo_init ( magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Initializes space for eigensolvers.

Parameters:
[in,out] solver_par magma_s_solver_par* structure containing all solver information
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sfrobenius ( magma_s_matrix  A,
magma_s_matrix  B,
real_Double_t *  res,
magma_queue_t  queue 
)

Computes the Frobenius norm of the difference between the CSR matrices A and B.

They need to share the same sparsity pattern!

Parameters:
[in] A magma_s_matrix sparse matrix in CSR
[in] B magma_s_matrix sparse matrix in CSR
[out] res real_Double_t* Frobenius norm of difference
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sicres ( magma_s_matrix  A,
magma_s_matrix  C,
magma_s_matrix  CT,
magma_s_matrix *  LU,
real_Double_t *  res,
real_Double_t *  nonlinres,
magma_queue_t  queue 
)

Computes the IC residual A - CC^T and returns the difference as well es the Frobenius norm of the difference.

Parameters:
[in] A magma_s_matrix input sparse matrix in CSR
[in] C magma_s_matrix input sparse matrix in CSR
[in] CT magma_s_matrix input sparse matrix in CSR
[in] LU magma_s_matrix* output sparse matrix in A-LU in CSR
[out] res real_Double_t* IC residual
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_silures ( magma_s_matrix  A,
magma_s_matrix  L,
magma_s_matrix  U,
magma_s_matrix *  LU,
real_Double_t *  res,
real_Double_t *  nonlinres,
magma_queue_t  queue 
)

Computes the ILU residual A - LU and returns the difference as well es the Frobenius norm of the difference.

Parameters:
[in] A magma_s_matrix input sparse matrix in CSR
[in] L magma_s_matrix input sparse matrix in CSR
[in] U magma_s_matrix input sparse matrix in CSR
[out] LU magma_s_matrix* output sparse matrix in A-LU in CSR
[out] res real_Double_t* Frobenius norm of difference
[out] nonlinres real_Double_t* Frobenius norm of difference
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sindexsort ( magma_index_t *  x,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of integers.

Parameters:
[in,out] x magma_index_t* array to sort
[in] first magma_int_t pointer to first element
[in] last magma_int_t pointer to last element
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sinitguess ( magma_s_matrix  A,
magma_s_matrix *  L,
magma_s_matrix *  U,
magma_queue_t  queue 
)

Computes an initial guess for the iterative ILU/IC.

Parameters:
[in] A magma_s_matrix sparse matrix in CSR
[out] L magma_s_matrix* sparse matrix in CSR
[out] U magma_s_matrix* sparse matrix in CSR
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sinitrecursiveLU ( magma_s_matrix  A,
magma_s_matrix *  B,
magma_queue_t  queue 
)

Using the iterative approach of computing ILU factorizations with increasing fill-in, it takes the input matrix A, containing the approximate factors, ( L and U as well ) computes a matrix with one higher level of fill-in, inserts the original approximation as initial guess, and provides the factors L and U also filled with the scaled initial guess.

Parameters:
[in] A magma_s_matrix* sparse matrix in CSR
[out] L magma_s_matrix* sparse matrix in CSR
[out] U magma_s_matrix* sparse matrix in CSR
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_slobpcg_res ( magma_int_t  num_rows,
magma_int_t  num_vecs,
magmaFloat_ptr  evalues,
magmaFloat_ptr  X,
magmaFloat_ptr  R,
magmaFloat_ptr  res,
magma_queue_t  queue 
)

This routine computes for Block-LOBPCG, the set of residuals.

R = Ax - x evalues It replaces: for(int i=0; i < n; i++){ magma_saxpy(m, MAGMA_S_MAKE(-evalues[i],0),blockX+i*m,1,blockR+i*m,1); } The memory layout of x is:

/ x1[0] x2[0] x3[0] \ | x1[1] x2[1] x3[1] | x = | x1[2] x2[2] x3[2] | = x1[0] x1[1] x1[2] x1[3] x1[4] x2[0] x2[1] . | x1[3] x2[3] x3[3] | \ x1[4] x2[4] x3[4] /

Parameters:
[in] num_rows magma_int_t number of rows
[in] num_vecs magma_int_t number of vectors
[in] evalues magmaFloat_ptr array of eigenvalues/approximations
[in] X magmaFloat_ptr block of eigenvector approximations
[in] R magmaFloat_ptr block of residuals
[in] res magmaFloat_ptr array of residuals
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_slobpcg_shift ( magma_int_t  num_rows,
magma_int_t  num_vecs,
magma_int_t  shift,
magmaFloat_ptr  x,
magma_queue_t  queue 
)

For a Block-LOBPCG, the set of residuals (entries consecutive in memory) shrinks and the vectors are shifted in case shift residuals drop below threshold.

The memory layout of x is:

/ x1[0] x2[0] x3[0] \ | x1[1] x2[1] x3[1] | x = | x1[2] x2[2] x3[2] | = x1[0] x2[0] x3[0] x1[1] x2[1] x3[1] x1[2] . | x1[3] x2[3] x3[3] | \ x1[4] x2[4] x3[4] /

Parameters:
[in] num_rows magma_int_t number of rows
[in] num_vecs magma_int_t number of vectors
[in] shift magma_int_t shift number
in/out] x magmaFloat_ptr input/output vector x
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sm_27stencil ( magma_int_t  n,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Generate a 27-point stencil for a 3D FD discretization.

Parameters:
[in] n magma_int_t number of rows
[out] A magma_s_matrix* matrix to generate
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sm_5stencil ( magma_int_t  n,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Generate a 5-point stencil for a 2D FD discretization.

Parameters:
[in] n magma_int_t number of rows
[out] A magma_s_matrix* matrix to generate
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smconvert ( magma_s_matrix  A,
magma_s_matrix *  B,
magma_storage_t  old_format,
magma_storage_t  new_format,
magma_queue_t  queue 
)

Converter between different sparse storage formats.

Parameters:
[in] A magma_s_matrix sparse matrix A
[out] B magma_s_matrix* copy of A in new format
[in] old_format magma_storage_t original storage format
[in] new_format magma_storage_t new storage format
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smcsrcompressor ( magma_s_matrix *  A,
magma_queue_t  queue 
)

Removes zeros in a CSR matrix.

Parameters:
[in,out] A magma_s_matrix* input/output matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smcsrcompressor_gpu ( magma_s_matrix *  A,
magma_queue_t  queue 
)

Removes zeros in a CSR matrix.

This is a GPU implementation of the CSR compressor.

Parameters:
A magma_s_matrix* input/output matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smdiagadd ( magma_s_matrix *  A,
float  add,
magma_queue_t  queue 
)

Adds a multiple of the Identity matrix to a matrix: A = A+add * I.

Parameters:
[in,out] A magma_s_matrix* input/output matrix
[in] add float scaling for the identity matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smdiff ( magma_s_matrix  A,
magma_s_matrix  B,
real_Double_t *  res,
magma_queue_t  queue 
)

Computes the Frobenius norm of the difference between the CSR matrices A and B.

They do not need to share the same sparsity pattern!

res = ||A-B||_F = sqrt( sum_ij (A_ij-B_ij)^2 )

Parameters:
[in] A magma_s_matrix sparse matrix in CSR
[in] B magma_s_matrix sparse matrix in CSR
[out] res real_Double_t* residual
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smfree ( magma_s_matrix *  A,
magma_queue_t  queue 
)

Free the memory of a magma_s_matrix.

Parameters:
[in,out] A magma_s_matrix* matrix to free
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smgenerator ( magma_int_t  n,
magma_int_t  offdiags,
magma_index_t *  diag_offset,
float *  diag_vals,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Generate a symmetric n x n CSR matrix for a stencil.

Parameters:
[in] n magma_int_t number of rows
[in] offdiags magma_int_t number of offdiagonals
[in] diag_offset magma_int_t* array containing the offsets

(length offsets+1)

Parameters:
[in] diag_vals float* array containing the values

(length offsets+1)

Parameters:
[out] A magma_s_matrix* matrix to generate
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smLdiagadd ( magma_s_matrix *  L,
magma_queue_t  queue 
)

Checks for a lower triangular matrix whether it is strictly lower triangular and in the negative case adds a unit diagonal.

It does this in-place.

Parameters:
[in,out] L magma_s_matrix* sparse matrix in CSR
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smlumerge ( magma_s_matrix  L,
magma_s_matrix  U,
magma_s_matrix *  A,
magma_queue_t  queue 
)

Takes an strictly lower triangular matrix L and an upper triangular matrix U and merges them into a matrix A containing the upper and lower triangular parts.

Parameters:
[in] L magma_s_matrix input strictly lower triangular matrix L
[in] U magma_s_matrix input upper triangular matrix U
[out] A magma_s_matrix* output matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smscale ( magma_s_matrix *  A,
magma_scale_t  scaling,
magma_queue_t  queue 
)

Scales a matrix.

Parameters:
[in,out] A magma_s_matrix* input/output matrix
[in] scaling magma_scale_t scaling type (unit rownorm / unit diagonal)
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smtransfer ( magma_s_matrix  A,
magma_s_matrix *  B,
magma_location_t  src,
magma_location_t  dst,
magma_queue_t  queue 
)

Copies a matrix from memory location src to memory location dst.

Parameters:
[in] A magma_s_matrix sparse matrix A
[out] B magma_s_matrix* copy of A
[in] src magma_location_t original location A
[in] dst magma_location_t location of the copy of A
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_smtranspose ( magma_s_matrix  A,
magma_s_matrix *  B,
magma_queue_t  queue 
)

Interface to cuSPARSE transpose.

Parameters:
[in] A magma_s_matrix input matrix (CSR)
[out] B magma_s_matrix* output matrix (CSR)
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_snonlinres ( magma_s_matrix  A,
magma_s_matrix  L,
magma_s_matrix  U,
magma_s_matrix *  LU,
real_Double_t *  res,
magma_queue_t  queue 
)

Computes the nonlinear residual A - LU and returns the difference as well es the Frobenius norm of the difference.

Parameters:
[in] A magma_s_matrix input sparse matrix in CSR
[in] L magma_s_matrix input sparse matrix in CSR
[in] U magma_s_matrix input sparse matrix in CSR
[out] LU magma_s_matrix* output sparse matrix in A-LU in CSR
[out] res real_Double_t* Frobenius norm of difference
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sparse_opts ( int  argc,
char **  argv,
magma_sopts *  opts,
int *  matrices,
magma_queue_t  queue 
)

Parses input options for a solver.

Parameters:
[in] argc int command line input
[in] argv char** command line input
[in,out] opts magma_sopts * magma solver options
[out] matrices int counter how many linear systems to process
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sprint_csr ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
float **  val,
magma_index_t **  row,
magma_index_t **  col,
magma_queue_t  queue 
)

Prints a CSR matrix in CSR format.

Parameters:
[in] n_row magma_int_t* number of rows in matrix
[in] n_col magma_int_t* number of columns in matrix
[in] nnz magma_int_t* number of nonzeros in matrix
[in] val float** value array of CSR
[in] row magma_index_t** row pointer of CSR
[in] col magma_index_t** column indices of CSR
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sprint_csr_mtx ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
float **  val,
magma_index_t **  row,
magma_index_t **  col,
magma_order_t  MajorType,
magma_queue_t  queue 
)

Prints a CSR matrix in Matrix Market format.

Parameters:
[in] n_row magma_int_t* number of rows in matrix
[in] n_col magma_int_t* number of columns in matrix
[in] nnz magma_int_t* number of nonzeros in matrix
[in] val float** value array of CSR
[in] row magma_index_t** row pointer of CSR
[in] col magma_index_t** column indices of CSR
[in] MajorType magma_index_t Row or Column sort default: 0 = RowMajor, 1 = ColMajor
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sprint_matrix ( magma_s_matrix  A,
magma_queue_t  queue 
)

Prints a sparse matrix in CSR format.

Parameters:
[in] A magma_s_matrix sparse matrix in Magma_CSR format
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sprint_vector ( magma_s_matrix  x,
magma_int_t  offset,
magma_int_t  visulen,
magma_queue_t  queue 
)

Visualizes part of a vector of type magma_s_matrix.

With input vector x , offset, visulen, the entries offset - (offset + visulen) of x are visualized.

Parameters:
[in] x magma_s_matrix vector to visualize
[in] offset magma_int_t start inex of visualization
[in] visulen magma_int_t number of entries to visualize
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sresidual ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix  x,
float *  res,
magma_queue_t  queue 
)

Computes the residual ||b-Ax|| for a solution approximation x.

Parameters:
[in] A magma_s_matrix input matrix A
[in] b magma_s_matrix RHS b
[in] x magma_s_matrix solution approximation
[out] res float* return residual
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_sresidualvec ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix  x,
magma_s_matrix *  r,
float *  res,
magma_queue_t  queue 
)

Computes the residual r = b-Ax for a solution approximation x.

It returns both, the actual residual and the residual vector

Parameters:
[in] A magma_s_matrix input matrix A
[in] b magma_s_matrix RHS b
[in] x magma_s_matrix solution approximation
[in,out] r magma_s_matrix* residual vector
[out] res float* return residual
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_srowentries ( magma_s_matrix *  A,
magma_queue_t  queue 
)

Checks the maximal number of nonzeros in a row of matrix A.

Inserts the data into max_nnz_row.

Parameters:
[in,out] A magma_s_matrix* sparse matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_ssolverinfo ( magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Prints information about a previously called solver.

Parameters:
[in] solver_par magma_s_solver_par* structure containing all solver information
[in,out] precond_par magma_s_preconditioner* structure containing all preconditioner information
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_ssolverinfo_free ( magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Frees any memory assocoiated with the verbose mode of solver_par.

The other values are set to default.

Parameters:
[in,out] solver_par magma_s_solver_par* structure containing all solver information
[in,out] precond_par magma_s_preconditioner* structure containing all preconditioner information
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_ssolverinfo_init ( magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Initializes all solver and preconditioner parameters.

Parameters:
[in,out] solver_par magma_s_solver_par* structure containing all solver information
[in,out] precond_par magma_s_preconditioner* structure containing all preconditioner information
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_ssymbilu ( magma_s_matrix *  A,
magma_int_t  levels,
magma_s_matrix *  L,
magma_s_matrix *  U,
magma_queue_t  queue 
)

This routine performs a symbolic ILU factorization.

The algorithm is taken from an implementation written by Edmond Chow.

Parameters:
[in,out] A magma_s_matrix* matrix in magma sparse matrix format containing the original matrix on input, and L,U on output
[in] levels magma_magma_int_t_t fill in level
[out] L magma_s_matrix* output lower triangular matrix in magma sparse matrix format empty on function call
[out] U magma_s_matrix* output upper triangular matrix in magma sparse matrix format empty on function call
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svget ( magma_s_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
float **  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back.

Parameters:
[in] v magma_s_matrix magma vector
[out] m magma_int_t number of rows
[out] n magma_int_t number of columns
[out] val float* array containing vector entries
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svget_dev ( magma_s_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
magmaFloat_ptr *  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back (located on DEV).

Parameters:
[in] v magma_s_matrix magma vector
[out] m magma_int_t number of rows
[out] n magma_int_t number of columns
[out] val magmaFloat_ptr array containing vector entries
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svinit ( magma_s_matrix *  x,
magma_location_t  mem_loc,
magma_int_t  num_rows,
magma_int_t  num_cols,
float  values,
magma_queue_t  queue 
)

Allocates memory for magma_s_matrix and initializes it with the passed value.

Parameters:
[out] x magma_s_matrix* vector to initialize
[in] mem_loc magma_location_t memory for vector
[in] num_rows magma_int_t desired length of vector
[in] values float entries in vector
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svread ( magma_s_matrix *  x,
magma_int_t  length,
char *  filename,
magma_queue_t  queue 
)

Reads in a float vector of length "length".

Parameters:
[out] x magma_s_matrix * vector to read in
[in] length magma_int_t length of vector
[in] filename char* file where vector is stored
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svset ( magma_int_t  m,
magma_int_t  n,
float *  val,
magma_s_matrix *  v,
magma_queue_t  queue 
)

Passes a vector to MAGMA.

Parameters:
[in] m magma_int_t number of rows
[in] n magma_int_t number of columns
[in] val float* array containing vector entries
[out] v magma_s_matrix* magma vector
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svset_dev ( magma_int_t  m,
magma_int_t  n,
magmaFloat_ptr  val,
magma_s_matrix *  v,
magma_queue_t  queue 
)

Passes a vector to MAGMA (located on DEV).

Parameters:
[in] m magma_int_t number of rows
[in] n magma_int_t number of columns
[in] val magmaFloat_ptr array containing vector entries
[out] v magma_s_matrix* magma vector
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svspread ( magma_s_matrix *  x,
const char *  filename,
magma_queue_t  queue 
)

Reads in a sparse vector-block stored in COO format.

Parameters:
[out] x magma_s_matrix * vector to read in
[in] filename char* file where vector is stored
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_svtranspose ( magma_s_matrix  x,
magma_s_matrix *  y,
magma_queue_t  queue 
)

Transposes a vector from col to row major and vice versa.

Parameters:
[in] x magma_s_matrix input vector
[out] y magma_s_matrix* output vector
[in] queue magma_queue_t Queue to execute in.
magma_int_t magma_swrite_csr_mtx ( magma_s_matrix  A,
magma_order_t  MajorType,
const char *  filename,
magma_queue_t  queue 
)

Writes a CSR matrix to a file using Matrix Market format.

Parameters:
[in] n_row magma_int_t number of rows in matrix
[in] n_col magma_int_t number of columns in matrix
[in] nnz magma_int_t number of nonzeros in matrix
[in] val float** value array of CSR TODO: why are these ** pointers? Wouldn't * pointers work?
[in] row magma_index_t** row pointer of CSR
[in] col magma_index_t** column indices of CSR
[in] MajorType magma_index_t Row or Column sort default: 0 = RowMajor, 1 = ColMajor TODO: use named constants (e.g., MagmaRowMajor), not numbers.
[in] filename const char* output-filname of the mtx matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t read_s_csr_from_mtx ( magma_storage_t *  type,
magma_location_t *  location,
magma_int_t *  n_row,
magma_int_t *  n_col,
magma_int_t *  nnz,
float **  val,
magma_index_t **  row,
magma_index_t **  col,
const char *  filename,
magma_queue_t  queue 
)

Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.

It duplicates the off-diagonal entries in the symmetric case.

Parameters:
[out] type magma_storage_t* storage type of matrix
[out] location magma_location_t* location of matrix
[out] n_row magma_int_t* number of rows in matrix
[out] n_col magma_int_t* number of columns in matrix
[out] nnz magma_int_t* number of nonzeros in matrix
[out] val float** value array of CSR output
[out] row magma_index_t** row pointer of CSR output
[out] col magma_index_t** column indices of CSR output
[in] filename const char* filname of the mtx matrix
[in] queue magma_queue_t Queue to execute in.
magma_int_t s_transpose_csr ( magma_int_t  n_rows,
magma_int_t  n_cols,
magma_int_t  nnz,
float *  values,
magma_index_t *  rowptr,
magma_index_t *  colind,
magma_int_t *  new_n_rows,
magma_int_t *  new_n_cols,
magma_int_t *  new_nnz,
float **  new_values,
magma_index_t **  new_rowptr,
magma_index_t **  new_colind,
magma_queue_t  queue 
)

Transposes a matrix stored in CSR format on the CPU host.

Parameters:
[in] n_rows magma_int_t number of rows in input matrix
[in] n_cols magma_int_t number of columns in input matrix
[in] nnz magma_int_t number of nonzeros in input matrix
[in] values float* value array of input matrix
[in] rowptr magma_index_t* row pointer of input matrix
[in] colind magma_index_t* column indices of input matrix
[in] new_n_rows magma_index_t* number of rows in transposed matrix
[in] new_n_cols magma_index_t* number of columns in transposed matrix
[in] new_nnz magma_index_t* number of nonzeros in transposed matrix
[in] new_values float** value array of transposed matrix
[in] new_rowptr magma_index_t** row pointer of transposed matrix
[in] new_colind magma_index_t** column indices of transposed matrix
[in] queue magma_queue_t Queue to execute in.

Generated on 3 May 2015 for MAGMA by  doxygen 1.6.1