![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_zvalinit_gpu (magma_int_t num_el, magmaDoubleComplex_ptr dval, magma_queue_t queue) |
Initializes a device array with zero. | |
magma_int_t | magma_zindexinit_gpu (magma_int_t num_el, magmaIndex_ptr dind, magma_queue_t queue) |
Initializes a device array with zero. | |
magma_int_t | magma_zbajac_csr (magma_int_t localiters, magma_z_matrix D, magma_z_matrix R, magma_z_matrix b, magma_z_matrix *x, magma_queue_t queue) |
This routine is a block-asynchronous Jacobi iteration performing s local Jacobi-updates within the block. | |
magma_int_t | magma_zbajac_csr_overlap (magma_int_t localiters, magma_int_t matrices, magma_int_t overlap, magma_z_matrix *D, magma_z_matrix *R, magma_z_matrix b, magma_z_matrix *x, magma_queue_t queue) |
This routine is a block-asynchronous Jacobi iteration with directed restricted additive Schwarz overlap (top-down) performing s local Jacobi-updates within the block. | |
magma_int_t | magma_zbcsrblockinfo5 (magma_int_t lustep, magma_int_t num_blocks, magma_int_t c_blocks, magma_int_t size_b, magma_index_t *blockinfo, magmaDoubleComplex_ptr val, magmaDoubleComplex_ptr *AII, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine copies the filled blocks from the original matrix A and initializes the blocks that will later be filled in the factorization process with zeros. | |
magma_int_t | magma_zbcsrvalcpy (magma_int_t size_b, magma_int_t num_blocks, magma_int_t num_zblocks, magmaDoubleComplex_ptr *Aval, magmaDoubleComplex_ptr *Bval, magmaDoubleComplex_ptr *Bval2, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine copies the filled blocks from the original matrix A and initializes the blocks that will later be filled in the factorization process with zeros. | |
magma_int_t | magma_zbcsrluegemm (magma_int_t size_b, magma_int_t num_brows, magma_int_t kblocks, magmaDoubleComplex_ptr *dA, magmaDoubleComplex_ptr *dB, magmaDoubleComplex_ptr *dC, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine updates all blocks in the trailing matrix. | |
magma_int_t | magma_zbcsrlupivloc (magma_int_t size_b, magma_int_t kblocks, magmaDoubleComplex_ptr *dA, magmaInt_ptr ipiv, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine updates all blocks in the trailing matrix. | |
magma_int_t | magma_zbcsrswp (magma_int_t r_blocks, magma_int_t size_b, magmaInt_ptr ipiv, magmaDoubleComplex_ptr x, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine swaps rows in the vector *x according to the pivoting in *ipiv. | |
magma_int_t | magma_zbcsrtrsv (magma_uplo_t uplo, magma_int_t r_blocks, magma_int_t c_blocks, magma_int_t size_b, magmaDoubleComplex_ptr A, magma_index_t *blockinfo, magmaDoubleComplex_ptr x, magma_queue_t queue) |
For a Block-CSR ILU factorization, this routine performs the triangular solves. | |
magma_int_t | magma_zcompact (magma_int_t m, magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magmaDouble_ptr dnorms, double tol, magmaInt_ptr active, magmaInt_ptr cBlock, magma_queue_t queue) |
ZCOMPACT takes a set of n vectors of size m (in dA) and their norms and compacts them into the cBlock size<=n vectors that have norms > tol. | |
magma_int_t | magma_zmprepare_batched_gpu (magma_uplo_t uplotype, magma_trans_t transtype, magma_diag_t diagtype, magma_z_matrix L, magma_z_matrix LC, magma_index_t *sizes, magma_index_t *locations, magmaDoubleComplex *trisystems, magmaDoubleComplex *rhs, magma_queue_t queue) |
This routine prepares the batch of small triangular systems that need to be solved for computing the ISAI preconditioner. | |
magma_int_t | magma_zjacobisetup_vector_gpu (magma_int_t num_rows, magma_z_matrix b, magma_z_matrix d, magma_z_matrix c, magma_z_matrix *x, magma_queue_t queue) |
Prepares the Jacobi Iteration according to x^(k+1) = D^(-1) * b - D^(-1) * (L+U) * x^k x^(k+1) = c - M * x^k. | |
magma_int_t | magma_zlobpcg_maxpy (magma_int_t num_rows, magma_int_t num_vecs, magmaDoubleComplex_ptr X, magmaDoubleComplex_ptr Y, magma_queue_t queue) |
This routine computes a axpy for a mxn matrix: | |
magma_int_t | magma_zbicgstab_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex omega, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr p, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgstab_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr s, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgstab_3 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex omega, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr t, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgstab_4 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex omega, magmaDoubleComplex_ptr y, magmaDoubleComplex_ptr z, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr t, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgmerge_spmv1 (magma_z_matrix A, magmaDoubleComplex_ptr d1, magmaDoubleComplex_ptr d2, magmaDoubleComplex_ptr dp, magmaDoubleComplex_ptr dr, magmaDoubleComplex_ptr dv, magmaDoubleComplex_ptr skp, magma_queue_t queue) |
Merges the first SpmV using CSR with the dot product and the computation of alpha. | |
magma_int_t | magma_zbicgmerge_spmv2 (magma_z_matrix A, magmaDoubleComplex_ptr d1, magmaDoubleComplex_ptr d2, magmaDoubleComplex_ptr ds, magmaDoubleComplex_ptr dt, magmaDoubleComplex_ptr skp, magma_queue_t queue) |
Merges the second SpmV using CSR with the dot product and the computation of omega. | |
magma_int_t | magma_zbicgmerge_xrbeta (magma_int_t n, magmaDoubleComplex_ptr d1, magmaDoubleComplex_ptr d2, magmaDoubleComplex_ptr rr, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr t, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr skp, magma_queue_t queue) |
Merges the second SpmV using CSR with the dot product and the computation of omega. | |
magma_int_t | magma_zbicgmerge1 (magma_int_t n, magmaDoubleComplex_ptr skp, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr p, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgmerge2 (magma_int_t n, magmaDoubleComplex_ptr skp, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr s, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgmerge3 (magma_int_t n, magmaDoubleComplex_ptr skp, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr t, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zbicgmerge4 (magma_int_t type, magmaDoubleComplex_ptr skp, magma_queue_t queue) |
Performs some parameter operations for the BiCGSTAB with scalars on GPU. | |
magma_int_t | magma_zmergeblockkrylov (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex_ptr alpha, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr x, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zcgmerge_spmv1 (magma_z_matrix A, magmaDoubleComplex_ptr d1, magmaDoubleComplex_ptr d2, magmaDoubleComplex_ptr dd, magmaDoubleComplex_ptr dz, magmaDoubleComplex_ptr skp, magma_queue_t queue) |
Merges the first SpmV using different formats with the dot product and the computation of rho. | |
magma_int_t | magma_zcgs_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr q, magmaDoubleComplex_ptr u, magmaDoubleComplex_ptr p, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zcgs_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex_ptr r, magmaDoubleComplex_ptr u, magmaDoubleComplex_ptr p, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zcgs_3 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex_ptr v_hat, magmaDoubleComplex_ptr u, magmaDoubleComplex_ptr q, magmaDoubleComplex_ptr t, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zcgs_4 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex_ptr u_hat, magmaDoubleComplex_ptr t, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zidr_smoothing_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex_ptr drs, magmaDoubleComplex_ptr dr, magmaDoubleComplex_ptr dt, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zidr_smoothing_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex omega, magmaDoubleComplex_ptr dx, magmaDoubleComplex_ptr dxs, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex rho, magmaDoubleComplex psi, magmaDoubleComplex_ptr y, magmaDoubleComplex_ptr z, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr w, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex pde, magmaDoubleComplex rde, magmaDoubleComplex_ptr y, magmaDoubleComplex_ptr z, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr q, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_3 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex_ptr pt, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr y, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_4 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex eta, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr pt, magmaDoubleComplex_ptr d, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_5 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex eta, magmaDoubleComplex pds, magmaDoubleComplex_ptr p, magmaDoubleComplex_ptr pt, magmaDoubleComplex_ptr d, magmaDoubleComplex_ptr s, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_6 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex rho, magmaDoubleComplex psi, magmaDoubleComplex_ptr y, magmaDoubleComplex_ptr z, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr w, magmaDoubleComplex_ptr wt, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_7 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex_ptr pt, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr vt, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zqmr_8 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex rho, magmaDoubleComplex psi, magmaDoubleComplex_ptr vt, magmaDoubleComplex_ptr wt, magmaDoubleComplex_ptr y, magmaDoubleComplex_ptr z, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr w, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_ztfqmr_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex sigma, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr Au, magmaDoubleComplex_ptr u_m, magmaDoubleComplex_ptr pu_m, magmaDoubleComplex_ptr u_mp1, magmaDoubleComplex_ptr w, magmaDoubleComplex_ptr d, magmaDoubleComplex_ptr Ad, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_ztfqmr_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex eta, magmaDoubleComplex_ptr d, magmaDoubleComplex_ptr Ad, magmaDoubleComplex_ptr x, magmaDoubleComplex_ptr r, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_ztfqmr_3 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex_ptr w, magmaDoubleComplex_ptr u_m, magmaDoubleComplex_ptr u_mp1, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_ztfqmr_4 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex beta, magmaDoubleComplex_ptr Au_new, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr Au, magma_queue_t queue) |
Merges multiple operations into one kernel: | |
magma_int_t | magma_ztfqmr_5 (magma_int_t num_rows, magma_int_t num_cols, magmaDoubleComplex alpha, magmaDoubleComplex sigma, magmaDoubleComplex_ptr v, magmaDoubleComplex_ptr Au, magmaDoubleComplex_ptr u_mp1, magmaDoubleComplex_ptr w, magmaDoubleComplex_ptr d, magmaDoubleComplex_ptr Ad, magma_queue_t queue) |
Mergels multiple operations into one kernel: | |
magma_int_t | magma_zparic_csr (magma_z_matrix A, magma_z_matrix A_CSR, magma_queue_t queue) |
This routine iteratively computes an incomplete LU factorization. | |
magma_int_t | magma_zparilu_csr (magma_z_matrix A, magma_z_matrix L, magma_z_matrix U, magma_queue_t queue) |
This routine iteratively computes an incomplete LU factorization. | |
magma_int_t | magma_zparilut_sweep_gpu (magma_z_matrix *A, magma_z_matrix *L, magma_z_matrix *U, magma_queue_t queue) |
This function does an ParILUT sweep. | |
magma_int_t | magma_zparilut_residuals_gpu (magma_z_matrix A, magma_z_matrix L, magma_z_matrix U, magma_z_matrix *R, magma_queue_t queue) |
This function computes the ILU residual in the locations included in the sparsity pattern of R. | |
magma_int_t magma_zvalinit_gpu | ( | magma_int_t | num_el, |
magmaDoubleComplex_ptr | dval, | ||
magma_queue_t | queue ) |
Initializes a device array with zero.
[in] | num_el | magma_int_t size of array |
[in,out] | dval | magmaDoubleComplex_ptr array to initialize |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zindexinit_gpu | ( | magma_int_t | num_el, |
magmaIndex_ptr | dind, | ||
magma_queue_t | queue ) |
Initializes a device array with zero.
[in] | num_el | magma_int_t size of array |
[in,out] | dind | magmaIndex_ptr array to initialize |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbajac_csr | ( | magma_int_t | localiters, |
magma_z_matrix | D, | ||
magma_z_matrix | R, | ||
magma_z_matrix | b, | ||
magma_z_matrix * | x, | ||
magma_queue_t | queue ) |
This routine is a block-asynchronous Jacobi iteration performing s local Jacobi-updates within the block.
Input format is two CSR matrices, one containing the diagonal blocks, one containing the rest.
[in] | localiters | magma_int_t number of local Jacobi-like updates |
[in] | D | magma_z_matrix input matrix with diagonal blocks |
[in] | R | magma_z_matrix input matrix with non-diagonal parts |
[in] | b | magma_z_matrix RHS |
[in] | x | magma_z_matrix* iterate/solution |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbajac_csr_overlap | ( | magma_int_t | localiters, |
magma_int_t | matrices, | ||
magma_int_t | overlap, | ||
magma_z_matrix * | D, | ||
magma_z_matrix * | R, | ||
magma_z_matrix | b, | ||
magma_z_matrix * | x, | ||
magma_queue_t | queue ) |
This routine is a block-asynchronous Jacobi iteration with directed restricted additive Schwarz overlap (top-down) performing s local Jacobi-updates within the block.
Input format is two CSR matrices, one containing the diagonal blocks, one containing the rest.
[in] | localiters | magma_int_t number of local Jacobi-like updates |
[in] | matrices | magma_int_t number of sub-matrices |
[in] | overlap | magma_int_t size of the overlap |
[in] | D | magma_z_matrix* set of matrices with diagonal blocks |
[in] | R | magma_z_matrix* set of matrices with non-diagonal parts |
[in] | b | magma_z_matrix RHS |
[in] | x | magma_z_matrix* iterate/solution |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrblockinfo5 | ( | magma_int_t | lustep, |
magma_int_t | num_blocks, | ||
magma_int_t | c_blocks, | ||
magma_int_t | size_b, | ||
magma_index_t * | blockinfo, | ||
magmaDoubleComplex_ptr | val, | ||
magmaDoubleComplex_ptr * | AII, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine copies the filled blocks from the original matrix A and initializes the blocks that will later be filled in the factorization process with zeros.
[in] | lustep | magma_int_t lustep |
[in] | num_blocks | magma_int_t number of nonzero blocks |
[in] | c_blocks | magma_int_t number of column-blocks |
[in] | size_b | magma_int_t blocksize |
[in] | blockinfo | magma_int_t* block filled? location? |
[in] | val | magmaDoubleComplex* pointers to the nonzero blocks in A |
[in] | AII | magmaDoubleComplex** pointers to the respective nonzero blocks in B |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrvalcpy | ( | magma_int_t | size_b, |
magma_int_t | num_blocks, | ||
magma_int_t | num_zblocks, | ||
magmaDoubleComplex_ptr * | Aval, | ||
magmaDoubleComplex_ptr * | Bval, | ||
magmaDoubleComplex_ptr * | Bval2, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine copies the filled blocks from the original matrix A and initializes the blocks that will later be filled in the factorization process with zeros.
[in] | size_b | magma_int_t blocksize in BCSR |
[in] | num_blocks | magma_int_t number of nonzero blocks |
[in] | num_zblocks | magma_int_t number of zero-blocks (will later be filled) |
[in] | Aval | magmaDoubleComplex_ptr * pointers to the nonzero blocks in A |
[in] | Bval | magmaDoubleComplex_ptr * pointers to the nonzero blocks in B |
[in] | Bval2 | magmaDoubleComplex_ptr * pointers to the zero blocks in B |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrluegemm | ( | magma_int_t | size_b, |
magma_int_t | num_brows, | ||
magma_int_t | kblocks, | ||
magmaDoubleComplex_ptr * | dA, | ||
magmaDoubleComplex_ptr * | dB, | ||
magmaDoubleComplex_ptr * | dC, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine updates all blocks in the trailing matrix.
[in] | size_b | magma_int_t blocksize in BCSR |
[in] | num_brows | magma_int_t number of block rows |
[in] | kblocks | magma_int_t number of blocks in row |
[in] | dA | magmaDoubleComplex** input blocks of matrix A |
[in] | dB | magmaDoubleComplex** input blocks of matrix B |
[in] | dC | magmaDoubleComplex** output blocks of matrix C |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrlupivloc | ( | magma_int_t | size_b, |
magma_int_t | kblocks, | ||
magmaDoubleComplex_ptr * | dA, | ||
magmaInt_ptr | ipiv, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine updates all blocks in the trailing matrix.
[in] | size_b | magma_int_t blocksize in BCSR |
[in] | kblocks | magma_int_t number of blocks |
[in] | dA | magmaDoubleComplex_ptr * matrix in BCSR |
[in] | ipiv | magmaInt_ptr array containing pivots |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrswp | ( | magma_int_t | r_blocks, |
magma_int_t | size_b, | ||
magmaInt_ptr | ipiv, | ||
magmaDoubleComplex_ptr | x, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine swaps rows in the vector *x according to the pivoting in *ipiv.
[in] | r_blocks | magma_int_t number of blocks |
[in] | size_b | magma_int_t blocksize in BCSR |
[in] | ipiv | magma_int_t* array containing pivots |
[in] | x | magmaDoubleComplex_ptr input/output vector x |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbcsrtrsv | ( | magma_uplo_t | uplo, |
magma_int_t | r_blocks, | ||
magma_int_t | c_blocks, | ||
magma_int_t | size_b, | ||
magmaDoubleComplex_ptr | A, | ||
magma_index_t * | blockinfo, | ||
magmaDoubleComplex_ptr | x, | ||
magma_queue_t | queue ) |
For a Block-CSR ILU factorization, this routine performs the triangular solves.
[in] | uplo | magma_uplo_t upper/lower fill structure |
[in] | r_blocks | magma_int_t number of blocks in row |
[in] | c_blocks | magma_int_t number of blocks in column |
[in] | size_b | magma_int_t blocksize in BCSR |
[in] | A | magmaDoubleComplex_ptr upper/lower factor |
[in] | blockinfo | magma_int_t* array containing matrix information |
[in] | x | magmaDoubleComplex_ptr input/output vector x |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcompact | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDoubleComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dnorms, | ||
double | tol, | ||
magmaInt_ptr | active, | ||
magmaInt_ptr | cBlock, | ||
magma_queue_t | queue ) |
ZCOMPACT takes a set of n vectors of size m (in dA) and their norms and compacts them into the cBlock size<=n vectors that have norms > tol.
The active mask array has 1 or 0, showing if a vector remained or not in the compacted resulting set of vectors.
[in] | m | INTEGER The number of rows of the matrix dA. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix dA. N >= 0. |
[in,out] | dA | COMPLEX DOUBLE PRECISION array, dimension (LDDA,N) The m by n matrix dA. |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). |
[in] | dnorms | DOUBLE PRECISION array, dimension N The norms of the N vectors in dA |
[in] | tol | DOUBLE PRECISON The tolerance value used in the criteria to compact or not. |
[in,out] | active | INTEGER array, dimension N A mask of 1s and 0s showing if a vector remains or has been removed |
[in,out] | cBlock | magmaInt_ptr The number of vectors that remain in dA (i.e., with norms > tol). |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zmprepare_batched_gpu | ( | magma_uplo_t | uplotype, |
magma_trans_t | transtype, | ||
magma_diag_t | diagtype, | ||
magma_z_matrix | L, | ||
magma_z_matrix | LC, | ||
magma_index_t * | sizes, | ||
magma_index_t * | locations, | ||
magmaDoubleComplex * | trisystems, | ||
magmaDoubleComplex * | rhs, | ||
magma_queue_t | queue ) |
This routine prepares the batch of small triangular systems that need to be solved for computing the ISAI preconditioner.
[in] | uplotype | magma_uplo_t input matrix |
[in] | transtype | magma_trans_t input matrix |
[in] | diagtype | magma_diag_t input matrix |
[in] | L | magma_z_matrix triangular factor for which the ISAI matrix is computed. Col-Major CSR storage. |
[in] | LC | magma_z_matrix sparsity pattern of the ISAI matrix. Col-Major CSR storage. |
[in,out] | sizes | magma_index_t* array containing the sizes of the small triangular systems |
[in,out] | locations | magma_index_t* array containing the locations in the respective column of L |
[in,out] | trisystems | magmaDoubleComplex* batch of generated small triangular systems. All systems are embedded in uniform memory blocks of size BLOCKSIZE x BLOCKSIZE |
[in,out] | rhs | magmaDoubleComplex* RHS of the small triangular systems |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zjacobisetup_vector_gpu | ( | magma_int_t | num_rows, |
magma_z_matrix | b, | ||
magma_z_matrix | d, | ||
magma_z_matrix | c, | ||
magma_z_matrix * | x, | ||
magma_queue_t | queue ) |
Prepares the Jacobi Iteration according to x^(k+1) = D^(-1) * b - D^(-1) * (L+U) * x^k x^(k+1) = c - M * x^k.
Returns the vector c. It calls a GPU kernel
[in] | num_rows | magma_int_t number of rows |
[in] | b | magma_z_matrix RHS b |
[in] | d | magma_z_matrix vector with diagonal entries |
[out] | c | magma_z_matrix* c = D^(-1) * b |
[out] | x | magma_z_matrix* iteration vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zlobpcg_maxpy | ( | magma_int_t | num_rows, |
magma_int_t | num_vecs, | ||
magmaDoubleComplex_ptr | X, | ||
magmaDoubleComplex_ptr | Y, | ||
magma_queue_t | queue ) |
This routine computes a axpy for a mxn matrix:
Y = X + Y
It replaces: magma_zaxpy(m*n, c_one, Y, 1, X, 1);
/ 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] /
[in] | num_rows | magma_int_t number of rows |
[in] | num_vecs | magma_int_t number of vectors |
[in] | X | magmaDoubleComplex_ptr input vector X |
[in,out] | Y | magmaDoubleComplex_ptr input/output vector Y |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgstab_1 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex | omega, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | p, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
p = r + beta * ( p - omega * v )
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | omega | magmaDoubleComplex scalar |
[in] | r | magmaDoubleComplex_ptr vector |
[in] | v | magmaDoubleComplex_ptr vector |
[in,out] | p | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgstab_2 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | s, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
s = r - alpha v
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | r | magmaDoubleComplex_ptr vector |
[in] | v | magmaDoubleComplex_ptr vector |
[in,out] | s | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgstab_3 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex | omega, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | t, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
x = x + alpha * p + omega * s r = s - omega * t
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | omega | magmaDoubleComplex scalar |
[in] | p | magmaDoubleComplex_ptr vector |
[in] | s | magmaDoubleComplex_ptr vector |
[in] | t | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgstab_4 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex | omega, | ||
magmaDoubleComplex_ptr | y, | ||
magmaDoubleComplex_ptr | z, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | t, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
x = x + alpha * y + omega * z r = s - omega * t
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | omega | magmaDoubleComplex scalar |
[in] | y | magmaDoubleComplex_ptr vector |
[in] | z | magmaDoubleComplex_ptr vector |
[in] | s | magmaDoubleComplex_ptr vector |
[in] | t | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgmerge_spmv1 | ( | magma_z_matrix | A, |
magmaDoubleComplex_ptr | d1, | ||
magmaDoubleComplex_ptr | d2, | ||
magmaDoubleComplex_ptr | dp, | ||
magmaDoubleComplex_ptr | dr, | ||
magmaDoubleComplex_ptr | dv, | ||
magmaDoubleComplex_ptr | skp, | ||
magma_queue_t | queue ) |
Merges the first SpmV using CSR with the dot product and the computation of alpha.
[in] | A | magma_z_matrix system matrix |
[in] | d1 | magmaDoubleComplex_ptr temporary vector |
[in] | d2 | magmaDoubleComplex_ptr temporary vector |
[in] | dp | magmaDoubleComplex_ptr input vector p |
[in] | dr | magmaDoubleComplex_ptr input vector r |
[in] | dv | magmaDoubleComplex_ptr output vector v |
[in,out] | skp | magmaDoubleComplex_ptr array for parameters ( skp[0]=alpha ) |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgmerge_spmv2 | ( | magma_z_matrix | A, |
magmaDoubleComplex_ptr | d1, | ||
magmaDoubleComplex_ptr | d2, | ||
magmaDoubleComplex_ptr | ds, | ||
magmaDoubleComplex_ptr | dt, | ||
magmaDoubleComplex_ptr | skp, | ||
magma_queue_t | queue ) |
Merges the second SpmV using CSR with the dot product and the computation of omega.
[in] | A | magma_z_matrix input matrix |
[in] | d1 | magmaDoubleComplex_ptr temporary vector |
[in] | d2 | magmaDoubleComplex_ptr temporary vector |
[in] | ds | magmaDoubleComplex_ptr input vector s |
[in] | dt | magmaDoubleComplex_ptr output vector t |
[in,out] | skp | magmaDoubleComplex_ptr array for parameters |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgmerge_xrbeta | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | d1, | ||
magmaDoubleComplex_ptr | d2, | ||
magmaDoubleComplex_ptr | rr, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | t, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | skp, | ||
magma_queue_t | queue ) |
Merges the second SpmV using CSR with the dot product and the computation of omega.
[in] | n | int dimension n |
[in] | d1 | magmaDoubleComplex_ptr temporary vector |
[in] | d2 | magmaDoubleComplex_ptr temporary vector |
[in] | rr | magmaDoubleComplex_ptr input vector rr |
[in] | r | magmaDoubleComplex_ptr input/output vector r |
[in] | p | magmaDoubleComplex_ptr input vector p |
[in] | s | magmaDoubleComplex_ptr input vector s |
[in] | t | magmaDoubleComplex_ptr input vector t |
[out] | x | magmaDoubleComplex_ptr output vector x |
[in] | skp | magmaDoubleComplex_ptr array for parameters |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgmerge1 | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | skp, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | p, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
p = beta*p p = p-omega*beta*v p = p+r
-> p = r + beta * ( p - omega * v )
[in] | n | int dimension n |
[in] | skp | magmaDoubleComplex_ptr set of scalar parameters |
[in] | v | magmaDoubleComplex_ptr input vector v |
[in] | r | magmaDoubleComplex_ptr input vector r |
[in,out] | p | magmaDoubleComplex_ptr input/output vector p |
[in] | queue | magma_queue_t queue to execute in. |
magma_int_t magma_zbicgmerge2 | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | skp, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | s, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
s=r s=s-alpha*v
-> s = r - alpha * v
[in] | n | int dimension n |
[in] | skp | magmaDoubleComplex_ptr set of scalar parameters |
[in] | r | magmaDoubleComplex_ptr input vector r |
[in] | v | magmaDoubleComplex_ptr input vector v |
[out] | s | magmaDoubleComplex_ptr output vector s |
[in] | queue | magma_queue_t queue to execute in. |
magma_int_t magma_zbicgmerge3 | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | skp, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | t, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
x=x+alpha*p x=x+omega*s r=s r=r-omega*t
-> x = x + alpha * p + omega * s -> r = s - omega * t
[in] | n | int dimension n |
[in] | skp | magmaDoubleComplex_ptr set of scalar parameters |
[in] | p | magmaDoubleComplex_ptr input p |
[in] | s | magmaDoubleComplex_ptr input s |
[in] | t | magmaDoubleComplex_ptr input t |
[in,out] | x | magmaDoubleComplex_ptr input/output x |
[in,out] | r | magmaDoubleComplex_ptr input/output r |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zbicgmerge4 | ( | magma_int_t | type, |
magmaDoubleComplex_ptr | skp, | ||
magma_queue_t | queue ) |
Performs some parameter operations for the BiCGSTAB with scalars on GPU.
[in] | type | int kernel type |
[in,out] | skp | magmaDoubleComplex_ptr vector with parameters |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zmergeblockkrylov | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex_ptr | alpha, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | x, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
v = y / rho y = y / rho w = wt / psi z = z / psi
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex_ptr matrix containing all SKP |
[in] | p | magmaDoubleComplex_ptr search directions |
[in,out] | x | magmaDoubleComplex_ptr approximation vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcgmerge_spmv1 | ( | magma_z_matrix | A, |
magmaDoubleComplex_ptr | d1, | ||
magmaDoubleComplex_ptr | d2, | ||
magmaDoubleComplex_ptr | dd, | ||
magmaDoubleComplex_ptr | dz, | ||
magmaDoubleComplex_ptr | skp, | ||
magma_queue_t | queue ) |
Merges the first SpmV using different formats with the dot product and the computation of rho.
[in] | A | magma_z_matrix input matrix |
[in] | d1 | magmaDoubleComplex_ptr temporary vector |
[in] | d2 | magmaDoubleComplex_ptr temporary vector |
[in] | dd | magmaDoubleComplex_ptr input vector d |
[out] | dz | magmaDoubleComplex_ptr input vector z |
[out] | skp | magmaDoubleComplex_ptr array for parameters ( skp[3]=rho ) |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcgs_1 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | q, | ||
magmaDoubleComplex_ptr | u, | ||
magmaDoubleComplex_ptr | p, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
u = r + beta q p = u + beta*(q + beta*p)
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | r | magmaDoubleComplex_ptr vector |
[in] | q | magmaDoubleComplex_ptr vector |
[in,out] | u | magmaDoubleComplex_ptr vector |
[in,out] | p | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcgs_2 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex_ptr | r, | ||
magmaDoubleComplex_ptr | u, | ||
magmaDoubleComplex_ptr | p, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
u = r p = r
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | r | magmaDoubleComplex_ptr vector |
[in,out] | u | magmaDoubleComplex_ptr vector |
[in,out] | p | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcgs_3 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex_ptr | v_hat, | ||
magmaDoubleComplex_ptr | u, | ||
magmaDoubleComplex_ptr | q, | ||
magmaDoubleComplex_ptr | t, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
q = u - alpha v_hat t = u + q
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | v_hat | magmaDoubleComplex_ptr vector |
[in] | u | magmaDoubleComplex_ptr vector |
[in,out] | q | magmaDoubleComplex_ptr vector |
[in,out] | t | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zcgs_4 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex_ptr | u_hat, | ||
magmaDoubleComplex_ptr | t, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
x = x + alpha u_hat r = r -alpha*A u_hat = r -alpha*t
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | u_hat | magmaDoubleComplex_ptr vector |
[in] | t | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zidr_smoothing_1 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex_ptr | drs, | ||
magmaDoubleComplex_ptr | dr, | ||
magmaDoubleComplex_ptr | dt, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
dt = drs - dr
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | drs | magmaDoubleComplex_ptr vector |
[in] | dr | magmaDoubleComplex_ptr vector |
[in,out] | dt | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zidr_smoothing_2 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | omega, | ||
magmaDoubleComplex_ptr | dx, | ||
magmaDoubleComplex_ptr | dxs, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
dxs = dxs - gamma*(dxs-dx)
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | omega | magmaDoubleComplex scalar |
[in] | dx | magmaDoubleComplex_ptr vector |
[in,out] | dxs | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_1 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | rho, | ||
magmaDoubleComplex | psi, | ||
magmaDoubleComplex_ptr | y, | ||
magmaDoubleComplex_ptr | z, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | w, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
v = y / rho y = y / rho w = wt / psi z = z / psi
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | rho | magmaDoubleComplex scalar |
[in] | psi | magmaDoubleComplex scalar |
[in,out] | y | magmaDoubleComplex_ptr vector |
[in,out] | z | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | w | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_2 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | pde, | ||
magmaDoubleComplex | rde, | ||
magmaDoubleComplex_ptr | y, | ||
magmaDoubleComplex_ptr | z, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | q, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
p = y - pde * p q = z - rde * q
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | pde | magmaDoubleComplex scalar |
[in] | rde | magmaDoubleComplex scalar |
[in] | y | magmaDoubleComplex_ptr vector |
[in] | z | magmaDoubleComplex_ptr vector |
[in,out] | p | magmaDoubleComplex_ptr vector |
[in,out] | q | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_3 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex_ptr | pt, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | y, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
v = pt - beta * v y = v
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | pt | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | y | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_4 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | eta, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | pt, | ||
magmaDoubleComplex_ptr | d, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
d = eta * p; s = eta * pt; x = x + d; r = r - s;
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | eta | magmaDoubleComplex scalar |
[in] | p | magmaDoubleComplex_ptr vector |
[in] | pt | magmaDoubleComplex_ptr vector |
[in,out] | d | magmaDoubleComplex_ptr vector |
[in,out] | s | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_5 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | eta, | ||
magmaDoubleComplex | pds, | ||
magmaDoubleComplex_ptr | p, | ||
magmaDoubleComplex_ptr | pt, | ||
magmaDoubleComplex_ptr | d, | ||
magmaDoubleComplex_ptr | s, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
d = eta * p + pds * d; s = eta * pt + pds * s; x = x + d; r = r - s;
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | eta | magmaDoubleComplex scalar |
[in] | pds | magmaDoubleComplex scalar |
[in] | p | magmaDoubleComplex_ptr vector |
[in] | pt | magmaDoubleComplex_ptr vector |
[in,out] | d | magmaDoubleComplex_ptr vector |
[in,out] | s | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_6 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex | rho, | ||
magmaDoubleComplex | psi, | ||
magmaDoubleComplex_ptr | y, | ||
magmaDoubleComplex_ptr | z, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | w, | ||
magmaDoubleComplex_ptr | wt, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
wt = wt - conj(beta) * w v = y / rho y = y / rho w = wt / psi z = wt / psi
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | rho | magmaDoubleComplex scalar |
[in] | psi | magmaDoubleComplex scalar |
[in,out] | y | magmaDoubleComplex_ptr vector |
[in,out] | z | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | w | magmaDoubleComplex_ptr vector |
[in,out] | wt | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_7 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex_ptr | pt, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | vt, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
vt = pt - beta * v
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | pt | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | vt | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zqmr_8 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | rho, | ||
magmaDoubleComplex | psi, | ||
magmaDoubleComplex_ptr | vt, | ||
magmaDoubleComplex_ptr | wt, | ||
magmaDoubleComplex_ptr | y, | ||
magmaDoubleComplex_ptr | z, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | w, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
v = y / rho y = y / rho w = wt / psi z = z / psi
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | rho | magmaDoubleComplex scalar |
[in] | psi | magmaDoubleComplex scalar |
[in] | vt | magmaDoubleComplex_ptr vector |
[in] | wt | magmaDoubleComplex_ptr vector |
[in,out] | y | magmaDoubleComplex_ptr vector |
[in,out] | z | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | w | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_ztfqmr_1 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex | sigma, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | Au, | ||
magmaDoubleComplex_ptr | u_m, | ||
magmaDoubleComplex_ptr | pu_m, | ||
magmaDoubleComplex_ptr | u_mp1, | ||
magmaDoubleComplex_ptr | w, | ||
magmaDoubleComplex_ptr | d, | ||
magmaDoubleComplex_ptr | Ad, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
u_mp1 = u_mp1 - alpha*v; w = w - alpha*Au; d = pu_m + sigma*d; Ad = Au + sigma*Ad;
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | sigma | magmaDoubleComplex scalar |
[in] | v | magmaDoubleComplex_ptr vector |
[in] | Au | magmaDoubleComplex_ptr vector |
[in,out] | u_m | magmaDoubleComplex_ptr vector |
[in,out] | pu_m | magmaDoubleComplex_ptr vector |
[in,out] | u_mp1 | magmaDoubleComplex_ptr vector |
[in,out] | w | magmaDoubleComplex_ptr vector |
[in,out] | d | magmaDoubleComplex_ptr vector |
[in,out] | Ad | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_ztfqmr_2 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | eta, | ||
magmaDoubleComplex_ptr | d, | ||
magmaDoubleComplex_ptr | Ad, | ||
magmaDoubleComplex_ptr | x, | ||
magmaDoubleComplex_ptr | r, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
x = x + eta * d r = r - eta * Ad
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | eta | magmaDoubleComplex scalar |
[in] | d | magmaDoubleComplex_ptr vector |
[in] | Ad | magmaDoubleComplex_ptr vector |
[in,out] | x | magmaDoubleComplex_ptr vector |
[in,out] | r | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_ztfqmr_3 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex_ptr | w, | ||
magmaDoubleComplex_ptr | u_m, | ||
magmaDoubleComplex_ptr | u_mp1, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
u_mp1 = w + beta*u_mp1
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | w | magmaDoubleComplex_ptr vector |
[in] | u_m | magmaDoubleComplex_ptr vector |
[in,out] | u_mp1 | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_ztfqmr_4 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | beta, | ||
magmaDoubleComplex_ptr | Au_new, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | Au, | ||
magma_queue_t | queue ) |
Merges multiple operations into one kernel:
v = Au_new + beta*(Au+beta*v); Au = Au_new
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | beta | magmaDoubleComplex scalar |
[in] | Au_new | magmaDoubleComplex_ptr vector |
[in,out] | v | magmaDoubleComplex_ptr vector |
[in,out] | Au | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_ztfqmr_5 | ( | magma_int_t | num_rows, |
magma_int_t | num_cols, | ||
magmaDoubleComplex | alpha, | ||
magmaDoubleComplex | sigma, | ||
magmaDoubleComplex_ptr | v, | ||
magmaDoubleComplex_ptr | Au, | ||
magmaDoubleComplex_ptr | u_mp1, | ||
magmaDoubleComplex_ptr | w, | ||
magmaDoubleComplex_ptr | d, | ||
magmaDoubleComplex_ptr | Ad, | ||
magma_queue_t | queue ) |
Mergels multiple operations into one kernel:
w = w - alpha*Au; d = pu_m + sigma*d; Ad = Au + sigma*Ad;
[in] | num_rows | magma_int_t dimension m |
[in] | num_cols | magma_int_t dimension n |
[in] | alpha | magmaDoubleComplex scalar |
[in] | sigma | magmaDoubleComplex scalar |
[in] | v | magmaDoubleComplex_ptr vector |
[in] | Au | magmaDoubleComplex_ptr vector |
[in,out] | u_mp1 | magmaDoubleComplex_ptr vector |
[in,out] | w | magmaDoubleComplex_ptr vector |
[in,out] | d | magmaDoubleComplex_ptr vector |
[in,out] | Ad | magmaDoubleComplex_ptr vector |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zparic_csr | ( | magma_z_matrix | A, |
magma_z_matrix | A_CSR, | ||
magma_queue_t | queue ) |
This routine iteratively computes an incomplete LU factorization.
For reference, see: E. Chow and A. Patel: "Fine-grained Parallel Incomplete LU Factorization", SIAM Journal on Scientific Computing, 37, C169-C193 (2015). This routine was used in the ISC 2015 paper: E. Chow et al.: "Asynchronous Iterative Algorithm for Computing Incomplete Factorizations on GPUs", ISC High Performance 2015, LNCS 9137, pp. 1-16, 2015.
The input format of the initial guess matrix A is Magma_CSRCOO, A_CSR is CSR or CSRCOO format.
[in] | A | magma_z_matrix input matrix A - initial guess (lower triangular) |
[in,out] | A_CSR | magma_z_matrix input/output matrix containing the IC approximation |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zparilu_csr | ( | magma_z_matrix | A, |
magma_z_matrix | L, | ||
magma_z_matrix | U, | ||
magma_queue_t | queue ) |
This routine iteratively computes an incomplete LU factorization.
For reference, see: E. Chow and A. Patel: "Fine-grained Parallel Incomplete LU Factorization", SIAM Journal on Scientific Computing, 37, C169-C193 (2015). This routine was used in the ISC 2015 paper: E. Chow et al.: "Asynchronous Iterative Algorithm for Computing Incomplete Factorizations on GPUs", ISC High Performance 2015, LNCS 9137, pp. 1-16, 2015.
The input format of the system matrix is COO, the lower triangular factor L is stored in CSR, the upper triangular factor U is transposed, then also stored in CSR (equivalent to CSC format for the non-transposed U). Every component of L and U is handled by one thread.
[in] | A | magma_z_matrix input matrix A determing initial guess & processing order |
[in,out] | L | magma_z_matrix input/output matrix L containing the lower triangular factor |
[in,out] | U | magma_z_matrix input/output matrix U containing the upper triangular factor |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zparilut_sweep_gpu | ( | magma_z_matrix * | A, |
magma_z_matrix * | L, | ||
magma_z_matrix * | U, | ||
magma_queue_t | queue ) |
This function does an ParILUT sweep.
The difference to the ParILU sweep is that the nonzero pattern of A and the incomplete factors L and U can be different. The pattern determing which elements are iterated are hence the pattern of L and U, not A. L has a unit diagonal.
This is the GPU version of the asynchronous ParILUT sweep.
[in] | A | magma_z_matrix* System matrix. The format is sorted CSR. |
[in,out] | L | magma_z_matrix* Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored. |
[in,out] | U | magma_z_matrix* Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored. |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zparilut_residuals_gpu | ( | magma_z_matrix | A, |
magma_z_matrix | L, | ||
magma_z_matrix | U, | ||
magma_z_matrix * | R, | ||
magma_queue_t | queue ) |
This function computes the ILU residual in the locations included in the sparsity pattern of R.
[in] | A | magma_z_matrix System matrix. The format is sorted CSR. |
[in] | L | magma_z_matrix Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored. |
[in] | U | magma_z_matrix Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored. |
[in,out] | R | magma_z_matrix* Sparsity pattern on which the ILU residual is computed. R is in COO format. On output, R contains the ILU residual. |
[in] | queue | magma_queue_t Queue to execute in. |