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

Functions

magma_int_t magma_sqr (magma_int_t m, magma_int_t n, magma_s_matrix A, magma_int_t lda, magma_s_matrix *Q, magma_s_matrix *R, magma_queue_t queue)
 This is a wrapper to call MAGMA QR on the data structure of sparse matrices. More...
 
magma_int_t magma_sbpcg (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_scg (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_scg_merge (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real matrix A. More...
 
magma_int_t magma_scg_res (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_scgs_merge (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_sidr_merge (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general real N-by-N matrix A. More...
 
magma_int_t magma_sidr_strms (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general real N-by-N matrix A. More...
 
magma_int_t magma_spcg_merge (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real matrix A. More...
 
magma_int_t magma_spidr_merge (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_spidr_strms (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 
magma_int_t magma_sptfqmr (magma_s_matrix A, magma_s_matrix b, magma_s_matrix *x, magma_s_solver_par *solver_par, magma_s_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A. More...
 

Detailed Description

Function Documentation

magma_int_t magma_sqr ( magma_int_t  m,
magma_int_t  n,
magma_s_matrix  A,
magma_int_t  lda,
magma_s_matrix *  Q,
magma_s_matrix *  R,
magma_queue_t  queue 
)

This is a wrapper to call MAGMA QR on the data structure of sparse matrices.

Output matrices Q and R reside on the same memory location as matrix A. On exit, Q is a M-by-N matrix in lda-by-N space. On exit, R is a min(M,N)-by-N upper trapezoidal matrix.

Parameters
[in]mmagma_int_t dimension m
[in]nmagma_int_t dimension n
[in]Amagma_s_matrix input matrix A
[in]ldamagma_int_t leading dimension matrix A
[in,out]Qmagma_s_matrix* input matrix Q
[in,out]Rmagma_s_matrix* input matrix R
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_sbpcg ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the block preconditioned Conjugate Gradient method.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b - can be a block
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]precond_parmagma_s_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_scg ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the Conjugate Gradient method.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_scg_merge ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real matrix A.

This is a GPU implementation of the Conjugate Gradient method in variant, where multiple operations are merged into one compute kernel.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_scg_res ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the preconditioned Conjugate Gradient method.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_scgs_merge ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the Conjugate Gradient Squared (CGS) method.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_sidr_merge ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a general real N-by-N matrix A.

This is a GPU implementation of the Induced Dimension Reduction method applying kernel fusion.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_sidr_strms ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a general real N-by-N matrix A.

This is a GPU implementation of the Induced Dimension Reduction method applying kernel fusion and kernel overlap.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_spcg_merge ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real matrix A.

This is a GPU implementation of the Conjugate Gradient method in variant, where multiple operations are merged into one compute kernel.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]precond_parmagma_s_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_spidr_merge ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the preconditioned Induced Dimension Reduction method applying kernel fusion.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]precond_parmagma_s_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_spidr_strms ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the preconditioned Induced Dimension Reduction method applying kernel fusion and kernel overlap.

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in]precond_parmagma_s_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_sptfqmr ( magma_s_matrix  A,
magma_s_matrix  b,
magma_s_matrix *  x,
magma_s_solver_par *  solver_par,
magma_s_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a real symmetric N-by-N positive definite matrix A.

This is a GPU implementation of the preconditioned transpose-free Quasi-Minimal Residual method (TFQMR).

Parameters
[in]Amagma_s_matrix input matrix A
[in]bmagma_s_matrix RHS b
[in,out]xmagma_s_matrix* solution approximation
[in,out]solver_parmagma_s_solver_par* solver parameters
[in,out]precond_parmagma_s_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.