MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages

Functions

magma_int_t magma_dbajac_csr (magma_int_t localiters, magma_d_matrix D, magma_d_matrix R, magma_d_matrix b, magma_d_matrix *x, magma_queue_t queue)
 This routine is a block-asynchronous Jacobi iteration performing s local Jacobi-updates within the block. More...
 
magma_int_t magma_dbajac_csr_overlap (magma_int_t localiters, magma_int_t matrices, magma_int_t overlap, magma_d_matrix *D, magma_d_matrix *R, magma_d_matrix b, magma_d_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. More...
 
magma_int_t magma_dcompact (magma_int_t m, magma_int_t n, magmaDouble_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. More...
 
magma_int_t magma_djacobisetup_vector_gpu (magma_int_t num_rows, magma_d_matrix b, magma_d_matrix d, magma_d_matrix c, magma_d_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. More...
 
magma_int_t magma_dlobpcg_maxpy (magma_int_t num_rows, magma_int_t num_vecs, magmaDouble_ptr X, magmaDouble_ptr Y, magma_queue_t queue)
 This routine computes a axpy for a mxn matrix: More...
 
magma_int_t magma_dbicgstab_1 (magma_int_t num_rows, magma_int_t num_cols, double beta, double omega, magmaDouble_ptr r, magmaDouble_ptr v, magmaDouble_ptr p, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgstab_2 (magma_int_t num_rows, magma_int_t num_cols, double alpha, magmaDouble_ptr r, magmaDouble_ptr v, magmaDouble_ptr s, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgstab_3 (magma_int_t num_rows, magma_int_t num_cols, double alpha, double omega, magmaDouble_ptr p, magmaDouble_ptr s, magmaDouble_ptr t, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgstab_4 (magma_int_t num_rows, magma_int_t num_cols, double alpha, double omega, magmaDouble_ptr y, magmaDouble_ptr z, magmaDouble_ptr s, magmaDouble_ptr t, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgmerge_spmv1 (magma_d_matrix A, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr dp, magmaDouble_ptr dr, magmaDouble_ptr dv, magmaDouble_ptr skp, magma_queue_t queue)
 Merges the first SpmV using CSR with the dot product and the computation of alpha. More...
 
magma_int_t magma_dbicgmerge_spmv2 (magma_d_matrix A, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr ds, magmaDouble_ptr dt, magmaDouble_ptr skp, magma_queue_t queue)
 Merges the second SpmV using CSR with the dot product and the computation of omega. More...
 
magma_int_t magma_dbicgmerge_xrbeta (magma_int_t n, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr rr, magmaDouble_ptr r, magmaDouble_ptr p, magmaDouble_ptr s, magmaDouble_ptr t, magmaDouble_ptr x, magmaDouble_ptr skp, magma_queue_t queue)
 Merges the second SpmV using CSR with the dot product and the computation of omega. More...
 
magma_int_t magma_dbicgmerge1 (magma_int_t n, magmaDouble_ptr skp, magmaDouble_ptr v, magmaDouble_ptr r, magmaDouble_ptr p, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgmerge2 (magma_int_t n, magmaDouble_ptr skp, magmaDouble_ptr r, magmaDouble_ptr v, magmaDouble_ptr s, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgmerge3 (magma_int_t n, magmaDouble_ptr skp, magmaDouble_ptr p, magmaDouble_ptr s, magmaDouble_ptr t, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dbicgmerge4 (magma_int_t type, magmaDouble_ptr skp, magma_queue_t queue)
 Performs some parameter operations for the BiCGSTAB with scalars on GPU. More...
 
magma_int_t magma_dmergeblockkrylov (magma_int_t num_rows, magma_int_t num_cols, magmaDouble_ptr alpha, magmaDouble_ptr p, magmaDouble_ptr x, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dcgmerge_spmv1 (magma_d_matrix A, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr dd, magmaDouble_ptr dz, magmaDouble_ptr skp, magma_queue_t queue)
 Merges the first SpmV using different formats with the dot product and the computation of rho. More...
 
magma_int_t magma_dcgs_1 (magma_int_t num_rows, magma_int_t num_cols, double beta, magmaDouble_ptr r, magmaDouble_ptr q, magmaDouble_ptr u, magmaDouble_ptr p, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dcgs_2 (magma_int_t num_rows, magma_int_t num_cols, magmaDouble_ptr r, magmaDouble_ptr u, magmaDouble_ptr p, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dcgs_3 (magma_int_t num_rows, magma_int_t num_cols, double alpha, magmaDouble_ptr v_hat, magmaDouble_ptr u, magmaDouble_ptr q, magmaDouble_ptr t, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dcgs_4 (magma_int_t num_rows, magma_int_t num_cols, double alpha, magmaDouble_ptr u_hat, magmaDouble_ptr t, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_didr_smoothing_1 (magma_int_t num_rows, magma_int_t num_cols, magmaDouble_ptr drs, magmaDouble_ptr dr, magmaDouble_ptr dt, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_didr_smoothing_2 (magma_int_t num_rows, magma_int_t num_cols, double omega, magmaDouble_ptr dx, magmaDouble_ptr dxs, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_1 (magma_int_t num_rows, magma_int_t num_cols, double rho, double psi, magmaDouble_ptr y, magmaDouble_ptr z, magmaDouble_ptr v, magmaDouble_ptr w, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_2 (magma_int_t num_rows, magma_int_t num_cols, double pde, double rde, magmaDouble_ptr y, magmaDouble_ptr z, magmaDouble_ptr p, magmaDouble_ptr q, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_3 (magma_int_t num_rows, magma_int_t num_cols, double beta, magmaDouble_ptr pt, magmaDouble_ptr v, magmaDouble_ptr y, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_4 (magma_int_t num_rows, magma_int_t num_cols, double eta, magmaDouble_ptr p, magmaDouble_ptr pt, magmaDouble_ptr d, magmaDouble_ptr s, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_5 (magma_int_t num_rows, magma_int_t num_cols, double eta, double pds, magmaDouble_ptr p, magmaDouble_ptr pt, magmaDouble_ptr d, magmaDouble_ptr s, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_6 (magma_int_t num_rows, magma_int_t num_cols, double beta, double rho, double psi, magmaDouble_ptr y, magmaDouble_ptr z, magmaDouble_ptr v, magmaDouble_ptr w, magmaDouble_ptr wt, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_7 (magma_int_t num_rows, magma_int_t num_cols, double beta, magmaDouble_ptr pt, magmaDouble_ptr v, magmaDouble_ptr vt, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dqmr_8 (magma_int_t num_rows, magma_int_t num_cols, double rho, double psi, magmaDouble_ptr vt, magmaDouble_ptr wt, magmaDouble_ptr y, magmaDouble_ptr z, magmaDouble_ptr v, magmaDouble_ptr w, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dtfqmr_1 (magma_int_t num_rows, magma_int_t num_cols, double alpha, double sigma, magmaDouble_ptr v, magmaDouble_ptr Au, magmaDouble_ptr u_m, magmaDouble_ptr pu_m, magmaDouble_ptr u_mp1, magmaDouble_ptr w, magmaDouble_ptr d, magmaDouble_ptr Ad, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dtfqmr_2 (magma_int_t num_rows, magma_int_t num_cols, double eta, magmaDouble_ptr d, magmaDouble_ptr Ad, magmaDouble_ptr x, magmaDouble_ptr r, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dtfqmr_3 (magma_int_t num_rows, magma_int_t num_cols, double beta, magmaDouble_ptr w, magmaDouble_ptr u_m, magmaDouble_ptr u_mp1, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dtfqmr_4 (magma_int_t num_rows, magma_int_t num_cols, double beta, magmaDouble_ptr Au_new, magmaDouble_ptr v, magmaDouble_ptr Au, magma_queue_t queue)
 Merges multiple operations into one kernel: More...
 
magma_int_t magma_dtfqmr_5 (magma_int_t num_rows, magma_int_t num_cols, double alpha, double sigma, magmaDouble_ptr v, magmaDouble_ptr Au, magmaDouble_ptr u_mp1, magmaDouble_ptr w, magmaDouble_ptr d, magmaDouble_ptr Ad, magma_queue_t queue)
 Mergels multiple operations into one kernel: More...
 
magma_int_t magma_dparic_csr (magma_d_matrix A, magma_d_matrix A_CSR, magma_queue_t queue)
 This routine iteratively computes an incomplete LU factorization. More...
 
magma_int_t magma_dparilu_csr (magma_d_matrix A, magma_d_matrix L, magma_d_matrix U, magma_queue_t queue)
 This routine iteratively computes an incomplete LU factorization. More...
 
magma_int_t magma_dparilut_sweep_gpu (magma_d_matrix *A, magma_d_matrix *L, magma_d_matrix *U, magma_queue_t queue)
 This function does an ParILUT sweep. More...
 
magma_int_t magma_dparilut_residuals_gpu (magma_d_matrix A, magma_d_matrix L, magma_d_matrix U, magma_d_matrix *R, magma_queue_t queue)
 This function computes the ILU residual in the locations included in the sparsity pattern of R. More...
 
magma_int_t magma_dvalinit_gpu (magma_int_t num_el, magmaDouble_ptr dval, magma_queue_t queue)
 Initializes a device array with zero. More...
 
magma_int_t magma_dindexinit_gpu (magma_int_t num_el, magmaIndex_ptr dind, magma_queue_t queue)
 Initializes a device array with zero. More...
 

Detailed Description

Function Documentation

magma_int_t magma_dbajac_csr ( magma_int_t  localiters,
magma_d_matrix  D,
magma_d_matrix  R,
magma_d_matrix  b,
magma_d_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.

Parameters
[in]localitersmagma_int_t number of local Jacobi-like updates
[in]Dmagma_d_matrix input matrix with diagonal blocks
[in]Rmagma_d_matrix input matrix with non-diagonal parts
[in]bmagma_d_matrix RHS
[in]xmagma_d_matrix* iterate/solution
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbajac_csr_overlap ( magma_int_t  localiters,
magma_int_t  matrices,
magma_int_t  overlap,
magma_d_matrix *  D,
magma_d_matrix *  R,
magma_d_matrix  b,
magma_d_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.

Parameters
[in]localitersmagma_int_t number of local Jacobi-like updates
[in]matricesmagma_int_t number of sub-matrices
[in]overlapmagma_int_t size of the overlap
[in]Dmagma_d_matrix* set of matrices with diagonal blocks
[in]Rmagma_d_matrix* set of matrices with non-diagonal parts
[in]bmagma_d_matrix RHS
[in]xmagma_d_matrix* iterate/solution
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcompact ( magma_int_t  m,
magma_int_t  n,
magmaDouble_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.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in,out]dACOMPLEX DOUBLE PRECISION array, dimension (LDDA,N) The m by n matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[in]dnormsDOUBLE PRECISION array, dimension N The norms of the N vectors in dA
[in]tolDOUBLE PRECISON The tolerance value used in the criteria to compact or not.
[in,out]activeINTEGER array, dimension N A mask of 1s and 0s showing if a vector remains or has been removed
[in,out]cBlockmagmaInt_ptr The number of vectors that remain in dA (i.e., with norms > tol).
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_djacobisetup_vector_gpu ( magma_int_t  num_rows,
magma_d_matrix  b,
magma_d_matrix  d,
magma_d_matrix  c,
magma_d_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

Parameters
[in]num_rowsmagma_int_t number of rows
[in]bmagma_d_matrix RHS b
[in]dmagma_d_matrix vector with diagonal entries
[out]cmagma_d_matrix* c = D^(-1) * b
[out]xmagma_d_matrix* iteration vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dlobpcg_maxpy ( magma_int_t  num_rows,
magma_int_t  num_vecs,
magmaDouble_ptr  X,
magmaDouble_ptr  Y,
magma_queue_t  queue 
)

This routine computes a axpy for a mxn matrix:

Y = X + Y

It replaces: magma_daxpy(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] /

Parameters
[in]num_rowsmagma_int_t number of rows
[in]num_vecsmagma_int_t number of vectors
[in]XmagmaDouble_ptr input vector X
[in,out]YmagmaDouble_ptr input/output vector Y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgstab_1 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
double  omega,
magmaDouble_ptr  r,
magmaDouble_ptr  v,
magmaDouble_ptr  p,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

p = r + beta * ( p - omega * v )

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]omegadouble scalar
[in]rmagmaDouble_ptr vector
[in]vmagmaDouble_ptr vector
[in,out]pmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgstab_2 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
magmaDouble_ptr  r,
magmaDouble_ptr  v,
magmaDouble_ptr  s,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

s = r - alpha v

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]rmagmaDouble_ptr vector
[in]vmagmaDouble_ptr vector
[in,out]smagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgstab_3 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
double  omega,
magmaDouble_ptr  p,
magmaDouble_ptr  s,
magmaDouble_ptr  t,
magmaDouble_ptr  x,
magmaDouble_ptr  r,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

x = x + alpha * p + omega * s r = s - omega * t

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]omegadouble scalar
[in]pmagmaDouble_ptr vector
[in]smagmaDouble_ptr vector
[in]tmagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgstab_4 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
double  omega,
magmaDouble_ptr  y,
magmaDouble_ptr  z,
magmaDouble_ptr  s,
magmaDouble_ptr  t,
magmaDouble_ptr  x,
magmaDouble_ptr  r,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

x = x + alpha * y + omega * z r = s - omega * t

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]omegadouble scalar
[in]ymagmaDouble_ptr vector
[in]zmagmaDouble_ptr vector
[in]smagmaDouble_ptr vector
[in]tmagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgmerge_spmv1 ( magma_d_matrix  A,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  dp,
magmaDouble_ptr  dr,
magmaDouble_ptr  dv,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Merges the first SpmV using CSR with the dot product and the computation of alpha.

Parameters
[in]Amagma_d_matrix system matrix
[in]d1magmaDouble_ptr temporary vector
[in]d2magmaDouble_ptr temporary vector
[in]dpmagmaDouble_ptr input vector p
[in]drmagmaDouble_ptr input vector r
[in]dvmagmaDouble_ptr output vector v
[in,out]skpmagmaDouble_ptr array for parameters ( skp[0]=alpha )
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgmerge_spmv2 ( magma_d_matrix  A,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  ds,
magmaDouble_ptr  dt,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Merges the second SpmV using CSR with the dot product and the computation of omega.

Parameters
[in]Amagma_d_matrix input matrix
[in]d1magmaDouble_ptr temporary vector
[in]d2magmaDouble_ptr temporary vector
[in]dsmagmaDouble_ptr input vector s
[in]dtmagmaDouble_ptr output vector t
[in,out]skpmagmaDouble_ptr array for parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgmerge_xrbeta ( magma_int_t  n,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  rr,
magmaDouble_ptr  r,
magmaDouble_ptr  p,
magmaDouble_ptr  s,
magmaDouble_ptr  t,
magmaDouble_ptr  x,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Merges the second SpmV using CSR with the dot product and the computation of omega.

Parameters
[in]nint dimension n
[in]d1magmaDouble_ptr temporary vector
[in]d2magmaDouble_ptr temporary vector
[in]rrmagmaDouble_ptr input vector rr
[in]rmagmaDouble_ptr input/output vector r
[in]pmagmaDouble_ptr input vector p
[in]smagmaDouble_ptr input vector s
[in]tmagmaDouble_ptr input vector t
[out]xmagmaDouble_ptr output vector x
[in]skpmagmaDouble_ptr array for parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgmerge1 ( magma_int_t  n,
magmaDouble_ptr  skp,
magmaDouble_ptr  v,
magmaDouble_ptr  r,
magmaDouble_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 )

Parameters
[in]nint dimension n
[in]skpmagmaDouble_ptr set of scalar parameters
[in]vmagmaDouble_ptr input vector v
[in]rmagmaDouble_ptr input vector r
[in,out]pmagmaDouble_ptr input/output vector p
[in]queuemagma_queue_t queue to execute in.
magma_int_t magma_dbicgmerge2 ( magma_int_t  n,
magmaDouble_ptr  skp,
magmaDouble_ptr  r,
magmaDouble_ptr  v,
magmaDouble_ptr  s,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

s=r s=s-alpha*v

-> s = r - alpha * v

Parameters
[in]nint dimension n
[in]skpmagmaDouble_ptr set of scalar parameters
[in]rmagmaDouble_ptr input vector r
[in]vmagmaDouble_ptr input vector v
[out]smagmaDouble_ptr output vector s
[in]queuemagma_queue_t queue to execute in.
magma_int_t magma_dbicgmerge3 ( magma_int_t  n,
magmaDouble_ptr  skp,
magmaDouble_ptr  p,
magmaDouble_ptr  s,
magmaDouble_ptr  t,
magmaDouble_ptr  x,
magmaDouble_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

Parameters
[in]nint dimension n
[in]skpmagmaDouble_ptr set of scalar parameters
[in]pmagmaDouble_ptr input p
[in]smagmaDouble_ptr input s
[in]tmagmaDouble_ptr input t
[in,out]xmagmaDouble_ptr input/output x
[in,out]rmagmaDouble_ptr input/output r
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dbicgmerge4 ( magma_int_t  type,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Performs some parameter operations for the BiCGSTAB with scalars on GPU.

Parameters
[in]typeint kernel type
[in,out]skpmagmaDouble_ptr vector with parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmergeblockkrylov ( magma_int_t  num_rows,
magma_int_t  num_cols,
magmaDouble_ptr  alpha,
magmaDouble_ptr  p,
magmaDouble_ptr  x,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

v = y / rho y = y / rho w = wt / psi z = z / psi

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphamagmaDouble_ptr matrix containing all SKP
[in]pmagmaDouble_ptr search directions
[in,out]xmagmaDouble_ptr approximation vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcgmerge_spmv1 ( magma_d_matrix  A,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  dd,
magmaDouble_ptr  dz,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Merges the first SpmV using different formats with the dot product and the computation of rho.

Parameters
[in]Amagma_d_matrix input matrix
[in]d1magmaDouble_ptr temporary vector
[in]d2magmaDouble_ptr temporary vector
[in]ddmagmaDouble_ptr input vector d
[out]dzmagmaDouble_ptr input vector z
[out]skpmagmaDouble_ptr array for parameters ( skp[3]=rho )
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcgs_1 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
magmaDouble_ptr  r,
magmaDouble_ptr  q,
magmaDouble_ptr  u,
magmaDouble_ptr  p,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

u = r + beta q p = u + beta*(q + beta*p)

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]rmagmaDouble_ptr vector
[in]qmagmaDouble_ptr vector
[in,out]umagmaDouble_ptr vector
[in,out]pmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcgs_2 ( magma_int_t  num_rows,
magma_int_t  num_cols,
magmaDouble_ptr  r,
magmaDouble_ptr  u,
magmaDouble_ptr  p,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

u = r p = r

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]rmagmaDouble_ptr vector
[in,out]umagmaDouble_ptr vector
[in,out]pmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcgs_3 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
magmaDouble_ptr  v_hat,
magmaDouble_ptr  u,
magmaDouble_ptr  q,
magmaDouble_ptr  t,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

q = u - alpha v_hat t = u + q

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]v_hatmagmaDouble_ptr vector
[in]umagmaDouble_ptr vector
[in,out]qmagmaDouble_ptr vector
[in,out]tmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcgs_4 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
magmaDouble_ptr  u_hat,
magmaDouble_ptr  t,
magmaDouble_ptr  x,
magmaDouble_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

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]u_hatmagmaDouble_ptr vector
[in]tmagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_didr_smoothing_1 ( magma_int_t  num_rows,
magma_int_t  num_cols,
magmaDouble_ptr  drs,
magmaDouble_ptr  dr,
magmaDouble_ptr  dt,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

dt = drs - dr

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]drsmagmaDouble_ptr vector
[in]drmagmaDouble_ptr vector
[in,out]dtmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_didr_smoothing_2 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  omega,
magmaDouble_ptr  dx,
magmaDouble_ptr  dxs,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

dxs = dxs - gamma*(dxs-dx)

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]omegadouble scalar
[in]dxmagmaDouble_ptr vector
[in,out]dxsmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_1 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  rho,
double  psi,
magmaDouble_ptr  y,
magmaDouble_ptr  z,
magmaDouble_ptr  v,
magmaDouble_ptr  w,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

v = y / rho y = y / rho w = wt / psi z = z / psi

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]rhodouble scalar
[in]psidouble scalar
[in,out]ymagmaDouble_ptr vector
[in,out]zmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]wmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_2 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  pde,
double  rde,
magmaDouble_ptr  y,
magmaDouble_ptr  z,
magmaDouble_ptr  p,
magmaDouble_ptr  q,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

p = y - pde * p q = z - rde * q

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]pdedouble scalar
[in]rdedouble scalar
[in]ymagmaDouble_ptr vector
[in]zmagmaDouble_ptr vector
[in,out]pmagmaDouble_ptr vector
[in,out]qmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_3 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
magmaDouble_ptr  pt,
magmaDouble_ptr  v,
magmaDouble_ptr  y,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

v = pt - beta * v y = v

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]ptmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]ymagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_4 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  eta,
magmaDouble_ptr  p,
magmaDouble_ptr  pt,
magmaDouble_ptr  d,
magmaDouble_ptr  s,
magmaDouble_ptr  x,
magmaDouble_ptr  r,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

d = eta * p; s = eta * pt; x = x + d; r = r - s;

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]etadouble scalar
[in]pmagmaDouble_ptr vector
[in]ptmagmaDouble_ptr vector
[in,out]dmagmaDouble_ptr vector
[in,out]smagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_5 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  eta,
double  pds,
magmaDouble_ptr  p,
magmaDouble_ptr  pt,
magmaDouble_ptr  d,
magmaDouble_ptr  s,
magmaDouble_ptr  x,
magmaDouble_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;

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]etadouble scalar
[in]pdsdouble scalar
[in]pmagmaDouble_ptr vector
[in]ptmagmaDouble_ptr vector
[in,out]dmagmaDouble_ptr vector
[in,out]smagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_6 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
double  rho,
double  psi,
magmaDouble_ptr  y,
magmaDouble_ptr  z,
magmaDouble_ptr  v,
magmaDouble_ptr  w,
magmaDouble_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

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]rhodouble scalar
[in]psidouble scalar
[in,out]ymagmaDouble_ptr vector
[in,out]zmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]wmagmaDouble_ptr vector
[in,out]wtmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_7 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
magmaDouble_ptr  pt,
magmaDouble_ptr  v,
magmaDouble_ptr  vt,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

vt = pt - beta * v

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]ptmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]vtmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dqmr_8 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  rho,
double  psi,
magmaDouble_ptr  vt,
magmaDouble_ptr  wt,
magmaDouble_ptr  y,
magmaDouble_ptr  z,
magmaDouble_ptr  v,
magmaDouble_ptr  w,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

v = y / rho y = y / rho w = wt / psi z = z / psi

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]rhodouble scalar
[in]psidouble scalar
[in]vtmagmaDouble_ptr vector
[in]wtmagmaDouble_ptr vector
[in,out]ymagmaDouble_ptr vector
[in,out]zmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]wmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dtfqmr_1 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
double  sigma,
magmaDouble_ptr  v,
magmaDouble_ptr  Au,
magmaDouble_ptr  u_m,
magmaDouble_ptr  pu_m,
magmaDouble_ptr  u_mp1,
magmaDouble_ptr  w,
magmaDouble_ptr  d,
magmaDouble_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;

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]sigmadouble scalar
[in]vmagmaDouble_ptr vector
[in]AumagmaDouble_ptr vector
[in,out]u_mmagmaDouble_ptr vector
[in,out]pu_mmagmaDouble_ptr vector
[in,out]u_mp1magmaDouble_ptr vector
[in,out]wmagmaDouble_ptr vector
[in,out]dmagmaDouble_ptr vector
[in,out]AdmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dtfqmr_2 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  eta,
magmaDouble_ptr  d,
magmaDouble_ptr  Ad,
magmaDouble_ptr  x,
magmaDouble_ptr  r,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

x = x + eta * d r = r - eta * Ad

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]etadouble scalar
[in]dmagmaDouble_ptr vector
[in]AdmagmaDouble_ptr vector
[in,out]xmagmaDouble_ptr vector
[in,out]rmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dtfqmr_3 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
magmaDouble_ptr  w,
magmaDouble_ptr  u_m,
magmaDouble_ptr  u_mp1,
magma_queue_t  queue 
)

Mergels multiple operations into one kernel:

u_mp1 = w + beta*u_mp1

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]wmagmaDouble_ptr vector
[in]u_mmagmaDouble_ptr vector
[in,out]u_mp1magmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dtfqmr_4 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  beta,
magmaDouble_ptr  Au_new,
magmaDouble_ptr  v,
magmaDouble_ptr  Au,
magma_queue_t  queue 
)

Merges multiple operations into one kernel:

v = Au_new + beta*(Au+beta*v); Au = Au_new

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]betadouble scalar
[in]Au_newmagmaDouble_ptr vector
[in,out]vmagmaDouble_ptr vector
[in,out]AumagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dtfqmr_5 ( magma_int_t  num_rows,
magma_int_t  num_cols,
double  alpha,
double  sigma,
magmaDouble_ptr  v,
magmaDouble_ptr  Au,
magmaDouble_ptr  u_mp1,
magmaDouble_ptr  w,
magmaDouble_ptr  d,
magmaDouble_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;

Parameters
[in]num_rowsmagma_int_t dimension m
[in]num_colsmagma_int_t dimension n
[in]alphadouble scalar
[in]sigmadouble scalar
[in]vmagmaDouble_ptr vector
[in]AumagmaDouble_ptr vector
[in,out]u_mp1magmaDouble_ptr vector
[in,out]wmagmaDouble_ptr vector
[in,out]dmagmaDouble_ptr vector
[in,out]AdmagmaDouble_ptr vector
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dparic_csr ( magma_d_matrix  A,
magma_d_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.

Parameters
[in]Amagma_d_matrix input matrix A - initial guess (lower triangular)
[in,out]A_CSRmagma_d_matrix input/output matrix containing the IC approximation
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dparilu_csr ( magma_d_matrix  A,
magma_d_matrix  L,
magma_d_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.

Parameters
[in]Amagma_d_matrix input matrix A determing initial guess & processing order
[in,out]Lmagma_d_matrix input/output matrix L containing the lower triangular factor
[in,out]Umagma_d_matrix input/output matrix U containing the upper triangular factor
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dparilut_sweep_gpu ( magma_d_matrix *  A,
magma_d_matrix *  L,
magma_d_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.

Parameters
[in]Amagma_d_matrix* System matrix. The format is sorted CSR.
[in,out]Lmagma_d_matrix* Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored.
[in,out]Umagma_d_matrix* Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dparilut_residuals_gpu ( magma_d_matrix  A,
magma_d_matrix  L,
magma_d_matrix  U,
magma_d_matrix *  R,
magma_queue_t  queue 
)

This function computes the ILU residual in the locations included in the sparsity pattern of R.

Parameters
[in]Amagma_d_matrix System matrix. The format is sorted CSR.
[in]Lmagma_d_matrix Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored.
[in]Umagma_d_matrix Current approximation for the lower triangular factor The format is MAGMA_CSRCOO. This is sorted CSR plus the rowindexes being stored.
[in,out]Rmagma_d_matrix* Sparsity pattern on which the ILU residual is computed. R is in COO format. On output, R contains the ILU residual.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dvalinit_gpu ( magma_int_t  num_el,
magmaDouble_ptr  dval,
magma_queue_t  queue 
)

Initializes a device array with zero.

Parameters
[in]num_elmagma_int_t size of array
[in,out]dvalmagmaDouble_ptr array to initialize
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dindexinit_gpu ( magma_int_t  num_el,
magmaIndex_ptr  dind,
magma_queue_t  queue 
)

Initializes a device array with zero.

Parameters
[in]num_elmagma_int_t size of array
[in,out]dindmagmaIndex_ptr array to initialize
[in]queuemagma_queue_t Queue to execute in.