MAGMA 2.10.0
Matrix Algebra for GPU and Multicore Architectures
Loading...
Searching...
No Matches
double-complex precision

Functions

magma_int_t magma_z_precond (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is chosen.
 
magma_int_t magma_z_precondsetup (magma_z_matrix A, magma_z_matrix b, magma_z_solver_par *solver, magma_z_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.
 
magma_int_t magma_z_applyprecond (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective preconditioner is applied.
 
magma_int_t magma_z_applyprecond_left (magma_trans_t trans, magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective left preconditioner is applied.
 
magma_int_t magma_z_applyprecond_right (magma_trans_t trans, magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 For a given input matrix A and vectors x, y and the preconditioner parameters, the respective right-preconditioner is applied.
 
magma_int_t magma_z_solver (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_zopts *zopts, magma_queue_t queue)
 This is an interface that allows to use any iterative solver on the linear system Ax = b.
 
magma_int_t magma_zapplycustomprecond_l (magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 This is an interface to the left solve for any custom preconditioner.
 
magma_int_t magma_zapplycustomprecond_r (magma_z_matrix b, magma_z_matrix *x, magma_z_preconditioner *precond, magma_queue_t queue)
 This is an interface to the right solve for any custom preconditioner.
 
magma_int_t magma_zresidual (magma_z_matrix A, magma_z_matrix b, magma_z_matrix x, double *res, magma_queue_t queue)
 Computes the residual ||b-Ax|| for a solution approximation x.
 
magma_int_t magma_zresidual_slice (magma_int_t start, magma_int_t end, magma_z_matrix A, magma_z_matrix b, magma_z_matrix x, double *res, magma_queue_t queue)
 Computes the residual r=||b-Ax|| for the slice r(start:end) for a solution approximation x.
 
magma_int_t magma_zresidualvec (magma_z_matrix A, magma_z_matrix b, magma_z_matrix x, magma_z_matrix *r, double *res, magma_queue_t queue)
 Computes the residual r = b-Ax for a solution approximation x.
 
magma_int_t magma_zcsrsplit (magma_int_t offset, magma_int_t bsize, magma_z_matrix A, magma_z_matrix *D, magma_z_matrix *R, magma_queue_t queue)
 Splits a CSR matrix into two matrices, one containing the diagonal blocks with the diagonal element stored first, one containing the rest of the original matrix.
 
magma_int_t magma_zmfree (magma_z_matrix *A, magma_queue_t queue)
 Free the memory of a magma_z_matrix.
 
magma_int_t magma_zprecondfree (magma_z_preconditioner *precond_par, magma_queue_t queue)
 Free a preconditioner.
 
magma_int_t magma_zfrobenius (magma_z_matrix A, magma_z_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B.
 
magma_int_t magma_zmatrix_addrowindex (magma_z_matrix *A, magma_queue_t queue)
 Adds to a CSR matrix an array containing the rowindexes.
 
magma_int_t magma_zcsrcoo_transpose (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Transposes a matrix that already contains rowidx.
 
magma_int_t magma_zmatrix_createrowptr (magma_int_t n, magma_index_t *row, magma_queue_t queue)
 This function generates a rowpointer out of a row-wise element count in parallel.
 
magma_int_t magma_zmatrix_swap (magma_z_matrix *A, magma_z_matrix *B, magma_queue_t queue)
 Swaps two matrices.
 
magma_int_t magma_zmatrix_tril (magma_z_matrix A, magma_z_matrix *L, magma_queue_t queue)
 Extracts the lower triangular of a matrix: L = tril(A).
 
magma_int_t magma_zmatrix_triu (magma_z_matrix A, magma_z_matrix *U, magma_queue_t queue)
 Extracts the lower triangular of a matrix: U = triu(A).
 
magma_int_t magma_zcsr_sort (magma_z_matrix *A, magma_queue_t queue)
 SOrts the elements in a CSR matrix for increasing column index.
 
magma_int_t magma_zrowentries (magma_z_matrix *A, magma_queue_t queue)
 Checks the maximal number of nonzeros in a row of matrix A.
 
magma_int_t magma_z_csr_compressor (magmaDoubleComplex **val, magma_index_t **row, magma_index_t **col, magmaDoubleComplex **valn, magma_index_t **rown, magma_index_t **coln, magma_int_t *n, magma_queue_t queue)
 Helper function to compress CSR containing zero-entries.
 
magma_int_t magma_zmconvert (magma_z_matrix A, magma_z_matrix *B, magma_storage_t old_format, magma_storage_t new_format, magma_queue_t queue)
 Converter between different sparse storage formats.
 
magma_int_t magma_zmcsrcompressor (magma_z_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix.
 
magma_int_t magma_zcsrset (magma_int_t m, magma_int_t n, magma_index_t *row, magma_index_t *col, magmaDoubleComplex *val, magma_z_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA.
 
magma_int_t magma_zcsrget (magma_z_matrix A, magma_int_t *m, magma_int_t *n, magma_index_t **row, magma_index_t **col, magmaDoubleComplex **val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure.
 
magma_int_t magma_zcsrset_gpu (magma_int_t m, magma_int_t n, magmaIndex_ptr row, magmaIndex_ptr col, magmaDoubleComplex_ptr val, magma_z_matrix *A, magma_queue_t queue)
 Passes a CSR matrix to MAGMA (located on DEV).
 
magma_int_t magma_zcsrget_gpu (magma_z_matrix A, magma_int_t *m, magma_int_t *n, magmaIndex_ptr *row, magmaIndex_ptr *col, magmaDoubleComplex_ptr *val, magma_queue_t queue)
 Passes a MAGMA matrix to CSR structure (located on DEV).
 
magma_int_t magma_zmdiagdom (magma_z_matrix M, double *min_dd, double *max_dd, double *avg_dd, magma_queue_t queue)
 This routine takes a CSR matrix and computes the average diagonal dominance.
 
magma_int_t magma_zmbdiagdom (magma_z_matrix M, magma_z_matrix blocksizes, double *min_dd, double *max_dd, double *avg_dd, magma_queue_t queue)
 This routine takes a CSR matrix and computes the average block-diagonal dominance.
 
magma_int_t magma_zmdiff (magma_z_matrix A, magma_z_matrix B, real_Double_t *res, magma_queue_t queue)
 Computes the Frobenius norm of the difference between the CSR matrices A and B.
 
magma_int_t magma_zmfrobenius (magma_z_matrix A, magma_z_matrix B, magma_z_matrix S, double *norm, magma_queue_t queue)
 Computes the Frobenius norm || A - B ||_S on the sparsity pattern of S.
 
magma_int_t magma_zmgenerator (magma_int_t n, magma_int_t offdiags, magma_index_t *diag_offset, magmaDoubleComplex *diag_vals, magma_z_matrix *A, magma_queue_t queue)
 Generate a symmetric n x n CSR matrix for a stencil.
 
magma_int_t magma_zm_27stencil (magma_int_t n, magma_z_matrix *A, magma_queue_t queue)
 Generate a 27-point stencil for a 3D FD discretization.
 
magma_int_t magma_zm_5stencil (magma_int_t n, magma_z_matrix *A, magma_queue_t queue)
 Generate a 5-point stencil for a 2D FD discretization.
 
magma_int_t magma_zsymbilu (magma_z_matrix *A, magma_int_t levels, magma_z_matrix *L, magma_z_matrix *U, magma_queue_t queue)
 This routine performs a symbolic ILU factorization.
 
magma_int_t read_z_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, magmaDoubleComplex **val, magma_index_t **row, magma_index_t **col, const char *filename, magma_queue_t queue)
 Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
 
magma_int_t magma_zwrite_csr_mtx (magma_z_matrix A, magma_order_t MajorType, const char *filename, magma_queue_t queue)
 Writes a CSR matrix to a file using Matrix Market format.
 
magma_int_t magma_zprint_csr_mtx (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, magmaDoubleComplex **val, magma_index_t **row, magma_index_t **col, magma_order_t MajorType, magma_queue_t queue)
 Prints a CSR matrix in Matrix Market format.
 
magma_int_t magma_zprint_csr (magma_int_t n_row, magma_int_t n_col, magma_int_t nnz, magmaDoubleComplex **val, magma_index_t **row, magma_index_t **col, magma_queue_t queue)
 Prints a CSR matrix in CSR format.
 
magma_int_t magma_zprint_matrix (magma_z_matrix A, magma_queue_t queue)
 Prints a sparse matrix in CSR format.
 
magma_int_t magma_z_csr_mtx (magma_z_matrix *A, const char *filename, magma_queue_t queue)
 Reads in a matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
 
magma_int_t magma_z_csr_mtxsymm (magma_z_matrix *A, const char *filename, magma_queue_t queue)
 Reads in a SYMMETRIC matrix stored in coo format from a Matrix Market (.mtx) file and converts it into CSR format.
 
magma_int_t magma_zmlumerge (magma_z_matrix L, magma_z_matrix U, magma_z_matrix *A, magma_queue_t queue)
 Takes an strictly lower triangular matrix L and an upper triangular matrix U and merges them into a matrix A containing the upper and lower triangular parts.
 
magma_int_t magma_zmscale (magma_z_matrix *A, magma_scale_t scaling, magma_queue_t queue)
 Scales a matrix.
 
magma_int_t magma_zmscale_matrix_rhs (magma_z_matrix *A, magma_z_matrix *b, magma_z_matrix *scaling_factors, magma_scale_t scaling, magma_queue_t queue)
 Scales a matrix and a right hand side vector of a Ax = b system.
 
magma_int_t magma_zmdiagadd (magma_z_matrix *A, magmaDoubleComplex add, magma_queue_t queue)
 Adds a multiple of the Identity matrix to a matrix: A = A+add * I.
 
magma_int_t magma_zmscale_generate (magma_int_t n, magma_scale_t *scaling, magma_side_t *side, magma_z_matrix *A, magma_z_matrix *scaling_factors, magma_queue_t queue)
 Generates n vectors of scaling factors from the A matrix and stores them in the factors matrix as column vectors in column major ordering.
 
magma_int_t magma_zmscale_apply (magma_int_t n, magma_side_t *side, magma_z_matrix *scaling_factors, magma_z_matrix *A, magma_queue_t queue)
 Applies n diagonal scaling matrices to a matrix A; n=[1,2], factor[i] is applied to side[i] of the matrix.
 
magma_int_t magma_zdimv (magma_z_matrix *vecA, magma_z_matrix *vecB, magma_queue_t queue)
 Multiplies a diagonal matrix (vecA) and a vector (vecB).
 
magma_int_t magma_zmshrink (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Shrinks a non-square matrix (m < n) to the smaller dimension.
 
magma_int_t magma_zmslice (magma_int_t num_slices, magma_int_t slice, magma_z_matrix A, magma_z_matrix *B, magma_z_matrix *ALOC, magma_z_matrix *ANLOC, magma_index_t *comm_i, magmaDoubleComplex *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:
 
magma_int_t magma_zmtransfer (magma_z_matrix A, magma_z_matrix *B, magma_location_t src, magma_location_t dst, magma_queue_t queue)
 Copies a matrix from memory location src to memory location dst.
 
magma_int_t z_transpose_csr (magma_int_t n_rows, magma_int_t n_cols, magma_int_t nnz, magmaDoubleComplex *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, magmaDoubleComplex **new_values, magma_index_t **new_rowptr, magma_index_t **new_colind, magma_queue_t queue)
 Transposes a matrix stored in CSR format on the CPU host.
 
magma_int_t magma_zmtranspose (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Interface to cuSPARSE transpose.
 
magma_int_t magma_z_cucsrtranspose (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Helper function to transpose CSR matrix.
 
magma_int_t magma_zmtransposeconjugate (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 This function forms the transpose conjugate of a matrix.
 
magma_int_t magma_zmtranspose_cpu (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Generates a transpose of A on the CPU.
 
magma_int_t magma_zmtransposeconj_cpu (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Generates a transpose conjugate of A on the CPU.
 
magma_int_t magma_zmtransposestruct_cpu (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Generates a transpose of the nonzero pattern of A on the CPU.
 
magma_int_t magma_zmtransposeabs_cpu (magma_z_matrix A, magma_z_matrix *B, magma_queue_t queue)
 Generates a transpose with absolute values of A on the CPU.
 
magma_int_t magma_zsolverinfo (magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Prints information about a previously called solver.
 
magma_int_t magma_zsolverinfo_free (magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Frees any memory assocoiated with the verbose mode of solver_par.
 
magma_int_t magma_zsolverinfo_init (magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Initializes all solver and preconditioner parameters.
 
magma_int_t magma_zeigensolverinfo_init (magma_z_solver_par *solver_par, magma_queue_t queue)
 Initializes space for eigensolvers.
 
magma_int_t magma_zsort (magmaDoubleComplex *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of values in increasing order.
 
magma_int_t magma_zmsort (magmaDoubleComplex *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.
 
magma_int_t magma_zindexsort (magma_index_t *x, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers in increasing order.
 
magma_int_t magma_zindexsortval (magma_index_t *x, magmaDoubleComplex *y, magma_int_t first, magma_int_t last, magma_queue_t queue)
 Sorts an array of integers, updates a respective array of values.
 
magma_int_t magma_zorderstatistics (magmaDoubleComplex *val, magma_int_t length, magma_int_t k, magma_int_t r, magmaDoubleComplex *element, magma_queue_t queue)
 Identifies the kth smallest/largest element in an array.
 
magma_int_t magma_zbitonic_sort (magma_int_t start, magma_int_t length, magmaDoubleComplex *seq, magma_int_t flag, magma_queue_t queue)
 Approximates the k-th smallest element in an array by using order-statistics with step-size inc.
 
magma_int_t magma_zparse_opts (int argc, char **argv, magma_zopts *opts, int *matrices, magma_queue_t queue)
 Parses input options for a solver.
 
magma_int_t magma_zvinit (magma_z_matrix *x, magma_location_t mem_loc, magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex values, magma_queue_t queue)
 Allocates memory for magma_z_matrix and initializes it with the passed value.
 
magma_int_t magma_zvinit_rand (magma_z_matrix *x, magma_location_t mem_loc, magma_int_t num_rows, magma_int_t num_cols, magma_queue_t queue)
 Allocates memory for magma_z_matrix and initializes it with random values.
 
magma_int_t magma_zprint_vector (magma_z_matrix x, magma_int_t offset, magma_int_t visulen, magma_queue_t queue)
 Visualizes part of a vector of type magma_z_matrix.
 
magma_int_t magma_zvread (magma_z_matrix *x, magma_int_t length, char *filename, magma_queue_t queue)
 Reads in a double vector of length "length".
 
magma_int_t magma_zwrite_vector (magma_z_matrix A, const char *filename, magma_queue_t queue)
 Writes a vector to a file.
 
magma_int_t magma_zvset (magma_int_t m, magma_int_t n, magmaDoubleComplex *val, magma_z_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA.
 
magma_int_t magma_zvcopy (magma_z_matrix v, magma_int_t *m, magma_int_t *n, magmaDoubleComplex *val, magma_queue_t queue)
 Passes a MAGMA vector back.
 
magma_int_t magma_zvset_dev (magma_int_t m, magma_int_t n, magmaDoubleComplex_ptr val, magma_z_matrix *v, magma_queue_t queue)
 Passes a vector to MAGMA (located on DEV).
 
magma_int_t magma_zvget (magma_z_matrix v, magma_int_t *m, magma_int_t *n, magmaDoubleComplex **val, magma_queue_t queue)
 Passes a MAGMA vector back.
 
magma_int_t magma_zvget_dev (magma_z_matrix v, magma_int_t *m, magma_int_t *n, magmaDoubleComplex_ptr *val, magma_queue_t queue)
 Passes a MAGMA vector back (located on DEV).
 
magma_int_t magma_zvcopy_dev (magma_z_matrix v, magma_int_t *m, magma_int_t *n, magmaDoubleComplex_ptr val, magma_queue_t queue)
 Passes a MAGMA vector back (located on DEV).
 
magma_int_t magma_zvtranspose (magma_z_matrix x, magma_z_matrix *y, magma_queue_t queue)
 Transposes a vector from col to row major and vice versa.
 
magma_int_t magma_zdiagcheck (magma_z_matrix dA, magma_queue_t queue)
 This routine checks for a CSR matrix whether there exists a zero on the diagonal.
 
magma_int_t magma_zcsr_sort_gpu (magma_z_matrix *A, magma_queue_t queue)
 Generates a matrix \(U = A \cup B\).
 
magma_int_t magma_zmconjugate (magma_z_matrix *A, magma_queue_t queue)
 This function conjugates a matrix.
 
magma_int_t magma_zmcsrcompressor_gpu (magma_z_matrix *A, magma_queue_t queue)
 Removes zeros in a CSR matrix.
 
magma_int_t magma_zlobpcg_res (magma_int_t num_rows, magma_int_t num_vecs, magmaDouble_ptr evalues, magmaDoubleComplex_ptr X, magmaDoubleComplex_ptr R, magmaDouble_ptr res, magma_queue_t queue)
 This routine computes for Block-LOBPCG, the set of residuals.
 
magma_int_t magma_zlobpcg_shift (magma_int_t num_rows, magma_int_t num_vecs, magma_int_t shift, magmaDoubleComplex_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.
 

Detailed Description

Function Documentation

◆ magma_z_precond()

magma_int_t magma_z_precond ( magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix * x,
magma_z_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_z_matrix sparse matrix A
[in]bmagma_z_matrix input vector b
[in]xmagma_z_matrix* output vector x
[in,out]precondmagma_z_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_precondsetup()

magma_int_t magma_z_precondsetup ( magma_z_matrix A,
magma_z_matrix b,
magma_z_solver_par * solver,
magma_z_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_z_matrix sparse matrix M
[in]bmagma_z_matrix input vector y
[in]solvermagma_z_solver_par solver structure using the preconditioner
[in,out]precondmagma_z_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_applyprecond()

magma_int_t magma_z_applyprecond ( magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix * x,
magma_z_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_z_matrix sparse matrix A
[in]bmagma_z_matrix input vector b
[in,out]xmagma_z_matrix* output vector x
[in]precondmagma_z_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_applyprecond_left()

magma_int_t magma_z_applyprecond_left ( magma_trans_t trans,
magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix * x,
magma_z_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]transmagma_trans_t mode of the preconditioner: MagmaTrans or MagmaNoTrans
[in]Amagma_z_matrix sparse matrix A
[in]bmagma_z_matrix input vector b
[in,out]xmagma_z_matrix* output vector x
[in]precondmagma_z_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_applyprecond_right()

magma_int_t magma_z_applyprecond_right ( magma_trans_t trans,
magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix * x,
magma_z_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]transmagma_trans_t mode of the preconditioner: MagmaTrans or MagmaNoTrans
[in]Amagma_z_matrix sparse matrix A
[in]bmagma_z_matrix input vector b
[in,out]xmagma_z_matrix* output vector x
[in]precondmagma_z_preconditioner preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_solver()

magma_int_t magma_z_solver ( magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix * x,
magma_zopts * zopts,
magma_queue_t queue )

This is an interface that allows to use any iterative solver on the linear system Ax = b.

All linear algebra objects are expected to be on the device, the linear algebra objects are MAGMA-sparse specific structures (dense matrix b, dense matrix x, sparse/dense matrix A). The additional parameter zopts contains information about the solver and the preconditioner. the type of solver the relative / absolute stopping criterion the maximum number of iterations the preconditioner type ... Please see magmasparse_types.h for details about the fields and magma_zutil_sparse.cpp for the possible options.

Parameters
[in]Amagma_z_matrix sparse matrix A
[in]bmagma_z_matrix input vector b
[in]xmagma_z_matrix* output vector x
[in]zoptsmagma_zopts options for solver and preconditioner
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zapplycustomprecond_l()

magma_int_t magma_zapplycustomprecond_l ( magma_z_matrix b,
magma_z_matrix * x,
magma_z_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_z_matrix RHS
[in,out]xmagma_z_matrix* vector to precondition
[in,out]precondmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zapplycustomprecond_r()

magma_int_t magma_zapplycustomprecond_r ( magma_z_matrix b,
magma_z_matrix * x,
magma_z_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_z_matrix RHS
[in,out]xmagma_z_matrix* vector to precondition
[in,out]precondmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zresidual()

magma_int_t magma_zresidual ( magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix x,
double * res,
magma_queue_t queue )

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in]xmagma_z_matrix solution approximation
[out]resmagmaDoubleComplex* return residual
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zresidual_slice()

magma_int_t magma_zresidual_slice ( magma_int_t start,
magma_int_t end,
magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix x,
double * 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_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in]xmagma_z_matrix solution approximation
[out]resmagmaDoubleComplex* return residual
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zresidualvec()

magma_int_t magma_zresidualvec ( magma_z_matrix A,
magma_z_matrix b,
magma_z_matrix x,
magma_z_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_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in]xmagma_z_matrix solution approximation
[in,out]rmagma_z_matrix* residual vector
[out]resmagmaDoubleComplex* return residual
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrsplit()

magma_int_t magma_zcsrsplit ( magma_int_t offset,
magma_int_t bsize,
magma_z_matrix A,
magma_z_matrix * D,
magma_z_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_z_matrix CSR input matrix
[out]Dmagma_z_matrix* CSR matrix containing diagonal blocks
[out]Rmagma_z_matrix* CSR matrix containing rest
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmfree()

magma_int_t magma_zmfree ( magma_z_matrix * A,
magma_queue_t queue )

Free the memory of a magma_z_matrix.

Note, this routine performs a magma_queue_sync on the queue passed to it prior to freeing any memory.

Parameters
[in,out]Amagma_z_matrix* matrix to free
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zprecondfree()

magma_int_t magma_zprecondfree ( magma_z_preconditioner * precond_par,
magma_queue_t queue )

Free a preconditioner.

Note, this routine performs a magma_queue_sync on the queue passed to it prior to freeing any memory.

Parameters
[in,out]precond_parmagma_z_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zfrobenius()

magma_int_t magma_zfrobenius ( magma_z_matrix A,
magma_z_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_z_matrix sparse matrix in CSR
[in]Bmagma_z_matrix sparse matrix in CSR
[out]resreal_Double_t* Frobenius norm of difference
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmatrix_addrowindex()

magma_int_t magma_zmatrix_addrowindex ( magma_z_matrix * A,
magma_queue_t queue )

Adds to a CSR matrix an array containing the rowindexes.

Parameters
[in,out]Amagma_z_matrix* Matrix where rowindexes should be added.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrcoo_transpose()

magma_int_t magma_zcsrcoo_transpose ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Transposes a matrix that already contains rowidx.

The idea is to use a linked list.

Parameters
[in]Amagma_z_matrix Matrix to transpose.
[out]Bmagma_z_matrix* Transposed matrix.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmatrix_createrowptr()

magma_int_t magma_zmatrix_createrowptr ( magma_int_t n,
magma_index_t * row,
magma_queue_t queue )

This function generates a rowpointer out of a row-wise element count in parallel.

Parameters
[in]nmagma_indnt_t row-count.
[in,out]rowmagma_index_t* Input: Vector of size n+1 containing the row-counts (offset by one). Output: Rowpointer.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmatrix_swap()

magma_int_t magma_zmatrix_swap ( magma_z_matrix * A,
magma_z_matrix * B,
magma_queue_t queue )

Swaps two matrices.

Useful if a loop modifies the name of a matrix.

Parameters
[in,out]Amagma_z_matrix* Matrix to be swapped with B.
[in,out]Bmagma_z_matrix* Matrix to be swapped with A.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmatrix_tril()

magma_int_t magma_zmatrix_tril ( magma_z_matrix A,
magma_z_matrix * L,
magma_queue_t queue )

Extracts the lower triangular of a matrix: L = tril(A).

The values of A are preserved.

Parameters
[in]Amagma_z_matrix Element part of this.
[out]Lmagma_z_matrix* Lower triangular part of A.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmatrix_triu()

magma_int_t magma_zmatrix_triu ( magma_z_matrix A,
magma_z_matrix * U,
magma_queue_t queue )

Extracts the lower triangular of a matrix: U = triu(A).

The values of A are preserved.

Parameters
[in]Amagma_z_matrix Element part of this.
[out]Umagma_z_matrix* Lower triangular part of A.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsr_sort()

magma_int_t magma_zcsr_sort ( magma_z_matrix * A,
magma_queue_t queue )

SOrts the elements in a CSR matrix for increasing column index.

Parameters
[in,out]Amagma_z_matrix* CSR matrix, sorted on output.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zrowentries()

magma_int_t magma_zrowentries ( magma_z_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_z_matrix* sparse matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_csr_compressor()

magma_int_t magma_z_csr_compressor ( magmaDoubleComplex ** val,
magma_index_t ** row,
magma_index_t ** col,
magmaDoubleComplex ** 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]valmagmaDoubleComplex** input val pointer to compress
[in]rowmagma_int_t** input row pointer to modify
[in]colmagma_int_t** input col pointer to compress
[in]valnmagmaDoubleComplex** 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_zmconvert()

magma_int_t magma_zmconvert ( magma_z_matrix A,
magma_z_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_z_matrix sparse matrix A
[out]Bmagma_z_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_zmcsrcompressor()

magma_int_t magma_zmcsrcompressor ( magma_z_matrix * A,
magma_queue_t queue )

Removes zeros in a CSR matrix.

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrset()

magma_int_t magma_zcsrset ( magma_int_t m,
magma_int_t n,
magma_index_t * row,
magma_index_t * col,
magmaDoubleComplex * val,
magma_z_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]valmagmaDoubleComplex* array containing matrix entries
[out]Amagma_z_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrget()

magma_int_t magma_zcsrget ( magma_z_matrix A,
magma_int_t * m,
magma_int_t * n,
magma_index_t ** row,
magma_index_t ** col,
magmaDoubleComplex ** val,
magma_queue_t queue )

Passes a MAGMA matrix to CSR structure.

Parameters
[in]Amagma_z_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]valmagmaDoubleComplex* array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrset_gpu()

magma_int_t magma_zcsrset_gpu ( magma_int_t m,
magma_int_t n,
magmaIndex_ptr row,
magmaIndex_ptr col,
magmaDoubleComplex_ptr val,
magma_z_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]valmagmaDoubleComplex_ptr array containing matrix entries
[out]Amagma_z_matrix* matrix in magma sparse matrix format
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsrget_gpu()

magma_int_t magma_zcsrget_gpu ( magma_z_matrix A,
magma_int_t * m,
magma_int_t * n,
magmaIndex_ptr * row,
magmaIndex_ptr * col,
magmaDoubleComplex_ptr * val,
magma_queue_t queue )

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

Parameters
[in]Amagma_z_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]valmagmaDoubleComplex_ptr array containing matrix entries
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmdiagdom()

magma_int_t magma_zmdiagdom ( magma_z_matrix M,
double * min_dd,
double * max_dd,
double * avg_dd,
magma_queue_t queue )

This routine takes a CSR matrix and computes the average diagonal dominance.

For each row i, it computes the abs(d_ii)/sum_j(abs(a_ij)). It returns max, min, and average.

Parameters
[in]Mmagma_z_matrix System matrix.
[out]min_dddouble Smallest diagonal dominance.
[out]max_dddouble Largest diagonal dominance.
[out]avg_dddouble Average diagonal dominance.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmbdiagdom()

magma_int_t magma_zmbdiagdom ( magma_z_matrix M,
magma_z_matrix blocksizes,
double * min_dd,
double * max_dd,
double * avg_dd,
magma_queue_t queue )

This routine takes a CSR matrix and computes the average block-diagonal dominance.

For each row i, it computes the abs( D_(i,:) ) / abs( A(i,:) \ D_(i,:) ). It returns max, min, and average. The input vector bsz contains the blocksizes.

Parameters
[in]Mmagma_z_matrix System matrix.
[in]blocksizesmagma_z_matrix Vector containing blocksizes (as DoubleComplex).
[out]min_dddouble Smallest diagonal dominance.
[out]max_dddouble Largest diagonal dominance.
[out]avg_dddouble Average diagonal dominance.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmdiff()

magma_int_t magma_zmdiff ( magma_z_matrix A,
magma_z_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_z_matrix sparse matrix in CSR
[in]Bmagma_z_matrix sparse matrix in CSR
[out]resreal_Double_t* residual
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmfrobenius()

magma_int_t magma_zmfrobenius ( magma_z_matrix A,
magma_z_matrix B,
magma_z_matrix S,
double * norm,
magma_queue_t queue )

Computes the Frobenius norm || A - B ||_S on the sparsity pattern of S.

Parameters
[in]Amagma_z_matrix input sparse matrix in CSR
[in]Bmagma_z_matrix input sparse matrix in CSR
[in]Smagma_z_matrix input sparsity pattern in CSR
[out]normdouble* Frobenius norm of difference on sparsity pattern S
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmgenerator()

magma_int_t magma_zmgenerator ( magma_int_t n,
magma_int_t offdiags,
magma_index_t * diag_offset,
magmaDoubleComplex * diag_vals,
magma_z_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_valsmagmaDoubleComplex* array containing the values
                            (length offsets+1)
[out]Amagma_z_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zm_27stencil()

magma_int_t magma_zm_27stencil ( magma_int_t n,
magma_z_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_z_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zm_5stencil()

magma_int_t magma_zm_5stencil ( magma_int_t n,
magma_z_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_z_matrix* matrix to generate
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zsymbilu()

magma_int_t magma_zsymbilu ( magma_z_matrix * A,
magma_int_t levels,
magma_z_matrix * L,
magma_z_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_z_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_z_matrix* output lower triangular matrix in magma sparse matrix format empty on function call
[out]Umagma_z_matrix* output upper triangular matrix in magma sparse matrix format empty on function call
[in]queuemagma_queue_t Queue to execute in.

◆ read_z_csr_from_mtx()

magma_int_t read_z_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,
magmaDoubleComplex ** 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]valmagmaDoubleComplex** 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.

◆ magma_zwrite_csr_mtx()

magma_int_t magma_zwrite_csr_mtx ( magma_z_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_z_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_zprint_csr_mtx()

magma_int_t magma_zprint_csr_mtx ( magma_int_t n_row,
magma_int_t n_col,
magma_int_t nnz,
magmaDoubleComplex ** 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]valmagmaDoubleComplex** 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_zprint_csr()

magma_int_t magma_zprint_csr ( magma_int_t n_row,
magma_int_t n_col,
magma_int_t nnz,
magmaDoubleComplex ** 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]valmagmaDoubleComplex** 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_zprint_matrix()

magma_int_t magma_zprint_matrix ( magma_z_matrix A,
magma_queue_t queue )

Prints a sparse matrix in CSR format.

Parameters
[in]Amagma_z_matrix sparse matrix in Magma_CSR format
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_csr_mtx()

magma_int_t magma_z_csr_mtx ( magma_z_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_z_matrix* matrix in magma sparse matrix format
[in]filenameconst char* filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_csr_mtxsymm()

magma_int_t magma_z_csr_mtxsymm ( magma_z_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_z_matrix* matrix in magma sparse matrix format
[in]filenameconst char* filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmlumerge()

magma_int_t magma_zmlumerge ( magma_z_matrix L,
magma_z_matrix U,
magma_z_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_z_matrix input strictly lower triangular matrix L
[in]Umagma_z_matrix input upper triangular matrix U
[out]Amagma_z_matrix* output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmscale()

magma_int_t magma_zmscale ( magma_z_matrix * A,
magma_scale_t scaling,
magma_queue_t queue )

Scales a matrix.

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in]scalingmagma_scale_t scaling type (unit rownorm / unit diagonal)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmscale_matrix_rhs()

magma_int_t magma_zmscale_matrix_rhs ( magma_z_matrix * A,
magma_z_matrix * b,
magma_z_matrix * scaling_factors,
magma_scale_t scaling,
magma_queue_t queue )

Scales a matrix and a right hand side vector of a Ax = b system.

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in,out]bmagma_z_matrix* input/output right hand side vector
[out]scaling_factorsmagma_z_matrix* output scaling factors vector
[in]scalingmagma_scale_t scaling type (unit rownorm / unit diagonal)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmdiagadd()

magma_int_t magma_zmdiagadd ( magma_z_matrix * A,
magmaDoubleComplex add,
magma_queue_t queue )

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

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in]addmagmaDoubleComplex scaling for the identity matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmscale_generate()

magma_int_t magma_zmscale_generate ( magma_int_t n,
magma_scale_t * scaling,
magma_side_t * side,
magma_z_matrix * A,
magma_z_matrix * scaling_factors,
magma_queue_t queue )

Generates n vectors of scaling factors from the A matrix and stores them in the factors matrix as column vectors in column major ordering.

Parameters
[in]nmagma_int_t number of diagonal scaling matrices
[in]scalingmagma_scale_t* array of scaling specifiers
[in]sidemagma_side_t* array of side specifiers
[in]Amagma_z_matrix* input matrix
[out]scaling_factorsmagma_z_matrix* array of diagonal matrices
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmscale_apply()

magma_int_t magma_zmscale_apply ( magma_int_t n,
magma_side_t * side,
magma_z_matrix * scaling_factors,
magma_z_matrix * A,
magma_queue_t queue )

Applies n diagonal scaling matrices to a matrix A; n=[1,2], factor[i] is applied to side[i] of the matrix.

Parameters
[in]nmagma_int_t number of diagonal scaling matrices
[in]sidemagma_side_t* array of side specifiers
[in]scaling_factorsmagma_z_matrix* array of diagonal matrices
[in,out]Amagma_z_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zdimv()

magma_int_t magma_zdimv ( magma_z_matrix * vecA,
magma_z_matrix * vecB,
magma_queue_t queue )

Multiplies a diagonal matrix (vecA) and a vector (vecB).

Parameters
[in]vecAmagma_z_matrix* input matrix
[in,out]vecBmagma_z_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmshrink()

magma_int_t magma_zmshrink ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

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

Parameters
[in]Amagma_z_matrix sparse matrix A
[out]Bmagma_z_matrix* sparse matrix A
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmslice()

magma_int_t magma_zmslice ( magma_int_t num_slices,
magma_int_t slice,
magma_z_matrix A,
magma_z_matrix * B,
magma_z_matrix * ALOC,
magma_z_matrix * ANLOC,
magma_index_t * comm_i,
magmaDoubleComplex * 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_z_matrix sparse matrix in CSR
[out]Bmagma_z_matrix* sparse matrix in CSR
[out]ALOCmagma_z_matrix* sparse matrix in CSR
[out]ANLOCmagma_z_matrix* sparse matrix in CSR
[in,out]comm_imagma_int_t* communication plan
[in,out]comm_vmagmaDoubleComplex* 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_zmtransfer()

magma_int_t magma_zmtransfer ( magma_z_matrix A,
magma_z_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_z_matrix sparse matrix A
[out]Bmagma_z_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.

◆ z_transpose_csr()

magma_int_t z_transpose_csr ( magma_int_t n_rows,
magma_int_t n_cols,
magma_int_t nnz,
magmaDoubleComplex * 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,
magmaDoubleComplex ** 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]valuesmagmaDoubleComplex* 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_valuesmagmaDoubleComplex** 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_zmtranspose()

magma_int_t magma_zmtranspose ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Interface to cuSPARSE transpose.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_z_cucsrtranspose()

magma_int_t magma_z_cucsrtranspose ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Helper function to transpose CSR matrix.

Using the CUSPARSE CSR2CSC function.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmtransposeconjugate()

magma_int_t magma_zmtransposeconjugate ( magma_z_matrix A,
magma_z_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_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmtranspose_cpu()

magma_int_t magma_zmtranspose_cpu ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Generates a transpose of A on the CPU.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmtransposeconj_cpu()

magma_int_t magma_zmtransposeconj_cpu ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Generates a transpose conjugate of A on the CPU.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmtransposestruct_cpu()

magma_int_t magma_zmtransposestruct_cpu ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Generates a transpose of the nonzero pattern of A on the CPU.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmtransposeabs_cpu()

magma_int_t magma_zmtransposeabs_cpu ( magma_z_matrix A,
magma_z_matrix * B,
magma_queue_t queue )

Generates a transpose with absolute values of A on the CPU.

Parameters
[in]Amagma_z_matrix input matrix (CSR)
[out]Bmagma_z_matrix* output matrix (CSR)
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zsolverinfo()

magma_int_t magma_zsolverinfo ( magma_z_solver_par * solver_par,
magma_z_preconditioner * precond_par,
magma_queue_t queue )

Prints information about a previously called solver.

Parameters
[in]solver_parmagma_z_solver_par* structure containing all solver information
[in,out]precond_parmagma_z_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zsolverinfo_free()

magma_int_t magma_zsolverinfo_free ( magma_z_solver_par * solver_par,
magma_z_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_z_solver_par* structure containing all solver information
[in,out]precond_parmagma_z_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zsolverinfo_init()

magma_int_t magma_zsolverinfo_init ( magma_z_solver_par * solver_par,
magma_z_preconditioner * precond_par,
magma_queue_t queue )

Initializes all solver and preconditioner parameters.

Parameters
[in,out]solver_parmagma_z_solver_par* structure containing all solver information
[in,out]precond_parmagma_z_preconditioner* structure containing all preconditioner information
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zeigensolverinfo_init()

magma_int_t magma_zeigensolverinfo_init ( magma_z_solver_par * solver_par,
magma_queue_t queue )

Initializes space for eigensolvers.

Parameters
[in,out]solver_parmagma_z_solver_par* structure containing all solver information
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zsort()

magma_int_t magma_zsort ( magmaDoubleComplex * x,
magma_int_t first,
magma_int_t last,
magma_queue_t queue )

Sorts an array of values in increasing order.

Parameters
[in,out]xmagmaDoubleComplex* 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_zmsort()

magma_int_t magma_zmsort ( magmaDoubleComplex * 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]xmagmaDoubleComplex* 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_zindexsort()

magma_int_t magma_zindexsort ( 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_zindexsortval()

magma_int_t magma_zindexsortval ( magma_index_t * x,
magmaDoubleComplex * 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]ymagmaDoubleComplex* 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_zorderstatistics()

magma_int_t magma_zorderstatistics ( magmaDoubleComplex * val,
magma_int_t length,
magma_int_t k,
magma_int_t r,
magmaDoubleComplex * element,
magma_queue_t queue )

Identifies the kth smallest/largest element in an array.

Parameters
[in,out]valmagmaDoubleComplex* 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]elementmagmaDoubleComplex* location of the respective element
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zbitonic_sort()

magma_int_t magma_zbitonic_sort ( magma_int_t start,
magma_int_t length,
magmaDoubleComplex * seq,
magma_int_t flag,
magma_queue_t queue )

Approximates the k-th smallest element in an array by using order-statistics with step-size inc.

Parameters
[in]startmagma_int_t Start position of the target array.
[in]lengthmagma_int_t Length of the target array.
[in,out]seqmagmaDoubleComplex* Target array, will be modified during operation.
[in]flagmagma_int_t ???
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zparse_opts()

magma_int_t magma_zparse_opts ( int argc,
char ** argv,
magma_zopts * 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_zopts * magma solver options
[out]matricesint counter how many linear systems to process
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvinit()

magma_int_t magma_zvinit ( magma_z_matrix * x,
magma_location_t mem_loc,
magma_int_t num_rows,
magma_int_t num_cols,
magmaDoubleComplex values,
magma_queue_t queue )

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

Parameters
[out]xmagma_z_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]valuesmagmaDoubleComplex entries in vector
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvinit_rand()

magma_int_t magma_zvinit_rand ( magma_z_matrix * x,
magma_location_t mem_loc,
magma_int_t num_rows,
magma_int_t num_cols,
magma_queue_t queue )

Allocates memory for magma_z_matrix and initializes it with random values.

Parameters
[out]xmagma_z_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]queuemagma_queue_t Queue to execute in.

◆ magma_zprint_vector()

magma_int_t magma_zprint_vector ( magma_z_matrix x,
magma_int_t offset,
magma_int_t visulen,
magma_queue_t queue )

Visualizes part of a vector of type magma_z_matrix.

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

Parameters
[in]xmagma_z_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_zvread()

magma_int_t magma_zvread ( magma_z_matrix * x,
magma_int_t length,
char * filename,
magma_queue_t queue )

Reads in a double vector of length "length".

Parameters
[out]xmagma_z_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_zwrite_vector()

magma_int_t magma_zwrite_vector ( magma_z_matrix A,
const char * filename,
magma_queue_t queue )

Writes a vector to a file.

Parameters
[in]Amagma_z_matrix matrix to write out
[in]filenameconst char* output-filname of the mtx matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvset()

magma_int_t magma_zvset ( magma_int_t m,
magma_int_t n,
magmaDoubleComplex * val,
magma_z_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]valmagmaDoubleComplex* array containing vector entries
[out]vmagma_z_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvcopy()

magma_int_t magma_zvcopy ( magma_z_matrix v,
magma_int_t * m,
magma_int_t * n,
magmaDoubleComplex * val,
magma_queue_t queue )

Passes a MAGMA vector back.

This function requires the array val to be already allocated (of size m x n).

Parameters
[in]vmagma_z_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaDoubleComplex* array of size m x n the vector entries are copied into
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvset_dev()

magma_int_t magma_zvset_dev ( magma_int_t m,
magma_int_t n,
magmaDoubleComplex_ptr val,
magma_z_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]valmagmaDoubleComplex* array containing vector entries
[out]vmagma_z_matrix* magma vector
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvget()

magma_int_t magma_zvget ( magma_z_matrix v,
magma_int_t * m,
magma_int_t * n,
magmaDoubleComplex ** val,
magma_queue_t queue )

Passes a MAGMA vector back.

Parameters
[in]vmagma_z_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaDoubleComplex* array containing vector entries
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvget_dev()

magma_int_t magma_zvget_dev ( magma_z_matrix v,
magma_int_t * m,
magma_int_t * n,
magmaDoubleComplex_ptr * val,
magma_queue_t queue )

Passes a MAGMA vector back (located on DEV).

Parameters
[in]vmagma_z_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaDoubleComplex_ptr array containing vector entries
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvcopy_dev()

magma_int_t magma_zvcopy_dev ( magma_z_matrix v,
magma_int_t * m,
magma_int_t * n,
magmaDoubleComplex_ptr val,
magma_queue_t queue )

Passes a MAGMA vector back (located on DEV).

This function requires the array val to be already allocated (of size m x n).

Parameters
[in]vmagma_z_matrix magma vector
[out]mmagma_int_t number of rows
[out]nmagma_int_t number of columns
[out]valmagmaDoubleComplex* array of size m x n on the device the vector entries are copied into
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zvtranspose()

magma_int_t magma_zvtranspose ( magma_z_matrix x,
magma_z_matrix * y,
magma_queue_t queue )

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

Parameters
[in]xmagma_z_matrix input vector
[out]ymagma_z_matrix* output vector
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zdiagcheck()

magma_int_t magma_zdiagcheck ( magma_z_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_z_matrix matrix in CSR format
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zcsr_sort_gpu()

magma_int_t magma_zcsr_sort_gpu ( magma_z_matrix * A,
magma_queue_t queue )

Generates a matrix \(U = A \cup B\).

If both matrices have a nonzero value in the same location, the value of A is used.

This is the GPU version of the operation.

Parameters
[in]Amagma_z_matrix Input matrix 1.
[in]Bmagma_z_matrix Input matrix 2.
[out]Umagma_z_matrix* \(U = A \cup B\). If both matrices have a nonzero value in the same location, the value of A is used.
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmconjugate()

magma_int_t magma_zmconjugate ( magma_z_matrix * A,
magma_queue_t queue )

This function conjugates a matrix.

For a real matrix, no value is changed.

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zmcsrcompressor_gpu()

magma_int_t magma_zmcsrcompressor_gpu ( magma_z_matrix * A,
magma_queue_t queue )

Removes zeros in a CSR matrix.

This is a GPU implementation of the CSR compressor.

Parameters
[in,out]Amagma_z_matrix* input/output matrix
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zlobpcg_res()

magma_int_t magma_zlobpcg_res ( magma_int_t num_rows,
magma_int_t num_vecs,
magmaDouble_ptr evalues,
magmaDoubleComplex_ptr X,
magmaDoubleComplex_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_zaxpy(m, MAGMA_Z_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]XmagmaDoubleComplex_ptr block of eigenvector approximations
[in]RmagmaDoubleComplex_ptr block of residuals
[in]resmagmaDouble_ptr array of residuals
[in]queuemagma_queue_t Queue to execute in.

◆ magma_zlobpcg_shift()

magma_int_t magma_zlobpcg_shift ( magma_int_t num_rows,
magma_int_t num_vecs,
magma_int_t shift,
magmaDoubleComplex_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]xmagmaDoubleComplex_ptr input/output vector x
[in]queuemagma_queue_t Queue to execute in.