Fork me on GitHub
 All Classes Files Functions Variables Groups Pages
Public Methods | List of all members
cusp::relaxation::jacobi< ValueType, MemorySpace > Class Template Reference

Detailed description

template<typename ValueType, typename MemorySpace>
class cusp::relaxation::jacobi< ValueType, MemorySpace >

Represents a Jacobi relaxation scheme.

Template Parameters
ValueTypevalue_type of the array
MemorySpacememory space of the array (cusp::host_memory or cusp::device_memory)
Overview
Extracts the matrix diagonal and performs weighted Jacobi relaxation
Example
#include <cusp/array1d.h>
#include <cusp/monitor.h>
#include <cusp/blas/blas.h>
// include cusp Jacobi header file
int main()
{
// Construct 5-pt Poisson example
// Initialize data
// Allocate temporaries
// Construct Jacobi relaxation class
// Compute initial residual
cusp::multiply(A, x, r);
cusp::blas::axpy(b, r, float(-1));
// Construct monitor with stopping criteria of 100 iterations or 1e-4 residual error
cusp::monitor<float> monitor(b, 100, 1e-4, 0, true);
// Iteratively solve system
while (!monitor.finished(r))
{
M(A, b, x);
cusp::multiply(A, x, r);
cusp::blas::axpy(b, r, float(-1));
++monitor;
}
}

Definition at line 95 of file jacobi.h.

#include <jacobi.h>

Inheritance diagram for cusp::relaxation::jacobi< ValueType, MemorySpace >:
cusp::linear_operator<ValueType, MemorySpace >

Public Methods

 jacobi (void)
 
template<typename MatrixType >
 jacobi (const MatrixType &A, ValueType omega=1.0)
 
template<typename MemorySpace2 >
 jacobi (const jacobi< ValueType, MemorySpace2 > &A)
 
template<typename MatrixType , typename VectorType1 , typename VectorType2 >
void operator() (const MatrixType &A, const VectorType1 &b, VectorType2 &x)
 
template<typename MatrixType , typename VectorType1 , typename VectorType2 >
void operator() (const MatrixType &A, const VectorType1 &b, VectorType2 &x, const ValueType omega)
 
- Public Methods inherited from cusp::linear_operator<ValueType, MemorySpace >
 linear_operator (void)
 
 linear_operator (intnum_rows, intnum_cols)
 
 linear_operator (intnum_rows, intnum_cols, intnum_entries)
 

Constructor & Destructor Documentation

template<typename ValueType, typename MemorySpace>
cusp::relaxation::jacobi< ValueType, MemorySpace >::jacobi ( void  )
inline

This constructor creates an empty jacobi smoother.

Definition at line 107 of file jacobi.h.

template<typename ValueType, typename MemorySpace>
template<typename MatrixType >
cusp::relaxation::jacobi< ValueType, MemorySpace >::jacobi ( const MatrixType &  A,
ValueType  omega = 1.0 
)

This constructor creates a jacobi smoother using a given matrix and sweeping strategy (FORWARD, BACKWARD, SYMMETRIC).

Template Parameters
MatrixTypeType of input matrix used to create this jacobi smoother.
Parameters
AInput matrix used to create smoother.
omegaDamping factor used in Jacobi smoother.
template<typename ValueType, typename MemorySpace>
template<typename MemorySpace2 >
cusp::relaxation::jacobi< ValueType, MemorySpace >::jacobi ( const jacobi< ValueType, MemorySpace2 > &  A)
inline

Copy constructor for jacobi smoother.

Template Parameters
MemorySpace2Memory space of input jacobi smoother.
Parameters
AInput jacobi smoother.

Definition at line 128 of file jacobi.h.

Member Function Documentation

template<typename ValueType, typename MemorySpace>
template<typename MatrixType , typename VectorType1 , typename VectorType2 >
void cusp::relaxation::jacobi< ValueType, MemorySpace >::operator() ( const MatrixType &  A,
const VectorType1 &  b,
VectorType2 &  x 
)

Perform Jacobi relaxation using default omega damping factor specified during construction of this jacobi smoother

Template Parameters
MatrixTypeType of input matrix.
VectorType1Type of input right-hand side vector.
VectorType2Type of input approximate solution vector.
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
template<typename ValueType, typename MemorySpace>
template<typename MatrixType , typename VectorType1 , typename VectorType2 >
void cusp::relaxation::jacobi< ValueType, MemorySpace >::operator() ( const MatrixType &  A,
const VectorType1 &  b,
VectorType2 &  x,
const ValueType  omega 
)

Perform Jacobi relaxation using specified omega damping factor.

Template Parameters
MatrixTypeType of input matrix.
VectorType1Type of input right-hand side vector.
VectorType2Type of input approximate solution vector.
Parameters
Amatrix of the linear system
xapproximate solution of the linear system
bright-hand side of the linear system
omegaDamping factor used in Jacobi smoother.