CUSP
Loading...
Searching...
No Matches
Krylov Methods

Iterative Krylov methods for hermitian and non-hermitian linear systems. More...

Functions

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.
 

Detailed Description

Iterative Krylov methods for hermitian and non-hermitian linear systems.

Function Documentation

◆ bicg()

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 Parameters
LinearOperatoris a matrix or subclass of linear_operator
VectorType1vector
Monitoris a monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
Atconjugate transpose of the matrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
Mtconjugate transpose of the preconditioner for A
Overview
Solves the linear system A x = b with preconditioner M.
Example

The following code snippet demonstrates how to use bicg to solve a 10x10 Poisson problem.

#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// set preconditioner (identity)
// solve the linear system A x = b
// because both A and M are hermitian we can use
// them for their own conjugate transpose
cusp::krylov::bicg(A, A, x, b, monitor, M, M);
return 0;
}
Biconjugate Gradient (BiCG) method.
The array1d class is a 1D vector container that may contain elements stored in "host" or "device" mem...
Definition array1d.h:99
Compressed sparse row (CSR) representation a sparse matrix.
Definition csr_matrix.h:108
Simple identity operator.
Implements standard convergence criteria and reporting for iterative solvers.
Definition monitor.h:102
Compressed Sparse Row matrix format.
void bicg(const LinearOperator &A, const LinearOperator &At, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M, Preconditioner &Mt)
Biconjugate Gradient method.
Monitor iterative solver convergence.
Poisson matrix generators.
See also
monitor

◆ bicgstab()

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 Parameters
LinearOperatoris a matrix or subclass of linear_operator
VectorType1vector
Monitoris a monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
Overview

Solves the linear system A x = b with preconditioner M.

Example

The following code snippet demonstrates how to use bicgstab to solve a 10x10 Poisson problem.

#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// set preconditioner (identity)
// solve the linear system A x = b
return 0;
}
Biconjugate Gradient Stabilized (BiCGstab) method.
void bicgstab(const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M)
Biconjugate Gradient Stabilized method.
See also
monitor

◆ bicgstab_m()

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
LinearOperatoris a matrix or subclass of linear_operator
Vectorvector
Monitoris a monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
sigmaarray of shifts
monitormonitors 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.

#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
size_t N_s = 4;
// set sigma values
sigma[0] = 0.1;
sigma[1] = 0.5;
sigma[2] = 1.0;
sigma[3] = 5.0;
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// solve the linear system A x = b
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

◆ cg()

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 Parameters
LinearOperatoris a matrix or subclass of linear_operator
VectorType1x input vector type
VectorType2b output vector type
Monitoris a monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
Overview
Solves the symmetric, positive-definite linear system A x = b with preconditioner M.
Note
A and M must be symmetric and positive-definite.
Example
The following code snippet demonstrates how to use cg to solve a 10x10 Poisson problem.
#include <cusp/monitor.h>
#include <cusp/krylov/cg.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// set preconditioner (identity)
// solve the linear system A x = b
cusp::krylov::cg(A, x, b, monitor, M);
return 0;
}
Conjugate Gradient (CG) method.
void cg(const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M)
Conjugate Gradient method.
See also
monitor

◆ cg_m()

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
Amatrix of the linear system
xsolutions of the system
bright-hand side of the linear system
sigmaarray of shifts
monitormonitors 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.

#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
size_t N_s = 4;
// set sigma values
sigma[0] = 0.1;
sigma[1] = 0.5;
sigma[2] = 1.0;
sigma[3] = 5.0;
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// solve the linear system A x = b
cusp::krylov::cg_m(A, x, b, sigma, monitor);
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

◆ cr()

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 Parameters
LinearOperatoris a matrix or subclass of linear_operator
VectorType1vector
Monitoris a monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
Overview
Solves a linear system using the conjugate residual method
Note
A and M must be symmetric and semi-definite.
Example
The following code snippet demonstrates how to use cr to solve a 10x10 Poisson problem.
#include <cusp/monitor.h>
#include <cusp/krylov/cr.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
// set preconditioner (identity)
// solve the linear system A x = b
cusp::krylov::cr(A, x, b, monitor, M);
return 0;
}
Conjugate Residual (CR) method.
void cr(const LinearOperator &A, VectorType1 &x, const VectorType2 &b, Monitor &monitor, Preconditioner &M)
Conjugate Residual method.
See also
monitor

◆ gmres()

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.

Template Parameters
LinearOperatoris a matrix or subclass of linear_operator
VectorType1vector
Monitoris a monitor such as default_monitor or verbose_monitor
Preconditioneris a matrix or subclass of linear_operator
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
restartthe method every restart inner iterations
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
Overview
Solves the nonsymmetric, linear system A x = b with preconditioner M.
Example

The following code snippet demonstrates how to use gmres to solve a 10x10 Poisson problem.

#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage for solution (x) and right hand side (b)
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<float> monitor(b, 100, 1e-6, 0, true);
int restart = 50;
// set preconditioner (identity)
// solve the linear system A x = b
cusp::krylov::gmres(A, x, b,restart, monitor, M);
return 0;
}
Generalized Minimum Residual (GMRES) method.
void gmres(const LinearOperator &A, VectorType1 &x, const VectorType2 &b, const size_t restart, Monitor &monitor, Preconditioner &M)
GMRES method.
See also
monitor