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

Functions

magma_int_t magma_zbaiter (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * x = b via the block asynchronous iteration method on GPU. More...
 
magma_int_t magma_zbaiter_overlap (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * x = b via the block asynchronous iteration method on GPU. More...
 
magma_int_t magma_zbicg (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zbicgstab (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zbicgstab_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zbicgstab_merge2 (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zbicgstab_merge3 (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zbombard (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex general matrix A. More...
 
magma_int_t magma_zbombard_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex general matrix A. More...
 
magma_int_t magma_zcgs (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_zfgmres (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex sparse matrix stored in the GPU memory. More...
 
magma_int_t magma_ziterref (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex Hermitian N-by-N positive definite matrix A. More...
 
magma_int_t magma_zjacobi (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex Hermitian N-by-N positive definite matrix A. More...
 
magma_int_t magma_zjacobidomainoverlap (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex Hermitian N-by-N positive definite matrix A. More...
 
magma_int_t magma_zlsqr (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A*X=B for X if A is consistent, otherwise it attempts to solve the least squares solution X that minimizes norm(B-A*X). More...
 
magma_int_t magma_zpbicg (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zpbicgstab (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex N-by-N general matrix. More...
 
magma_int_t magma_zpbicgstab_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general matrix. More...
 
magma_int_t magma_zpcgs (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_zpcgs_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_zpqmr (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general complex matrix A. More...
 
magma_int_t magma_zpqmr_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general complex matrix A. More...
 
magma_int_t magma_zptfqmr_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_z_preconditioner *precond_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_zqmr (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a general complex matrix A. More...
 
magma_int_t magma_zqmr_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex general matrix A. More...
 
magma_int_t magma_ztfqmr (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_ztfqmr_merge (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 
magma_int_t magma_ztfqmr_unrolled (magma_z_matrix A, magma_z_matrix b, magma_z_matrix *x, magma_z_solver_par *solver_par, magma_queue_t queue)
 Solves a system of linear equations A * X = B where A is a complex matrix A. More...
 

Detailed Description

Function Documentation

magma_int_t magma_zbaiter ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * x = b via the block asynchronous iteration method on GPU.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbaiter_overlap ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * x = b via the block asynchronous iteration method on GPU.

It used restricted additive Schwarz overlap in top-down direction.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbicg ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Biconjugate Gradient method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbicgstab ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Biconjugate Gradient Stabilized method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbicgstab_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Biconjugate Gradient Stabilized method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbicgstab_merge2 ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Biconjugate Gradient Stabilized method. The difference to magma_zbicgstab is that we use specifically designed kernels merging multiple operations into one kernel.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbicgstab_merge3 ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Biconjugate Gradient Stabilized method. The difference to magma_zbicgstab is that we use specifically designed kernels merging multiple operations into one kernel.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbombard ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the iterative bombardment suggested in Barrett et al. ''Algorithmic bombardment for the iterative solution of linear systems: A poly-iterative approach'' using BiCGSTAB, CGS and MQR. At this point, it only works for symmetric A.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zbombard_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the iterative bombardment suggested in Barrett et al. ''Algorithmic bombardment for the iterative solution of linear systems: A poly-iterative approach'' using BiCGSTAB, CGS and MQR. At this point, it only works for symmetric A.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zcgs ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zfgmres ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B where A is a complex sparse matrix stored in the GPU memory.

X and B are complex vectors stored on the GPU memory. This is a GPU implementation of the right-preconditioned flexible GMRES.

Parameters
[in]Amagma_z_matrix descriptor for matrix A
[in]bmagma_z_matrix RHS b vector
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ziterref ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Iterative Refinement method. The inner solver is passed via the preconditioner argument.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in,out]precond_parmagma_z_preconditioner* inner solver
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zjacobi ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Jacobi method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zjacobidomainoverlap ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Jacobi method allowing for domain overlap.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zlsqr ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

Solves a system of linear equations A*X=B for X if A is consistent, otherwise it attempts to solve the least squares solution X that minimizes norm(B-A*X).

The N-by-P coefficient matrix A need not be square but the right hand side column vector B must have length N.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in,out]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpbicg ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in,out]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpbicgstab ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the preconditioned Biconjugate Gradient Stabelized method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpbicgstab_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the preconditioned Biconjugate Gradient Stabelized method.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpcgs ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpcgs_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpqmr ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpqmr_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the preconditioned Quasi-Minimal Residual method (QMR) using custom-designed kernels.

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zptfqmr_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_z_preconditioner *  precond_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in,out]precond_parmagma_z_preconditioner* preconditioner
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zqmr ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Quasi-Minimal Residual method (QMR).

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zqmr_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

This is a GPU implementation of the Quasi-Minimal Residual method (QMR).

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ztfqmr ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ztfqmr_merge ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_ztfqmr_unrolled ( magma_z_matrix  A,
magma_z_matrix  b,
magma_z_matrix *  x,
magma_z_solver_par *  solver_par,
magma_queue_t  queue 
)

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

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

Parameters
[in]Amagma_z_matrix input matrix A
[in]bmagma_z_matrix RHS b
[in,out]xmagma_z_matrix* solution approximation
[in,out]solver_parmagma_z_solver_par* solver parameters
[in]queuemagma_queue_t Queue to execute in.