![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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.
[in] | m | magma_int_t dimension m |
[in] | n | magma_int_t dimension n |
[in] | A | magma_s_matrix input matrix A |
[in] | lda | magma_int_t leading dimension matrix A |
[in,out] | Q | magma_s_matrix* input matrix Q |
[in,out] | R | magma_s_matrix* input matrix R |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b - can be a block |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | precond_par | magma_s_preconditioner* preconditioner |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | precond_par | magma_s_preconditioner* preconditioner |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | precond_par | magma_s_preconditioner* preconditioner |
[in] | queue | magma_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.
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in] | precond_par | magma_s_preconditioner* preconditioner |
[in] | queue | magma_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).
[in] | A | magma_s_matrix input matrix A |
[in] | b | magma_s_matrix RHS b |
[in,out] | x | magma_s_matrix* solution approximation |
[in,out] | solver_par | magma_s_solver_par* solver parameters |
[in,out] | precond_par | magma_s_preconditioner* preconditioner |
[in] | queue | magma_queue_t Queue to execute in. |