MAGMA  2.0.2
Matrix Algebra for GPU and Multicore Architectures
single-complex precision

Functions

magma_int_t magma_cresidual (magma_c_matrix A, magma_c_matrix b, magma_c_matrix x, float *res, magma_queue_t queue)
 Computes the residual ||b-Ax|| for a solution approximation x. More...
 
magma_int_t magma_cresidual_slice (magma_int_t start, magma_int_t end, magma_c_matrix A, magma_c_matrix b, magma_c_matrix x, float *res, magma_queue_t queue)
 Computes the residual r=||b-Ax|| for the slice r(start:end) for a solution approximation x. More...
 
magma_int_t magma_cresidualvec (magma_c_matrix A, magma_c_matrix b, magma_c_matrix x, magma_c_matrix *r, float *res, magma_queue_t queue)
 Computes the residual r = b-Ax for a solution approximation x. More...
 
magma_int_t magma_c_precond (magma_c_matrix A, magma_c_matrix b, magma_c_matrix *x, magma_c_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. More...
 
magma_int_t magma_c_precondsetup (magma_c_matrix A, magma_c_matrix b, magma_c_solver_par *solver, magma_c_preconditioner *precond, magma_queue_t queue)
 For a given input matrix M and vectors x, y and the preconditioner parameters, the respective preconditioner is preprocessed. More...
 
magma_int_t magma_c_applyprecond (magma_c_matrix A, magma_c_matrix b, magma_c_matrix *x, magma_c_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. More...
 
magma_int_t magma_c_applyprecond_left (magma_trans_t trans, magma_c_matrix A, magma_c_matrix b, magma_c_matrix *x, magma_c_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. More...
 
magma_int_t magma_c_applyprecond_right (magma_trans_t trans, magma_c_matrix A, magma_c_matrix b, magma_c_matrix *x, magma_c_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. More...
 
magma_int_t magma_c_solver (magma_c_matrix A, magma_c_matrix b, magma_c_matrix *x, magma_copts *zopts, magma_queue_t queue)
 Allows the user to choose a solver. More...
 
magma_int_t magma_capplycustomprecond_l (magma_c_matrix b, magma_c_matrix *x, magma_c_preconditioner *precond, magma_queue_t queue)
 This is an interface to the left solve for any custom preconditioner. More...
 
magma_int_t magma_capplycustomprecond_r (magma_c_matrix b, magma_c_matrix *x, magma_c_preconditioner *precond, magma_queue_t queue)
 This is an interface to the right solve for any custom preconditioner. More...
 
magma_int_t magma_ccsrsplit (magma_int_t offset, magma_int_t bsize, magma_c_matrix A, magma_c_matrix *D, magma_c_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. More...
 
magma_int_t magma_cdomainoverlap (magma_index_t num_rows, magma_int_t *num_indices, magma_index_t *rowptr, magma_index_t *colidx, magma_index_t *x, magma_queue_t queue)
 Generates the update list. More...
 
magma_int_t magma_cmfree (magma_c_matrix *A, magma_queue_t queue)
 Free the memory of a magma_c_matrix. More...
 
magma_int_t magma_cfrobenius (magma_c_matrix A, magma_c_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B. More...
 
magma_int_t magma_cnonlinres (magma_c_matrix A, magma_c_matrix L, magma_c_matrix U, magma_c_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. More...
 
magma_int_t magma_cilures (magma_c_matrix A, magma_c_matrix L, magma_c_matrix U, magma_c_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. More...
 
magma_int_t magma_cicres (magma_c_matrix A, magma_c_matrix C, magma_c_matrix CT, magma_c_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. More...
 
magma_int_t magma_cinitguess (magma_c_matrix A, magma_c_matrix *L, magma_c_matrix *U, magma_queue_t queue)
 Computes an initial guess for the iterative ILU/IC. More...
 
magma_int_t magma_cinitrecursiveLU (magma_c_matrix A, magma_c_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. More...
 
magma_int_t magma_cmLdiagadd (magma_c_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. More...
 
magma_int_t magma_crowentries (magma_c_matrix *A, magma_queue_t queue)
 Checks the maximal number of nonzeros in a row of matrix A. More...
 
magma_int_t magma_cdiameter (magma_c_matrix *A, magma_queue_t queue)
 Computes the diameter of a sparse matrix and stores the value in diameter. More...
 
magma_int_t magma_c_csr_compressor (magmaFloatComplex **val, magma_index_t **row, magma_index_t **col, magmaFloatComplex **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. More...
 
magma_int_t magma_cmconvert (magma_c_matrix A, magma_c_matrix *B, magma_storage_t old_format, magma_storage_t new_format, magma_queue_t queue)
 Converter between different sparse storage formats. More...
 
magma_int_t magma_cmcsrcompressor (magma_c_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix. More...
 
magma_int_t magma_ccsrset (magma_int_t m, magma_int_t n, magma_index_t *row, magma_index_t *col, magmaFloatComplex *val, magma_c_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA. More...
 
magma_int_t magma_ccsrget (magma_c_matrix A, magma_int_t *m, magma_int_t *n, magma_index_t **row, magma_index_t **col, magmaFloatComplex **val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure. More...
 
magma_int_t magma_ccsrset_gpu (magma_int_t m, magma_int_t n, magmaIndex_ptr row, magmaIndex_ptr col, magmaFloatComplex_ptr val, magma_c_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA (located on DEV). More...
 
magma_int_t magma_ccsrget_gpu (magma_c_matrix A, magma_int_t *m, magma_int_t *n, magmaIndex_ptr *row, magmaIndex_ptr *col, magmaFloatComplex_ptr *val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure (located on DEV). More...
 
magma_int_t magma_cmdiff (magma_c_matrix A, magma_c_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B. More...
 
magma_int_t magma_cmgenerator (magma_int_t n, magma_int_t offdiags, magma_index_t *diag_offset, magmaFloatComplex *diag_vals, magma_c_matrix *A, magma_queue_t queue)
 Generate a symmetric n x n CSR matrix for a stencil. More...
 
magma_int_t magma_cm_27stencil (magma_int_t n, magma_c_matrix *A, magma_queue_t queue)
 Generate a 27-point stencil for a 3D FD discretization. More...
 
magma_int_t magma_cm_5stencil (magma_int_t n, magma_c_matrix *A, magma_queue_t queue)
 Generate a 5-point stencil for a 2D FD discretization. More...
 
magma_int_t magma_csymbilu (magma_c_matrix *A, magma_int_t levels, magma_c_matrix *L, magma_c_matrix *U, magma_queue_t queue)
 This routine performs a symbolic ILU factorization. More...
 
magma_int_t read_c_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, magmaFloatComplex **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. More...
 
magma_int_t magma_cwrite_csr_mtx (magma_c_matrix A, magma_order_t MajorType, const char *filename, magma_queue_t queue)
 Writes a CSR matrix to a file using Matrix Market format. More...
 
magma_int_t magma_cprint_csr_mtx (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, magmaFloatComplex **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. More...
 
magma_int_t magma_cprint_csr (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, magmaFloatComplex **val, magma_index_t **row, magma_index_t **col, magma_queue_t queue)
 Prints a CSR matrix in CSR format. More...
 
magma_int_t magma_cprint_matrix (magma_c_matrix A, magma_queue_t queue)
 Prints a sparse matrix in CSR format. More...
 
magma_int_t magma_c_csr_mtx (magma_c_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. More...
 
magma_int_t magma_c_csr_mtxsymm (magma_c_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. More...
 
magma_int_t magma_cmlumerge (magma_c_matrix L, magma_c_matrix U, magma_c_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. More...
 
magma_int_t magma_cmscale (magma_c_matrix *A, magma_scale_t scaling, magma_queue_t queue)
 Scales a matrix. More...
 
magma_int_t magma_cmdiagadd (magma_c_matrix *A, magmaFloatComplex add, magma_queue_t queue)
 Adds a multiple of the Identity matrix to a matrix: A = A+add * I. More...
 
magma_int_t magma_cmshrink (magma_c_matrix A, magma_c_matrix *B, magma_queue_t queue)
 Shrinks a non-square matrix (m < n) to the smaller dimension. More...
 
magma_int_t magma_cmslice (magma_int_t num_slices, magma_int_t slice, magma_c_matrix A, magma_c_matrix *B, magma_c_matrix *ALOC, magma_c_matrix *ANLOC, magma_index_t *comm_i, magmaFloatComplex *comm_v, magma_int_t *start, magma_int_t *end, magma_queue_t queue)
 Takes a matrix and extracts a slice for solving the system in parallel: More...
 
magma_int_t magma_cmtransfer (magma_c_matrix A, magma_c_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. More...
 
magma_int_t c_transpose_csr (magma_int_t n_rows, magma_int_t n_cols, magma_int_t nnz, magmaFloatComplex *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, magmaFloatComplex **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. More...
 
magma_int_t magma_cmtranspose (magma_c_matrix A, magma_c_matrix *B, magma_queue_t queue)
 Interface to cuSPARSE transpose. More...
 
magma_int_t magma_c_cucsrtranspose (magma_c_matrix A, magma_c_matrix *B, magma_queue_t queue)
 Helper function to transpose CSR matrix. More...
 
magma_int_t magma_cmtransposeconjugate (magma_c_matrix A, magma_c_matrix *B, magma_queue_t queue)
 This function forms the transpose conjugate of a matrix. More...
 
magma_int_t magma_csolverinfo (magma_c_solver_par *solver_par, magma_c_preconditioner *precond_par, magma_queue_t queue)
 Prints information about a previously called solver. More...
 
magma_int_t magma_csolverinfo_free (magma_c_solver_par *solver_par, magma_c_preconditioner *precond_par, magma_queue_t queue)
 Frees any memory assocoiated with the verbose mode of solver_par. More...
 
magma_int_t magma_csolverinfo_init (magma_c_solver_par *solver_par, magma_c_preconditioner *precond_par, magma_queue_t queue)
 Initializes all solver and preconditioner parameters. More...
 
magma_int_t magma_ceigensolverinfo_init (magma_c_solver_par *solver_par, magma_queue_t queue)
 Initializes space for eigensolvers. More...
 
magma_int_t magma_csort (magmaFloatComplex *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of values in increasing order. More...
 
magma_int_t magma_cmsort (magmaFloatComplex *x, magma_index_t *col, magma_index_t *row, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of values in increasing order. More...
 
magma_int_t magma_cindexsort (magma_index_t *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers in increasing order. More...
 
magma_int_t magma_cindexsortval (magma_index_t *x, magmaFloatComplex *y, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers, updates a respective array of values. More...
 
magma_int_t magma_cmorderstatistics (magmaFloatComplex *val, magma_index_t *col, magma_index_t *row, magma_int_t length, magma_int_t k, magma_int_t r, magmaFloatComplex *element, magma_queue_t queue)
 Identifies the kth smallest/largest element in an array and reorders such that these elements come to the front. More...
 
magma_int_t magma_corderstatistics (magmaFloatComplex *val, magma_int_t length, magma_int_t k, magma_int_t r, magmaFloatComplex *element, magma_queue_t queue)
 Identifies the kth smallest/largest element in an array. More...
 
magma_int_t magma_cparse_opts (int argc, char **argv, magma_copts *opts, int *matrices, magma_queue_t queue)
 Parses input options for a solver. More...
 
magma_int_t magma_cvinit (magma_c_matrix *x, magma_location_t mem_loc, magma_int_t num_rows, magma_int_t num_cols, magmaFloatComplex values, magma_queue_t queue)
 Allocates memory for magma_c_matrix and initializes it with the passed value. More...
 
magma_int_t magma_cprint_vector (magma_c_matrix x, magma_int_t offset, magma_int_t visulen, magma_queue_t queue)
 Visualizes part of a vector of type magma_c_matrix. More...
 
magma_int_t magma_cvread (magma_c_matrix *x, magma_int_t length, char *filename, magma_queue_t queue)
 Reads in a float vector of length "length". More...
 
magma_int_t magma_cvspread (magma_c_matrix *x, const char *filename, magma_queue_t queue)
 Reads in a sparse vector-block stored in COO format. More...
 
magma_int_t magma_cvset (magma_int_t m, magma_int_t n, magmaFloatComplex *val, magma_c_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA. More...
 
magma_int_t magma_cvget (magma_c_matrix v, magma_int_t *m, magma_int_t *n, magmaFloatComplex **val, magma_queue_t queue)
 Passes a MAGMA vector back. More...
 
magma_int_t magma_cvset_dev (magma_int_t m, magma_int_t n, magmaFloatComplex_ptr val, magma_c_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA (located on DEV). More...
 
magma_int_t magma_cvget_dev (magma_c_matrix v, magma_int_t *m, magma_int_t *n, magmaFloatComplex_ptr *val, magma_queue_t queue)
 Passes a MAGMA vector back (located on DEV). More...
 
magma_int_t magma_cvtranspose (magma_c_matrix x, magma_c_matrix *y, magma_queue_t queue)
 Transposes a vector from col to row major and vice versa. More...
 
magma_int_t magma_vector_clag2z (magma_c_vector x, magma_z_matrix *y, magma_queue_t queue)
 convertes magma_c_vector from C to Z More...
 
magma_int_t magma_sparse_matrix_clag2z (magma_c_sparse_matrix A, magma_z_matrix *B, magma_queue_t queue)
 convertes magma_c_sparse_matrix from C to Z More...
 

Detailed Description

Function Documentation

magma_int_t c_transpose_csr ( magma_int_t  n_rows,
magma_int_t  n_cols,
magma_int_t  nnz,
magmaFloatComplex *  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,
magmaFloatComplex **  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_rowsmagma_int_t number of rows in input matrix
[in]n_colsmagma_int_t number of columns in input matrix
[in]nnzmagma_int_t number of nonzeros in input matrix
[in]valuesmagmaFloatComplex* value array of input matrix
[in]rowptrmagma_index_t* row pointer of input matrix
[in]colindmagma_index_t* column indices of input matrix
[in]new_n_rowsmagma_index_t* number of rows in transposed matrix
[in]new_n_colsmagma_index_t* number of columns in transposed matrix
[in]new_nnzmagma_index_t* number of nonzeros in transposed matrix
[in]new_valuesmagmaFloatComplex** value array of transposed matrix
[in]new_rowptrmagma_index_t** row pointer of transposed matrix
[in]new_colindmagma_index_t** column indices of transposed matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_applyprecond ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]Amagma_c_matrix sparse matrix A
[in]bmagma_c_matrix input vector b
[in,out]xmagma_c_matrix* output vector x
[in]precondmagma_c_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_applyprecond_left ( magma_trans_t  trans,
magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]Amagma_c_matrix sparse matrix A
[in]bmagma_c_matrix input vector b
[in,out]xmagma_c_matrix* output vector x
[in]precondmagma_c_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_applyprecond_right ( magma_trans_t  trans,
magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]Amagma_c_matrix sparse matrix A
[in]bmagma_c_matrix input vector b
[in,out]xmagma_c_matrix* output vector x
[in]precondmagma_c_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_csr_compressor ( magmaFloatComplex **  val,
magma_index_t **  row,
magma_index_t **  col,
magmaFloatComplex **  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]valmagmaFloatComplex** input val pointer to compress
[in]rowmagma_int_t** input row pointer to modify
[in]colmagma_int_t** input col pointer to compress
[in]valnmagmaFloatComplex** output val pointer
[out]rownmagma_int_t** output row pointer
[out]colnmagma_int_t** output col pointer
[out]nmagma_int_t* number of rows in matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_csr_mtx ( magma_c_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]Amagma_c_matrix* matrix in magma sparse matrix format
[in]filenameconst char* filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_csr_mtxsymm ( magma_c_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]Amagma_c_matrix* matrix in magma sparse matrix format
[in]filenameconst char* filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_cucsrtranspose ( magma_c_matrix  A,
magma_c_matrix *  B,
magma_queue_t  queue 
)

Helper function to transpose CSR matrix.

Using the CUSPARSE CSR2CSC function.

Parameters
[in]Amagma_c_matrix input matrix (CSR)
[out]Bmagma_c_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_precond ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]Amagma_c_matrix sparse matrix A
[in]bmagma_c_matrix input vector b
[in]xmagma_c_matrix* output vector x
[in,out]precondmagma_c_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_precondsetup ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_solver_par *  solver,
magma_c_preconditioner *  precond,
magma_queue_t  queue 
)

For a given input matrix M 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]Amagma_c_matrix sparse matrix M
[in]bmagma_c_matrix input vector y
[in,out]precondmagma_c_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_c_solver ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix *  x,
magma_copts *  zopts,
magma_queue_t  queue 
)

Allows the user to choose a solver.

Parameters
[in]Amagma_c_matrix sparse matrix A
[in]bmagma_c_matrix input vector b
[in]xmagma_c_matrix* output vector x
[in]zoptsmagma_copts options for solver and preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_capplycustomprecond_l ( magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]bmagma_c_matrix RHS
[in,out]xmagma_c_matrix* vector to precondition
[in,out]precondmagma_c_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_capplycustomprecond_r ( magma_c_matrix  b,
magma_c_matrix *  x,
magma_c_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]bmagma_c_matrix RHS
[in,out]xmagma_c_matrix* vector to precondition
[in,out]precondmagma_c_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ccsrget ( magma_c_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magma_index_t **  row,
magma_index_t **  col,
magmaFloatComplex **  val,
magma_queue_t  queue 
)

Passes a MAGMA matrix to CSR structure.

Parameters
[in]Amagma_c_matrix magma sparse matrix in CSR format
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]rowmagma_index_t* row pointer
[out]colmagma_index_t* column indices
[out]valmagmaFloatComplex* array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ccsrget_gpu ( magma_c_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magmaIndex_ptr *  row,
magmaIndex_ptr *  col,
magmaFloatComplex_ptr *  val,
magma_queue_t  queue 
)

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

Parameters
[in]Amagma_c_matrix magma sparse matrix in CSR format
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]rowmagmaIndex_ptr row pointer
[out]colmagmaIndex_ptr column indices
[out]valmagmaFloatComplex_ptr array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ccsrset ( magma_int_t  m,
magma_int_t  n,
magma_index_t *  row,
magma_index_t *  col,
magmaFloatComplex *  val,
magma_c_matrix *  A,
magma_queue_t  queue 
)

Passes a CSR matrix to MAGMA.

Parameters
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]rowmagma_index_t* row pointer
[in]colmagma_index_t* column indices
[in]valmagmaFloatComplex* array containing matrix entries
[out]Amagma_c_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ccsrset_gpu ( magma_int_t  m,
magma_int_t  n,
magmaIndex_ptr  row,
magmaIndex_ptr  col,
magmaFloatComplex_ptr  val,
magma_c_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]rowmagmaIndex_ptr row pointer
[in]colmagmaIndex_ptr column indices
[in]valmagmaFloatComplex_ptr array containing matrix entries
[out]Amagma_c_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ccsrsplit ( magma_int_t  offset,
magma_int_t  bsize,
magma_c_matrix  A,
magma_c_matrix *  D,
magma_c_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]offsetmagma_int_t size of the first block
[in]bsizemagma_int_t size of the diagonal blocks
[in]Amagma_c_matrix CSR input matrix
[out]Dmagma_c_matrix* CSR matrix containing diagonal blocks
[out]Rmagma_c_matrix* CSR matrix containing rest
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cdiameter ( magma_c_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in,out]Amagma_c_matrix* sparse matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cdomainoverlap ( magma_index_t  num_rows,
magma_int_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]xmagma_index_t* array to sort
[in]num_rowsmagma_int_t number of rows in matrix
[out]num_indicesmagma_int_t* number of indices in array
[in]rowptrmagma_index_t* rowpointer of matrix
[in]colidxmagma_index_t* colindices of matrix
[in]xmagma_index_t* array containing indices for domain overlap
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ceigensolverinfo_init ( magma_c_solver_par *  solver_par,
magma_queue_t  queue 
)

Initializes space for eigensolvers.

Parameters
[in,out]solver_parmagma_c_solver_par* structure containing all solver information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cfrobenius ( magma_c_matrix  A,
magma_c_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]Amagma_c_matrix sparse matrix in CSR
[in]Bmagma_c_matrix sparse matrix in CSR
[out]resreal_Double_t* Frobenius norm of difference
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cicres ( magma_c_matrix  A,
magma_c_matrix  C,
magma_c_matrix  CT,
magma_c_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]Amagma_c_matrix input sparse matrix in CSR
[in]Cmagma_c_matrix input sparse matrix in CSR
[in]CTmagma_c_matrix input sparse matrix in CSR
[in]LUmagma_c_matrix* output sparse matrix in A-LU in CSR
[out]resreal_Double_t* IC residual
[out]nonlinresreal_Double_t* nonlinear residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cilures ( magma_c_matrix  A,
magma_c_matrix  L,
magma_c_matrix  U,
magma_c_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]Amagma_c_matrix input sparse matrix in CSR
[in]Lmagma_c_matrix input sparse matrix in CSR
[in]Umagma_c_matrix input sparse matrix in CSR
[out]LUmagma_c_matrix* output sparse matrix in A-LU in CSR
[out]resreal_Double_t* Frobenius norm of difference
[out]nonlinresreal_Double_t* Frobenius norm of difference
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cindexsort ( magma_index_t *  x,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of integers in increasing order.

Parameters
[in,out]xmagma_index_t* array to sort
[in]firstmagma_int_t pointer to first element
[in]lastmagma_int_t pointer to last element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cindexsortval ( magma_index_t *  x,
magmaFloatComplex *  y,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of integers, updates a respective array of values.

Parameters
[in,out]xmagma_index_t* array to sort
[in,out]ymagmaFloatComplex* array to sort
[in]firstmagma_int_t pointer to first element
[in]lastmagma_int_t pointer to last element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cinitguess ( magma_c_matrix  A,
magma_c_matrix *  L,
magma_c_matrix *  U,
magma_queue_t  queue 
)

Computes an initial guess for the iterative ILU/IC.

Parameters
[in]Amagma_c_matrix sparse matrix in CSR
[out]Lmagma_c_matrix* sparse matrix in CSR
[out]Umagma_c_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cinitrecursiveLU ( magma_c_matrix  A,
magma_c_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]Amagma_c_matrix* sparse matrix in CSR
[out]Bmagma_c_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cm_27stencil ( magma_int_t  n,
magma_c_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in]nmagma_int_t number of rows
[out]Amagma_c_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cm_5stencil ( magma_int_t  n,
magma_c_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in]nmagma_int_t number of rows
[out]Amagma_c_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmconvert ( magma_c_matrix  A,
magma_c_matrix *  B,
magma_storage_t  old_format,
magma_storage_t  new_format,
magma_queue_t  queue 
)

Converter between different sparse storage formats.

Parameters
[in]Amagma_c_matrix sparse matrix A
[out]Bmagma_c_matrix* copy of A in new format
[in]old_formatmagma_storage_t original storage format
[in]new_formatmagma_storage_t new storage format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmcsrcompressor ( magma_c_matrix *  A,
magma_queue_t  queue 
)

Removes zeros in a CSR matrix.

Parameters
[in,out]Amagma_c_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmdiagadd ( magma_c_matrix *  A,
magmaFloatComplex  add,
magma_queue_t  queue 
)

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

Parameters
[in,out]Amagma_c_matrix* input/output matrix
[in]addmagmaFloatComplex scaling for the identity matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmdiff ( magma_c_matrix  A,
magma_c_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]Amagma_c_matrix sparse matrix in CSR
[in]Bmagma_c_matrix sparse matrix in CSR
[out]resreal_Double_t* residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmfree ( magma_c_matrix *  A,
magma_queue_t  queue 
)

Free the memory of a magma_c_matrix.

Parameters
[in,out]Amagma_c_matrix* matrix to free
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmgenerator ( magma_int_t  n,
magma_int_t  offdiags,
magma_index_t *  diag_offset,
magmaFloatComplex *  diag_vals,
magma_c_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in]nmagma_int_t number of rows
[in]offdiagsmagma_int_t number of offdiagonals
[in]diag_offsetmagma_int_t* array containing the offsets
                            (length offsets+1)
[in]diag_valsmagmaFloatComplex* array containing the values
                            (length offsets+1)
[out]Amagma_c_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmLdiagadd ( magma_c_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]Lmagma_c_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmlumerge ( magma_c_matrix  L,
magma_c_matrix  U,
magma_c_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]Lmagma_c_matrix input strictly lower triangular matrix L
[in]Umagma_c_matrix input upper triangular matrix U
[out]Amagma_c_matrix* output matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmorderstatistics ( magmaFloatComplex *  val,
magma_index_t *  col,
magma_index_t *  row,
magma_int_t  length,
magma_int_t  k,
magma_int_t  r,
magmaFloatComplex *  element,
magma_queue_t  queue 
)

Identifies the kth smallest/largest element in an array and reorders such that these elements come to the front.

The related arrays col and row are also reordered.

Parameters
[in,out]valmagmaFloatComplex* Target array, will be modified during operation.
[in,out]colmagma_index_t* Target array, will be modified during operation.
[in,out]rowmagma_index_t* Target array, will be modified during operation.
[in]lengthmagma_int_t Length of the target array.
[in]kmagma_int_t Element to be identified (largest/smallest).
[in]rmagma_int_t rule how to sort: '1' -> largest, '0' -> smallest
[out]elementmagmaFloatComplex* location of the respective element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmscale ( magma_c_matrix *  A,
magma_scale_t  scaling,
magma_queue_t  queue 
)

Scales a matrix.

Parameters
[in,out]Amagma_c_matrix* input/output matrix
[in]scalingmagma_scale_t scaling type (unit rownorm / unit diagonal)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmshrink ( magma_c_matrix  A,
magma_c_matrix *  B,
magma_queue_t  queue 
)

Shrinks a non-square matrix (m < n) to the smaller dimension.

Parameters
[in]Amagma_c_matrix sparse matrix A
[out]Bmagma_c_matrix* sparse matrix A
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmslice ( magma_int_t  num_slices,
magma_int_t  slice,
magma_c_matrix  A,
magma_c_matrix *  B,
magma_c_matrix *  ALOC,
magma_c_matrix *  ANLOC,
magma_index_t *  comm_i,
magmaFloatComplex *  comm_v,
magma_int_t *  start,
magma_int_t *  end,
magma_queue_t  queue 
)

Takes a matrix and extracts a slice for solving the system in parallel:

B = A( i:i+n, : ) and ALOC = A(i:i+n,i:i+n) and ANLOCA(0:start - end:n,:)

B is of size n x n, ALOC of size end-start x end-start, ANLOC of size end-start x n

The last slice might be smaller. For the non-local parts, B is the identity. comm contains 1ess in the locations that are non-local but needed to solve local system.

Parameters
[in]num_slicesmagma_int_t number of slices
[in]slicemagma_int_t slice id (0.. num_slices-1)
[in]Amagma_c_matrix sparse matrix in CSR
[out]Bmagma_c_matrix* sparse matrix in CSR
[out]ALOCmagma_c_matrix* sparse matrix in CSR
[out]ANLOCmagma_c_matrix* sparse matrix in CSR
[in,out]comm_imagma_int_t* communication plan
[in,out]comm_vmagmaFloatComplex* communication plan
[in]queuemagma_queue_t Queue to execute in.
[out]startmagma_int_t* start of slice (row-index)
[out]endmagma_int_t* end of slice (row-index)
magma_int_t magma_cmsort ( magmaFloatComplex *  x,
magma_index_t *  col,
magma_index_t *  row,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of values in increasing order.

Parameters
[in,out]xmagmaFloatComplex* array to sort
[in,out]colmagma_index_t* Target array, will be modified during operation.
[in,out]rowmagma_index_t* Target array, will be modified during operation.
[in]firstmagma_int_t pointer to first element
[in]lastmagma_int_t pointer to last element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmtransfer ( magma_c_matrix  A,
magma_c_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]Amagma_c_matrix sparse matrix A
[out]Bmagma_c_matrix* copy of A
[in]srcmagma_location_t original location A
[in]dstmagma_location_t location of the copy of A
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmtranspose ( magma_c_matrix  A,
magma_c_matrix *  B,
magma_queue_t  queue 
)

Interface to cuSPARSE transpose.

Parameters
[in]Amagma_c_matrix input matrix (CSR)
[out]Bmagma_c_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cmtransposeconjugate ( magma_c_matrix  A,
magma_c_matrix *  B,
magma_queue_t  queue 
)

This function forms the transpose conjugate of a matrix.

For a real-value matrix, the output is the transpose.

Parameters
[in]Amagma_c_matrix input matrix (CSR)
[out]Bmagma_c_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cnonlinres ( magma_c_matrix  A,
magma_c_matrix  L,
magma_c_matrix  U,
magma_c_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]Amagma_c_matrix input sparse matrix in CSR
[in]Lmagma_c_matrix input sparse matrix in CSR
[in]Umagma_c_matrix input sparse matrix in CSR
[out]LUmagma_c_matrix* output sparse matrix in A-LU in CSR
[out]resreal_Double_t* Frobenius norm of difference
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_corderstatistics ( magmaFloatComplex *  val,
magma_int_t  length,
magma_int_t  k,
magma_int_t  r,
magmaFloatComplex *  element,
magma_queue_t  queue 
)

Identifies the kth smallest/largest element in an array.

Parameters
[in,out]valmagmaFloatComplex* Target array, will be modified during operation.
[in]lengthmagma_int_t Length of the target array.
[in]kmagma_int_t Element to be identified (largest/smallest).
[in]rmagma_int_t rule how to sort: '1' -> largest, '0' -> smallest
[out]elementmagmaFloatComplex* location of the respective element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cparse_opts ( int  argc,
char **  argv,
magma_copts *  opts,
int *  matrices,
magma_queue_t  queue 
)

Parses input options for a solver.

Parameters
[in]argcint command line input
[in]argvchar** command line input
[in,out]optsmagma_copts * magma solver options
[out]matricesint counter how many linear systems to process
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cprint_csr ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
magmaFloatComplex **  val,
magma_index_t **  row,
magma_index_t **  col,
magma_queue_t  queue 
)

Prints a CSR matrix in CSR format.

Parameters
[in]n_rowmagma_int_t* number of rows in matrix
[in]n_colmagma_int_t* number of columns in matrix
[in]nnzmagma_int_t* number of nonzeros in matrix
[in]valmagmaFloatComplex** value array of CSR
[in]rowmagma_index_t** row pointer of CSR
[in]colmagma_index_t** column indices of CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cprint_csr_mtx ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
magmaFloatComplex **  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_rowmagma_int_t* number of rows in matrix
[in]n_colmagma_int_t* number of columns in matrix
[in]nnzmagma_int_t* number of nonzeros in matrix
[in]valmagmaFloatComplex** value array of CSR
[in]rowmagma_index_t** row pointer of CSR
[in]colmagma_index_t** column indices of CSR
[in]MajorTypemagma_index_t Row or Column sort default: 0 = RowMajor, 1 = ColMajor
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cprint_matrix ( magma_c_matrix  A,
magma_queue_t  queue 
)

Prints a sparse matrix in CSR format.

Parameters
[in]Amagma_c_matrix sparse matrix in Magma_CSR format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cprint_vector ( magma_c_matrix  x,
magma_int_t  offset,
magma_int_t  visulen,
magma_queue_t  queue 
)

Visualizes part of a vector of type magma_c_matrix.

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

Parameters
[in]xmagma_c_matrix vector to visualize
[in]offsetmagma_int_t start inex of visualization
[in]visulenmagma_int_t number of entries to visualize
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cresidual ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix  x,
float *  res,
magma_queue_t  queue 
)

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

Parameters
[in]Amagma_c_matrix input matrix A
[in]bmagma_c_matrix RHS b
[in]xmagma_c_matrix solution approximation
[out]resmagmaFloatComplex* return residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cresidual_slice ( magma_int_t  start,
magma_int_t  end,
magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix  x,
float *  res,
magma_queue_t  queue 
)

Computes the residual r=||b-Ax|| for the slice r(start:end) for a solution approximation x.

Parameters
[in]startmagma_int_t start of slice (row-index)
[in]endmagma_int_t end of slice (row-index)
[in]Amagma_c_matrix input matrix A
[in]bmagma_c_matrix RHS b
[in]xmagma_c_matrix solution approximation
[out]resmagmaFloatComplex* return residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cresidualvec ( magma_c_matrix  A,
magma_c_matrix  b,
magma_c_matrix  x,
magma_c_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]Amagma_c_matrix input matrix A
[in]bmagma_c_matrix RHS b
[in]xmagma_c_matrix solution approximation
[in,out]rmagma_c_matrix* residual vector
[out]resmagmaFloatComplex* return residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_crowentries ( magma_c_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]Amagma_c_matrix* sparse matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_csolverinfo ( magma_c_solver_par *  solver_par,
magma_c_preconditioner *  precond_par,
magma_queue_t  queue 
)

Prints information about a previously called solver.

Parameters
[in]solver_parmagma_c_solver_par* structure containing all solver information
[in,out]precond_parmagma_c_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_csolverinfo_free ( magma_c_solver_par *  solver_par,
magma_c_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_parmagma_c_solver_par* structure containing all solver information
[in,out]precond_parmagma_c_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_csolverinfo_init ( magma_c_solver_par *  solver_par,
magma_c_preconditioner *  precond_par,
magma_queue_t  queue 
)

Initializes all solver and preconditioner parameters.

Parameters
[in,out]solver_parmagma_c_solver_par* structure containing all solver information
[in,out]precond_parmagma_c_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_csort ( magmaFloatComplex *  x,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of values in increasing order.

Parameters
[in,out]xmagmaFloatComplex* array to sort
[in]firstmagma_int_t pointer to first element
[in]lastmagma_int_t pointer to last element
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_csymbilu ( magma_c_matrix *  A,
magma_int_t  levels,
magma_c_matrix *  L,
magma_c_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]Amagma_c_matrix* matrix in magma sparse matrix format containing the original matrix on input, and L,U on output
[in]levelsmagma_magma_int_t_t fill in level
[out]Lmagma_c_matrix* output lower triangular matrix in magma sparse matrix format empty on function call
[out]Umagma_c_matrix* output upper triangular matrix in magma sparse matrix format empty on function call
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvget ( magma_c_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
magmaFloatComplex **  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back.

Parameters
[in]vmagma_c_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaFloatComplex* array containing vector entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvget_dev ( magma_c_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
magmaFloatComplex_ptr *  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back (located on DEV).

Parameters
[in]vmagma_c_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaFloatComplex_ptr array containing vector entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvinit ( magma_c_matrix *  x,
magma_location_t  mem_loc,
magma_int_t  num_rows,
magma_int_t  num_cols,
magmaFloatComplex  values,
magma_queue_t  queue 
)

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

Parameters
[out]xmagma_c_matrix* vector to initialize
[in]mem_locmagma_location_t memory for vector
[in]num_rowsmagma_int_t desired length of vector
[in]num_colsmagma_int_t desired width of vector-block (columns of dense matrix)
[in]valuesmagmaFloatComplex entries in vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvread ( magma_c_matrix *  x,
magma_int_t  length,
char *  filename,
magma_queue_t  queue 
)

Reads in a float vector of length "length".

Parameters
[out]xmagma_c_matrix * vector to read in
[in]lengthmagma_int_t length of vector
[in]filenamechar* file where vector is stored
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvset ( magma_int_t  m,
magma_int_t  n,
magmaFloatComplex *  val,
magma_c_matrix *  v,
magma_queue_t  queue 
)

Passes a vector to MAGMA.

Parameters
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]valmagmaFloatComplex* array containing vector entries
[out]vmagma_c_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvset_dev ( magma_int_t  m,
magma_int_t  n,
magmaFloatComplex_ptr  val,
magma_c_matrix *  v,
magma_queue_t  queue 
)

Passes a vector to MAGMA (located on DEV).

Parameters
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]valmagmaFloatComplex_ptr array containing vector entries
[out]vmagma_c_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvspread ( magma_c_matrix *  x,
const char *  filename,
magma_queue_t  queue 
)

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

Parameters
[out]xmagma_c_matrix * vector to read in
[in]filenamechar* file where vector is stored
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cvtranspose ( magma_c_matrix  x,
magma_c_matrix *  y,
magma_queue_t  queue 
)

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

Parameters
[in]xmagma_c_matrix input vector
[out]ymagma_c_matrix* output vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cwrite_csr_mtx ( magma_c_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]Amagma_c_matrix matrix to write out
[in]MajorTypemagma_index_t Row or Column sort default: 0 = RowMajor, 1 = ColMajor TODO: use named constants (e.g., MagmaRowMajor), not numbers.
[in]filenameconst char* output-filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_sparse_matrix_clag2z ( magma_c_sparse_matrix  A,
magma_z_matrix *  B,
magma_queue_t  queue 
)

convertes magma_c_sparse_matrix from C to Z

Parameters
Amagma_c_sparse_matrix input matrix descriptor
Bmagma_z_matrix* output matrix descriptor
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_vector_clag2z ( magma_c_vector  x,
magma_z_matrix *  y,
magma_queue_t  queue 
)

convertes magma_c_vector from C to Z

Parameters
[in]xmagma_c_vector input vector descriptor
[out]ymagma_z_matrix* output vector descriptor
[in]queuemagma_queue_t Queue to execute in.
magma_int_t read_c_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,
magmaFloatComplex **  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]typemagma_storage_t* storage type of matrix
[out]locationmagma_location_t* location of matrix
[out]n_rowmagma_int_t* number of rows in matrix
[out]n_colmagma_int_t* number of columns in matrix
[out]nnzmagma_int_t* number of nonzeros in matrix
[out]valmagmaFloatComplex** value array of CSR output
[out]rowmagma_index_t** row pointer of CSR output
[out]colmagma_index_t** column indices of CSR output
[in]filenameconst char* filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.