MAGMA  1.6.3
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
double precision

Functions

magma_int_t magma_dresidual (magma_d_matrix A, magma_d_matrix b, magma_d_matrix x, double *res, magma_queue_t queue)
 Computes the residual ||b-Ax|| for a solution approximation x. More...
 
magma_int_t magma_dresidualvec (magma_d_matrix A, magma_d_matrix b, magma_d_matrix x, magma_d_matrix *r, double *res, magma_queue_t queue)
 Computes the residual r = b-Ax for a solution approximation x. More...
 
magma_int_t magma_d_precond (magma_d_matrix A, magma_d_matrix b, magma_d_matrix *x, magma_d_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_d_precondsetup (magma_d_matrix A, magma_d_matrix b, magma_d_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. More...
 
magma_int_t magma_d_applyprecond (magma_d_matrix A, magma_d_matrix b, magma_d_matrix *x, magma_d_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_d_applyprecond_left (magma_d_matrix A, magma_d_matrix b, magma_d_matrix *x, magma_d_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_d_applyprecond_right (magma_d_matrix A, magma_d_matrix b, magma_d_matrix *x, magma_d_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_d_solver (magma_d_matrix A, magma_d_matrix b, magma_d_matrix *x, magma_dopts *zopts, magma_queue_t queue)
 ALlows the user to choose a solver. More...
 
magma_int_t magma_dapplycustomprecond_l (magma_d_matrix b, magma_d_matrix *x, magma_d_preconditioner *precond, magma_queue_t queue)
 This is an interface to the left solve for any custom preconditioner. More...
 
magma_int_t magma_dapplycustomprecond_r (magma_d_matrix b, magma_d_matrix *x, magma_d_preconditioner *precond, magma_queue_t queue)
 This is an interface to the right solve for any custom preconditioner. More...
 
magma_int_t magma_dcsrsplit (magma_int_t bsize, magma_d_matrix A, magma_d_matrix *D, magma_d_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_ddomainoverlap (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_dmfree (magma_d_matrix *A, magma_queue_t queue)
 Free the memory of a magma_d_matrix. More...
 
magma_int_t magma_dfrobenius (magma_d_matrix A, magma_d_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_dnonlinres (magma_d_matrix A, magma_d_matrix L, magma_d_matrix U, magma_d_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_dilures (magma_d_matrix A, magma_d_matrix L, magma_d_matrix U, magma_d_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_dicres (magma_d_matrix A, magma_d_matrix C, magma_d_matrix CT, magma_d_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_dinitguess (magma_d_matrix A, magma_d_matrix *L, magma_d_matrix *U, magma_queue_t queue)
 Computes an initial guess for the iterative ILU/IC. More...
 
magma_int_t magma_dinitrecursiveLU (magma_d_matrix A, magma_d_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_dmLdiagadd (magma_d_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_drowentries (magma_d_matrix *A, magma_queue_t queue)
 Checks the maximal number of nonzeros in a row of matrix A. More...
 
magma_int_t magma_ddiameter (magma_d_matrix *A, magma_queue_t queue)
 Computes the diameter of a sparse matrix and stores the value in diameter. More...
 
magma_int_t magma_d_csr_compressor (double **val, magma_index_t **row, magma_index_t **col, double **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_dmconvert (magma_d_matrix A, magma_d_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_dmcsrcompressor (magma_d_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix. More...
 
magma_int_t magma_dcsrset (magma_int_t m, magma_int_t n, magma_index_t *row, magma_index_t *col, double *val, magma_d_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA. More...
 
magma_int_t magma_dcsrget (magma_d_matrix A, magma_int_t *m, magma_int_t *n, magma_index_t **row, magma_index_t **col, double **val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure. More...
 
magma_int_t magma_dcsrset_gpu (magma_int_t m, magma_int_t n, magmaIndex_ptr row, magmaIndex_ptr col, magmaDouble_ptr val, magma_d_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA (located on DEV). More...
 
magma_int_t magma_dcsrget_gpu (magma_d_matrix A, magma_int_t *m, magma_int_t *n, magmaIndex_ptr *row, magmaIndex_ptr *col, magmaDouble_ptr *val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure (located on DEV). More...
 
magma_int_t magma_dmdiff (magma_d_matrix A, magma_d_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_dmgenerator (magma_int_t n, magma_int_t offdiags, magma_index_t *diag_offset, double *diag_vals, magma_d_matrix *A, magma_queue_t queue)
 Generate a symmetric n x n CSR matrix for a stencil. More...
 
magma_int_t magma_dm_27stencil (magma_int_t n, magma_d_matrix *A, magma_queue_t queue)
 Generate a 27-point stencil for a 3D FD discretization. More...
 
magma_int_t magma_dm_5stencil (magma_int_t n, magma_d_matrix *A, magma_queue_t queue)
 Generate a 5-point stencil for a 2D FD discretization. More...
 
magma_int_t magma_dsymbilu (magma_d_matrix *A, magma_int_t levels, magma_d_matrix *L, magma_d_matrix *U, magma_queue_t queue)
 This routine performs a symbolic ILU factorization. More...
 
magma_int_t read_d_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, double **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_dwrite_csr_mtx (magma_d_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_dprint_csr_mtx (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, double **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_dprint_csr (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, double **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_dprint_matrix (magma_d_matrix A, magma_queue_t queue)
 Prints a sparse matrix in CSR format. More...
 
magma_int_t magma_d_csr_mtx (magma_d_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_d_csr_mtxsymm (magma_d_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_dmlumerge (magma_d_matrix L, magma_d_matrix U, magma_d_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_dmscale (magma_d_matrix *A, magma_scale_t scaling, magma_queue_t queue)
 Scales a matrix. More...
 
magma_int_t magma_dmdiagadd (magma_d_matrix *A, double add, magma_queue_t queue)
 Adds a multiple of the Identity matrix to a matrix: A = A+add * I. More...
 
magma_int_t magma_dmtransfer (magma_d_matrix A, magma_d_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 d_transpose_csr (magma_int_t n_rows, magma_int_t n_cols, magma_int_t nnz, double *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, double **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_dmtranspose (magma_d_matrix A, magma_d_matrix *B, magma_queue_t queue)
 Interface to cuSPARSE transpose. More...
 
magma_int_t magma_d_cucsrtranspose (magma_d_matrix A, magma_d_matrix *B, magma_queue_t queue)
 Helper function to transpose CSR matrix. More...
 
magma_int_t magma_dsolverinfo (magma_d_solver_par *solver_par, magma_d_preconditioner *precond_par, magma_queue_t queue)
 Prints information about a previously called solver. More...
 
magma_int_t magma_dsolverinfo_free (magma_d_solver_par *solver_par, magma_d_preconditioner *precond_par, magma_queue_t queue)
 Frees any memory assocoiated with the verbose mode of solver_par. More...
 
magma_int_t magma_dsolverinfo_init (magma_d_solver_par *solver_par, magma_d_preconditioner *precond_par, magma_queue_t queue)
 Initializes all solver and preconditioner parameters. More...
 
magma_int_t magma_deigensolverinfo_init (magma_d_solver_par *solver_par, magma_queue_t queue)
 Initializes space for eigensolvers. More...
 
magma_int_t magma_dindexsort (magma_index_t *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers. More...
 
magma_int_t magma_dparse_opts (int argc, char **argv, magma_dopts *opts, int *matrices, magma_queue_t queue)
 Parses input options for a solver. More...
 
magma_int_t magma_dvinit (magma_d_matrix *x, magma_location_t mem_loc, magma_int_t num_rows, magma_int_t num_cols, double values, magma_queue_t queue)
 Allocates memory for magma_d_matrix and initializes it with the passed value. More...
 
magma_int_t magma_dprint_vector (magma_d_matrix x, magma_int_t offset, magma_int_t visulen, magma_queue_t queue)
 Visualizes part of a vector of type magma_d_matrix. More...
 
magma_int_t magma_dvread (magma_d_matrix *x, magma_int_t length, char *filename, magma_queue_t queue)
 Reads in a double vector of length "length". More...
 
magma_int_t magma_dvspread (magma_d_matrix *x, const char *filename, magma_queue_t queue)
 Reads in a sparse vector-block stored in COO format. More...
 
magma_int_t magma_dvset (magma_int_t m, magma_int_t n, double *val, magma_d_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA. More...
 
magma_int_t magma_dvget (magma_d_matrix v, magma_int_t *m, magma_int_t *n, double **val, magma_queue_t queue)
 Passes a MAGMA vector back. More...
 
magma_int_t magma_dvset_dev (magma_int_t m, magma_int_t n, magmaDouble_ptr val, magma_d_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA (located on DEV). More...
 
magma_int_t magma_dvget_dev (magma_d_matrix v, magma_int_t *m, magma_int_t *n, magmaDouble_ptr *val, magma_queue_t queue)
 Passes a MAGMA vector back (located on DEV). More...
 
magma_int_t magma_dvtranspose (magma_d_matrix x, magma_d_matrix *y, magma_queue_t queue)
 Transposes a vector from col to row major and vice versa. More...
 
magma_int_t magma_dlobpcg_res (magma_int_t num_rows, magma_int_t num_vecs, magmaDouble_ptr evalues, magmaDouble_ptr X, magmaDouble_ptr R, magmaDouble_ptr res, magma_queue_t queue)
 This routine computes for Block-LOBPCG, the set of residuals. More...
 
magma_int_t magma_dlobpcg_shift (magma_int_t num_rows, magma_int_t num_vecs, magma_int_t shift, magmaDouble_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. More...
 
magma_int_t magma_dcopyscale (int n, int k, magmaDouble_ptr r, magmaDouble_ptr v, magmaDouble_ptr skp, magma_queue_t queue)
 Computes the correction term of the pipelined GMRES according to P. More...
 
magma_int_t magma_ddiagcheck (magma_d_matrix dA, magma_queue_t queue)
 This routine checks for a CSR matrix whether there exists a zero on the diagonal. More...
 
magma_int_t magma_dmcsrcompressor_gpu (magma_d_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix. More...
 

Detailed Description

Function Documentation

magma_int_t d_transpose_csr ( magma_int_t  n_rows,
magma_int_t  n_cols,
magma_int_t  nnz,
double *  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,
double **  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]valuesdouble* 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_valuesdouble** 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_d_applyprecond ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector b
[in,out]xmagma_d_matrix* output vector x
[in]precondmagma_d_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_applyprecond_left ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector b
[in,out]xmagma_d_matrix* output vector x
[in]precondmagma_d_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_applyprecond_right ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector b
[in,out]xmagma_d_matrix* output vector x
[in]precondmagma_d_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_csr_compressor ( double **  val,
magma_index_t **  row,
magma_index_t **  col,
double **  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]valdouble** input val pointer to compress
[in]rowmagma_int_t** input row pointer to modify
[in]colmagma_int_t** input col pointer to compress
[in]valndouble** 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_d_csr_mtx ( magma_d_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_d_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_d_csr_mtxsymm ( magma_d_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_d_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_d_cucsrtranspose ( magma_d_matrix  A,
magma_d_matrix *  B,
magma_queue_t  queue 
)

Helper function to transpose CSR matrix.

Using the CUSPARSE CSR2CSC function.

Parameters
[in]Amagma_d_matrix input matrix (CSR)
[out]Bmagma_d_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_precond ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector b
[in]xmagma_d_matrix* output vector x
[in,out]precondmagma_d_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_precondsetup ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_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]Amagma_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector y
[in,out]precondmagma_d_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_solver ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix *  x,
magma_dopts *  zopts,
magma_queue_t  queue 
)

ALlows the user to choose a solver.

Parameters
[in]Amagma_d_matrix sparse matrix A
[in]bmagma_d_matrix input vector b
[in]xmagma_d_matrix* output vector x
[in]zoptsmagma_dopts options for solver and preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dapplycustomprecond_l ( magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix RHS
[in,out]xmagma_d_matrix* vector to precondition
[in,out]precondmagma_d_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dapplycustomprecond_r ( magma_d_matrix  b,
magma_d_matrix *  x,
magma_d_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_d_matrix RHS
[in,out]xmagma_d_matrix* vector to precondition
[in,out]precondmagma_d_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcopyscale ( int  n,
int  k,
magmaDouble_ptr  r,
magmaDouble_ptr  v,
magmaDouble_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]nint length of v_i
[in]kint

skp entries v_i^T * r ( without r )

Parameters
[in]rmagmaDouble_ptr vector of length n
[in]vmagmaDouble_ptr vector of length n
[in]skpmagmaDouble_ptr array of parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcsrget ( magma_d_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magma_index_t **  row,
magma_index_t **  col,
double **  val,
magma_queue_t  queue 
)

Passes a MAGMA matrix to CSR structure.

Parameters
[in]Amagma_d_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]valdouble* array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcsrget_gpu ( magma_d_matrix  A,
magma_int_t *  m,
magma_int_t *  n,
magmaIndex_ptr *  row,
magmaIndex_ptr *  col,
magmaDouble_ptr *  val,
magma_queue_t  queue 
)

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

Parameters
[in]Amagma_d_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]valmagmaDouble_ptr array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcsrset ( magma_int_t  m,
magma_int_t  n,
magma_index_t *  row,
magma_index_t *  col,
double *  val,
magma_d_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]valdouble* array containing matrix entries
[out]Amagma_d_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcsrset_gpu ( magma_int_t  m,
magma_int_t  n,
magmaIndex_ptr  row,
magmaIndex_ptr  col,
magmaDouble_ptr  val,
magma_d_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]valmagmaDouble_ptr array containing matrix entries
[out]Amagma_d_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcsrsplit ( magma_int_t  bsize,
magma_d_matrix  A,
magma_d_matrix *  D,
magma_d_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]bsizemagma_int_t size of the diagonal blocks
[in]Amagma_d_matrix CSR input matrix
[out]Dmagma_d_matrix* CSR matrix containing diagonal blocks
[out]Rmagma_d_matrix* CSR matrix containing rest
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ddiagcheck ( magma_d_matrix  dA,
magma_queue_t  queue 
)

This routine checks for a CSR matrix whether there exists a zero on the diagonal.

This can be the diagonal entry missing or an explicit zero.

Parameters
[in]dAmagma_d_matrix matrix in CSR format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ddiameter ( magma_d_matrix *  A,
magma_queue_t  queue 
)

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

Parameters
[in,out]Amagma_d_matrix* sparse matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ddomainoverlap ( 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_deigensolverinfo_init ( magma_d_solver_par *  solver_par,
magma_queue_t  queue 
)

Initializes space for eigensolvers.

Parameters
[in,out]solver_parmagma_d_solver_par* structure containing all solver information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dfrobenius ( magma_d_matrix  A,
magma_d_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_d_matrix sparse matrix in CSR
[in]Bmagma_d_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_dicres ( magma_d_matrix  A,
magma_d_matrix  C,
magma_d_matrix  CT,
magma_d_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_d_matrix input sparse matrix in CSR
[in]Cmagma_d_matrix input sparse matrix in CSR
[in]CTmagma_d_matrix input sparse matrix in CSR
[in]LUmagma_d_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_dilures ( magma_d_matrix  A,
magma_d_matrix  L,
magma_d_matrix  U,
magma_d_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_d_matrix input sparse matrix in CSR
[in]Lmagma_d_matrix input sparse matrix in CSR
[in]Umagma_d_matrix input sparse matrix in CSR
[out]LUmagma_d_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_dindexsort ( magma_index_t *  x,
magma_int_t  first,
magma_int_t  last,
magma_queue_t  queue 
)

Sorts an array of integers.

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_dinitguess ( magma_d_matrix  A,
magma_d_matrix *  L,
magma_d_matrix *  U,
magma_queue_t  queue 
)

Computes an initial guess for the iterative ILU/IC.

Parameters
[in]Amagma_d_matrix sparse matrix in CSR
[out]Lmagma_d_matrix* sparse matrix in CSR
[out]Umagma_d_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dinitrecursiveLU ( magma_d_matrix  A,
magma_d_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_d_matrix* sparse matrix in CSR
[out]Bmagma_d_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dlobpcg_res ( magma_int_t  num_rows,
magma_int_t  num_vecs,
magmaDouble_ptr  evalues,
magmaDouble_ptr  X,
magmaDouble_ptr  R,
magmaDouble_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_daxpy(m, MAGMA_D_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_rowsmagma_int_t number of rows
[in]num_vecsmagma_int_t number of vectors
[in]evaluesmagmaDouble_ptr array of eigenvalues/approximations
[in]XmagmaDouble_ptr block of eigenvector approximations
[in]RmagmaDouble_ptr block of residuals
[in]resmagmaDouble_ptr array of residuals
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dlobpcg_shift ( magma_int_t  num_rows,
magma_int_t  num_vecs,
magma_int_t  shift,
magmaDouble_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_rowsmagma_int_t number of rows
[in]num_vecsmagma_int_t number of vectors
[in]shiftmagma_int_t shift number
[in,out]xmagmaDouble_ptr input/output vector x
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dm_27stencil ( magma_int_t  n,
magma_d_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_d_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dm_5stencil ( magma_int_t  n,
magma_d_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_d_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmconvert ( magma_d_matrix  A,
magma_d_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_d_matrix sparse matrix A
[out]Bmagma_d_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_dmcsrcompressor ( magma_d_matrix *  A,
magma_queue_t  queue 
)

Removes zeros in a CSR matrix.

Parameters
[in,out]Amagma_d_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmcsrcompressor_gpu ( magma_d_matrix *  A,
magma_queue_t  queue 
)

Removes zeros in a CSR matrix.

This is a GPU implementation of the CSR compressor.

Parameters
Amagma_d_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmdiagadd ( magma_d_matrix *  A,
double  add,
magma_queue_t  queue 
)

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

Parameters
[in,out]Amagma_d_matrix* input/output matrix
[in]adddouble scaling for the identity matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmdiff ( magma_d_matrix  A,
magma_d_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_d_matrix sparse matrix in CSR
[in]Bmagma_d_matrix sparse matrix in CSR
[out]resreal_Double_t* residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmfree ( magma_d_matrix *  A,
magma_queue_t  queue 
)

Free the memory of a magma_d_matrix.

Parameters
[in,out]Amagma_d_matrix* matrix to free
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmgenerator ( magma_int_t  n,
magma_int_t  offdiags,
magma_index_t *  diag_offset,
double *  diag_vals,
magma_d_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_valsdouble* array containing the values
                            (length offsets+1)
[out]Amagma_d_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmLdiagadd ( magma_d_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_d_matrix* sparse matrix in CSR
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmlumerge ( magma_d_matrix  L,
magma_d_matrix  U,
magma_d_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_d_matrix input strictly lower triangular matrix L
[in]Umagma_d_matrix input upper triangular matrix U
[out]Amagma_d_matrix* output matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmscale ( magma_d_matrix *  A,
magma_scale_t  scaling,
magma_queue_t  queue 
)

Scales a matrix.

Parameters
[in,out]Amagma_d_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_dmtransfer ( magma_d_matrix  A,
magma_d_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_d_matrix sparse matrix A
[out]Bmagma_d_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_dmtranspose ( magma_d_matrix  A,
magma_d_matrix *  B,
magma_queue_t  queue 
)

Interface to cuSPARSE transpose.

Parameters
[in]Amagma_d_matrix input matrix (CSR)
[out]Bmagma_d_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dnonlinres ( magma_d_matrix  A,
magma_d_matrix  L,
magma_d_matrix  U,
magma_d_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_d_matrix input sparse matrix in CSR
[in]Lmagma_d_matrix input sparse matrix in CSR
[in]Umagma_d_matrix input sparse matrix in CSR
[out]LUmagma_d_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_dparse_opts ( int  argc,
char **  argv,
magma_dopts *  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_dopts * magma solver options
[out]matricesint counter how many linear systems to process
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dprint_csr ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
double **  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]valdouble** 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_dprint_csr_mtx ( magma_int_t  n_row,
magma_int_t  n_col,
magma_int_t  nnz,
double **  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]valdouble** 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_dprint_matrix ( magma_d_matrix  A,
magma_queue_t  queue 
)

Prints a sparse matrix in CSR format.

Parameters
[in]Amagma_d_matrix sparse matrix in Magma_CSR format
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dprint_vector ( magma_d_matrix  x,
magma_int_t  offset,
magma_int_t  visulen,
magma_queue_t  queue 
)

Visualizes part of a vector of type magma_d_matrix.

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

Parameters
[in]xmagma_d_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_dresidual ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix  x,
double *  res,
magma_queue_t  queue 
)

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

Parameters
[in]Amagma_d_matrix input matrix A
[in]bmagma_d_matrix RHS b
[in]xmagma_d_matrix solution approximation
[out]resdouble* return residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dresidualvec ( magma_d_matrix  A,
magma_d_matrix  b,
magma_d_matrix  x,
magma_d_matrix *  r,
double *  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_d_matrix input matrix A
[in]bmagma_d_matrix RHS b
[in]xmagma_d_matrix solution approximation
[in,out]rmagma_d_matrix* residual vector
[out]resdouble* return residual
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_drowentries ( magma_d_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_d_matrix* sparse matrix
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dsolverinfo ( magma_d_solver_par *  solver_par,
magma_d_preconditioner *  precond_par,
magma_queue_t  queue 
)

Prints information about a previously called solver.

Parameters
[in]solver_parmagma_d_solver_par* structure containing all solver information
[in,out]precond_parmagma_d_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dsolverinfo_free ( magma_d_solver_par *  solver_par,
magma_d_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_d_solver_par* structure containing all solver information
[in,out]precond_parmagma_d_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dsolverinfo_init ( magma_d_solver_par *  solver_par,
magma_d_preconditioner *  precond_par,
magma_queue_t  queue 
)

Initializes all solver and preconditioner parameters.

Parameters
[in,out]solver_parmagma_d_solver_par* structure containing all solver information
[in,out]precond_parmagma_d_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dsymbilu ( magma_d_matrix *  A,
magma_int_t  levels,
magma_d_matrix *  L,
magma_d_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_d_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_d_matrix* output lower triangular matrix in magma sparse matrix format empty on function call
[out]Umagma_d_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_dvget ( magma_d_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
double **  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back.

Parameters
[in]vmagma_d_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valdouble* array containing vector entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvget_dev ( magma_d_matrix  v,
magma_int_t *  m,
magma_int_t *  n,
magmaDouble_ptr *  val,
magma_queue_t  queue 
)

Passes a MAGMA vector back (located on DEV).

Parameters
[in]vmagma_d_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaDouble_ptr array containing vector entries
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvinit ( magma_d_matrix *  x,
magma_location_t  mem_loc,
magma_int_t  num_rows,
magma_int_t  num_cols,
double  values,
magma_queue_t  queue 
)

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

Parameters
[out]xmagma_d_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]valuesdouble entries in vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvread ( magma_d_matrix *  x,
magma_int_t  length,
char *  filename,
magma_queue_t  queue 
)

Reads in a double vector of length "length".

Parameters
[out]xmagma_d_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_dvset ( magma_int_t  m,
magma_int_t  n,
double *  val,
magma_d_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]valdouble* array containing vector entries
[out]vmagma_d_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvset_dev ( magma_int_t  m,
magma_int_t  n,
magmaDouble_ptr  val,
magma_d_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]valmagmaDouble_ptr array containing vector entries
[out]vmagma_d_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvspread ( magma_d_matrix *  x,
const char *  filename,
magma_queue_t  queue 
)

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

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

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

Parameters
[in]xmagma_d_matrix input vector
[out]ymagma_d_matrix* output vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dwrite_csr_mtx ( magma_d_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_d_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 read_d_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,
double **  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]valdouble** 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.