|
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename Monitor , typename Preconditioner > |
| void | cusp::krylov::bicg (const LinearOperator &A, const LinearOperator &At, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M, Preconditioner &Mt) |
| | Biconjugate Gradient method.
|
| |
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename Monitor , typename Preconditioner > |
| void | cusp::krylov::bicgstab (const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M) |
| | Biconjugate Gradient Stabilized method.
|
| |
| template<class LinearOperator , class VectorType1 , class VectorType2 , class VectorType3 , class Monitor > |
| thrust::detail::enable_if_convertible_t< typename LinearOperator::format, cusp::known_format > | cusp::krylov::bicgstab_m (LinearOperator &A, VectorType1 &x, VectorType2 &b, VectorType3 &sigma, Monitor &monitor) |
| | Multi-mass Biconjugate Gradient stabilized method.
|
| |
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename Monitor , typename Preconditioner > |
| void | cusp::krylov::cg (const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M) |
| | Conjugate Gradient method.
|
| |
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename VectorType3 , typename Monitor > |
| void | cusp::krylov::cg_m (const LinearOperator &A, VectorType1 &x, const VectorType2 &b, const VectorType3 &sigma, Monitor &monitor) |
| | Multi-mass Conjugate Gradient method.
|
| |
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename Monitor , typename Preconditioner > |
| void | cusp::krylov::cr (const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M) |
| | Conjugate Residual method.
|
| |
| template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename Monitor , typename Preconditioner > |
| void | cusp::krylov::gmres (const LinearOperator &A, VectorType1 &x, const VectorType2 &b, const size_t restart, Monitor &monitor, Preconditioner &M) |
| | GMRES method.
|
| |
Iterative Krylov methods for hermitian and non-hermitian linear systems.
template<class LinearOperator , class VectorType1 , class VectorType2 , class VectorType3 , class Monitor >
| thrust::detail::enable_if_convertible_t< typename LinearOperator::format, cusp::known_format > cusp::krylov::bicgstab_m |
( |
LinearOperator & |
A, |
|
|
VectorType1 & |
x, |
|
|
VectorType2 & |
b, |
|
|
VectorType3 & |
sigma, |
|
|
Monitor & |
monitor |
|
) |
| |
Multi-mass Biconjugate Gradient stabilized method.
- Template Parameters
-
- Parameters
-
| 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 | monitors iteration and determines stopping conditions |
- Overview
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.
- See also
- http://arxiv.org/abs/hep-lat/9612014
- Example
The following code snippet demonstrates how to use bicgstab to solve a 10x10 Poisson problem.
int main(void)
{
size_t N_s = 4;
sigma[0] = 0.1;
sigma[1] = 0.5;
sigma[2] = 1.0;
sigma[3] = 5.0;
return 0;
}
Multi-mass Biconjugate Gradient stabilized (BiCGstab-M) method.
thrust::detail::enable_if_convertible_t< typename LinearOperator::format, cusp::known_format > bicgstab_m(LinearOperator &A, VectorType1 &x, VectorType2 &b, VectorType3 &sigma, Monitor &monitor)
Multi-mass Biconjugate Gradient stabilized method.
- See also
monitor
template<typename LinearOperator , typename VectorType1 , typename VectorType2 , typename VectorType3 , typename Monitor >
| void cusp::krylov::cg_m |
( |
const LinearOperator & |
A, |
|
|
VectorType1 & |
x, |
|
|
const VectorType2 & |
b, |
|
|
const VectorType3 & |
sigma, |
|
|
Monitor & |
monitor |
|
) |
| |
Multi-mass Conjugate Gradient method.
- Parameters
-
| 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 iteration and determines stopping conditions |
- Overview
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.
- See also
- http://arxiv.org/abs/hep-lat/9612014
- Example
The following code snippet demonstrates how to use bicgstab to solve a 10x10 Poisson problem.
int main(void)
{
size_t N_s = 4;
sigma[0] = 0.1;
sigma[1] = 0.5;
sigma[2] = 1.0;
sigma[3] = 5.0;
return 0;
}
Multi-mass Conjugate Gradient (CG-M) method.
void cg_m(const LinearOperator &A, VectorType1 &x, const VectorType2 &b, const VectorType3 &sigma, Monitor &monitor)
Multi-mass Conjugate Gradient method.
- See also
monitor