Iterative Krylov methods for hermitian and non-hermitian linear systems. More...
Functions | |
template<class LinearOperator , class Vector , class Monitor , class Preconditioner > | |
void | cusp::krylov::bicg (LinearOperator &A, LinearOperator &At, Vector &x, Vector &b, Monitor &monitor, Preconditioner &M, Preconditioner &Mt) |
Biconjugate Gradient method. More... | |
template<class LinearOperator , class Vector , class Monitor , class Preconditioner > | |
void | cusp::krylov::bicgstab (LinearOperator &A, Vector &x, Vector &b, Monitor &monitor, Preconditioner &M) |
Biconjugate Gradient Stabilized method. More... | |
template<class LinearOperator , class VectorType1 , class VectorType2 , class VectorType3 , class Monitor > | |
thrust::detail::enable_if_convertible < typename LinearOperator::format, cusp::known_format >::type | cusp::krylov::bicgstab_m (LinearOperator &A, VectorType1 &x, VectorType2 &b, VectorType3 &sigma, Monitor &monitor) |
Multi-mass Biconjugate Gradient stabilized method. More... | |
template<class LinearOperator , class Vector , class Monitor , class Preconditioner > | |
void | cusp::krylov::cg (LinearOperator &A, Vector &x, Vector &b, Monitor &monitor, Preconditioner &M) |
Conjugate Gradient method. More... | |
template<class LinearOperator , class VectorType1 , class VectorType2 , class VectorType3 , class Monitor > | |
thrust::detail::enable_if_convertible < typename LinearOperator::format, cusp::known_format >::type | cusp::krylov::cg_m (LinearOperator &A, VectorType1 &x, VectorType2 &b, VectorType3 &sigma, Monitor &monitor) |
Multi-mass Conjugate Gradient method. More... | |
template<class LinearOperator , class Vector , class Monitor , class Preconditioner > | |
void | cusp::krylov::cr (LinearOperator &A, Vector &x, Vector &b, Monitor &monitor, Preconditioner &M) |
Conjugate Residual method. More... | |
template<class LinearOperator , class Vector , class Monitor , class Preconditioner > | |
void | cusp::krylov::gmres (LinearOperator &A, Vector &x, Vector &b, const size_t restart, Monitor &monitor, Preconditioner &M) |
GMRES method. More... | |
void cusp::krylov::bicg | ( | LinearOperator & | A, |
LinearOperator & | At, | ||
Vector & | x, | ||
Vector & | b, | ||
Monitor & | monitor, | ||
Preconditioner & | M, | ||
Preconditioner & | Mt | ||
) |
Biconjugate Gradient method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
At | conjugate tranpose of the matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
monitor | montiors iteration and determines stopping conditions |
M | preconditioner for A |
Mt | conjugate tranpose of the preconditioner for A |
M
.The following code snippet demonstrates how to use bicg
to solve a 10x10 Poisson problem.
monitor
void cusp::krylov::bicgstab | ( | LinearOperator & | A, |
Vector & | x, | ||
Vector & | b, | ||
Monitor & | monitor, | ||
Preconditioner & | M | ||
) |
Biconjugate Gradient Stabilized method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
monitor | montiors iteration and determines stopping conditions |
M | preconditioner for A |
Solves the linear system A x = b with preconditioner M
.
The following code snippet demonstrates how to use bicgstab
to solve a 10x10 Poisson problem.
monitor
thrust::detail::enable_if_convertible<typename LinearOperator::format,cusp::known_format>::type cusp::krylov::bicgstab_m | ( | LinearOperator & | A, |
VectorType1 & | x, | ||
VectorType2 & | b, | ||
VectorType3 & | sigma, | ||
Monitor & | monitor | ||
) |
Multi-mass Biconjugate Gradient stabilized method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
sigma | array of shifts |
monitor | montiors iteration and determines stopping conditions |
This routine solves systems of the type
(A+sigma Id)x = b
for a number of different sigma, iteratively, for sparse A, without additional matrix-vector multiplication.
The following code snippet demonstrates how to use bicgstab
to solve a 10x10 Poisson problem.
monitor
void cusp::krylov::cg | ( | LinearOperator & | A, |
Vector & | x, | ||
Vector & | b, | ||
Monitor & | monitor, | ||
Preconditioner & | M | ||
) |
Conjugate Gradient method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
monitor | monitors iteration and determines stopping conditions |
M | preconditioner for A |
M
.A
and M
must be symmetric and positive-definite.cg
to solve a 10x10 Poisson problem.monitor
thrust::detail::enable_if_convertible<typename LinearOperator::format,cusp::known_format>::type cusp::krylov::cg_m | ( | LinearOperator & | A, |
VectorType1 & | x, | ||
VectorType2 & | b, | ||
VectorType3 & | sigma, | ||
Monitor & | monitor | ||
) |
Multi-mass Conjugate Gradient method.
A | matrix of the linear system |
x | solutions of the system |
b | right-hand side of the linear system |
sigma | array of shifts |
monitor | monitors interation and determines stoppoing conditions |
Solves the symmetric, positive-definited linear system (A+sigma) x = b for some set of constant shifts sigma
for the price of the smallest shift iteratively, for sparse A, without additional matrix-vector multiplication.
The following code snippet demonstrates how to use bicgstab
to solve a 10x10 Poisson problem.
monitor
void cusp::krylov::cr | ( | LinearOperator & | A, |
Vector & | x, | ||
Vector & | b, | ||
Monitor & | monitor, | ||
Preconditioner & | M | ||
) |
Conjugate Residual method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
monitor | montiors iteration and determines stopping conditions |
M | preconditioner for A |
A
and M
must be symmetric and semi-definite.cr
to solve a 10x10 Poisson problem.monitor
void cusp::krylov::gmres | ( | LinearOperator & | A, |
Vector & | x, | ||
Vector & | b, | ||
const size_t | restart, | ||
Monitor & | monitor, | ||
Preconditioner & | M | ||
) |
GMRES method.
LinearOperator | is a matrix or subclass of linear_operator |
Vector | vector |
Monitor | is a monitor such as default_monitor or verbose_monitor |
Preconditioner | is a matrix or subclass of linear_operator |
A | matrix of the linear system |
x | approximate solution of the linear system |
b | right-hand side of the linear system |
restart | the method every restart inner iterations |
monitor | montiors iteration and determines stopping conditions |
M | preconditioner for A |
M
.The following code snippet demonstrates how to use gmres
to solve a 10x10 Poisson problem.
monitor