ug4
|
Jacobi-Iteration. More...
#include <jacobi.h>
Public Types | |
typedef TAlgebra | algebra_type |
Algebra type. | |
typedef IPreconditioner< TAlgebra > | base_type |
Base type. | |
typedef IPreconditioner< TAlgebra >::matrix_operator_type | matrix_operator_type |
Matrix Operator type. | |
typedef TAlgebra::matrix_type | matrix_type |
Matrix type. | |
typedef TAlgebra::vector_type | vector_type |
Vector type. | |
![]() | |
typedef TAlgebra | algebra_type |
Algebra type. | |
typedef MatrixOperator< matrix_type, vector_type > | matrix_operator_type |
Matrix Operator type. | |
typedef TAlgebra::matrix_type | matrix_type |
Matrix type. | |
typedef TAlgebra::vector_type | vector_type |
Vector type. | |
![]() | |
typedef TAlgebra::vector_type | codomain_function_type |
Range space. | |
typedef TAlgebra::vector_type | domain_function_type |
Domain space. | |
![]() | |
typedef TAlgebra | algebra_type |
type of algebra | |
typedef TAlgebra::matrix_type | matrix_type |
type of matrix | |
typedef TAlgebra::vector_type | vector_type |
type of vector | |
![]() | |
typedef TAlgebra::vector_type | vector_type |
type of vector | |
Public Member Functions | |
virtual SmartPtr< ILinearIterator< vector_type > > | clone () |
Clone. | |
Jacobi () | |
default constructor | |
Jacobi (const Jacobi< TAlgebra > &parent) | |
clone constructor | |
Jacobi (number damp) | |
constructor setting the damping parameter | |
void | set_block (bool b) |
sets if blocked jacobi is used (inverting block-diagonal), or plain (scalar) diagonal if false | |
virtual bool | supports_parallel () const |
returns if parallel solving is supported | |
virtual | ~Jacobi () |
Destructor. | |
![]() | |
virtual bool | apply_update_defect (vector_type &c, vector_type &d) |
compute new correction c = B*d and update defect d:= d - L*c | |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | approx_operator () |
underlying matrix based operator used for the preconditioner | |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | defect_operator () |
underlying matrix based operator for calculation of defect | |
virtual bool | init (SmartPtr< ILinearOperator< vector_type > > J, const vector_type &u) |
implements the ILinearIterator-interface for matrix based preconditioner | |
bool | init (SmartPtr< ILinearOperator< vector_type > > L) |
implements the ILinearIterator-interface for matrix based preconditioner | |
bool | init (SmartPtr< MatrixOperator< matrix_type, vector_type > > Op) |
initializes the preconditioner for a matrix based operator | |
IPreconditioner () | |
default constructor | |
IPreconditioner (const IPreconditioner< TAlgebra > &parent) | |
clone constructor | |
IPreconditioner (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter) | |
constructor setting debug writer | |
virtual void | set_approximation (SmartPtr< MatrixOperator< matrix_type, vector_type > > approx) |
virtual | ~IPreconditioner () |
virtual destructor | |
![]() | |
virtual std::string | config_string () const |
SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > | damping () |
returns the scaling | |
ILinearIterator () | |
constructor | |
ILinearIterator (const ILinearIterator< TAlgebra::vector_type, TAlgebra::vector_type > &parent) | |
copy constructor | |
virtual bool | init (SmartPtr< ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type > > J, const TAlgebra::vector_type &u)=0 |
initialize for operator J(u) and linearization point u | |
virtual bool | init (SmartPtr< ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type > > L)=0 |
initialize for linear operator L | |
void | set_damp (number factor) |
sets the damping to a constant factor | |
void | set_damp (SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > spScaling) |
sets a scaling for the correction | |
virtual | ~ILinearIterator () |
virtual destructor | |
![]() | |
SmartPtr< IDebugWriter< algebra_type > > | debug_writer () |
returns the debug writer | |
ConstSmartPtr< IDebugWriter< algebra_type > > | debug_writer () const |
bool | debug_writer_valid () const |
returns true if the debug writer is set | |
DebugWritingObject () | |
DebugWritingObject (const DebugWritingObject< algebra_type > &parent) | |
clone constructor | |
DebugWritingObject (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter) | |
virtual void | set_debug (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter) |
set debug writer | |
virtual | ~DebugWritingObject () |
virtual destructor | |
![]() | |
virtual void | set_debug (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter) |
set debug writer | |
SmartPtr< IVectorDebugWriter< vector_type > > | vector_debug_writer () |
returns the debug writer | |
ConstSmartPtr< IVectorDebugWriter< vector_type > > | vector_debug_writer () const |
bool | vector_debug_writer_valid () const |
returns true if the debug writer is set | |
VectorDebugWritingObject () | |
VectorDebugWritingObject (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter) | |
void | write_debug (const vector_type &vec, const char *filename) |
writing debug output for a vector (if debug writer set) | |
virtual | ~VectorDebugWritingObject () |
virtual destructor | |
Protected Types | |
typedef block_traits< typenamematrix_type::value_type >::inverse_type | inverse_type |
type of block-inverse | |
Protected Member Functions | |
virtual bool | apply (vector_type &c, const vector_type &d) |
compute new correction c = B*d | |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | approx_operator () |
underlying matrix based operator used for the preconditioner | |
virtual const char * | name () const |
Name of preconditioner. | |
virtual bool | postprocess () |
Postprocess routine. | |
virtual bool | preprocess (SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp) |
Preprocess routine. | |
virtual bool | step (SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp, vector_type &c, const vector_type &d) |
computes a new correction c = B*d | |
![]() | |
void | enter_debug_writer_section (const char *secDir) |
enters a debugging section | |
void | enter_debug_writer_section (std::string secDir) |
enters a debugging section | |
void | leave_debug_writer_section () |
leaves a debugging section | |
void | write_debug (const matrix_type &mat, const char *filename) |
write debug output for a matrix (if debug writer set) | |
void | write_debug (const matrix_type &mat, std::string name) |
write debug output for a matrix (if debug writer set) | |
![]() | |
void | enter_vector_debug_writer_section (const char *secDir) |
enters a debugging section | |
void | enter_vector_debug_writer_section (std::string secDir) |
enters a debugging section | |
void | leave_vector_debug_writer_section () |
leaves a debugging section | |
void | print_debugger_message (const char *msg) |
prints a debugger message (listing all the sections) | |
void | print_debugger_message (std::string msg) |
prints a debugger message (listing all the sections) | |
virtual void | write_debug (const vector_type &vec, std::string name) |
writing debug output for a vector (if debug writer set) | |
Protected Attributes | |
bool | m_bBlock |
std::vector< inverse_type > | m_diagInv |
storage of the inverse diagonal in parallel | |
![]() | |
bool | m_bInit |
init flag indicating if init has been called | |
bool | m_bOtherApproxOperator |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | m_spApproxOperator |
underlying matrix based operator used for the preconditioner | |
SmartPtr< ILinearOperator< vector_type > > | m_spDefectOperator |
underlying matrix based operator for calculation of defect | |
![]() | |
SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > | m_spDamping |
the scaling | |
![]() | |
SmartPtr< IDebugWriter< algebra_type > > | m_spDebugWriter |
Debug Writer. | |
![]() | |
SmartPtr< IVectorDebugWriter< vector_type > > | m_spVectorDebugWriter |
Debug Writer. | |
Jacobi-Iteration.
Here, the Jacobi-iteration is described for solving the linear equation
\f$ A * x = b. A \in R^{nxn}, x \in R^n, b \in R^n \f$.
Most of the common linear iteration-methods base on the decomposition of A into its diagonal (D) and strict-upper(-U) and strict-lower part (-L),
\( A = D - L - U \).
Among others, W. Hackbusch ('Iterative Loesung grosser Gleichungssysteme'), distinguishes three different forms for describing a linear iteration scheme. The general 'first normal-form' of a linear iteration scheme takes the form
\( x^{m+1} = M * x^m + N * b \),
with some Matrices \( M \) and \( N \in R^{nxn} \). m denotes the iteration index. The general 'second normal-form' of a linear iteration scheme takes the form
\( x^{m+1} = x^m - N * (A * x^m - b) \).
Those linear iteration schemes, which can be represented by the second normal-form are the linear, consistent iteration schemes.
Introducing the correction \( c{m+1} := x^{m+1} - x^m \) and the defect \( d^m := b - A * x^m \) the second normal-form can be rewritten as
\( c = N * d \).
The matrix of the second normal-form for the Jacobi-method takes the simple form
\( N = D^{-1} \). .
References:
Jacobi Preconditioner
typedef TAlgebra ug::Jacobi< TAlgebra >::algebra_type |
Algebra type.
typedef IPreconditioner<TAlgebra> ug::Jacobi< TAlgebra >::base_type |
Base type.
|
protected |
type of block-inverse
typedef IPreconditioner<TAlgebra>::matrix_operator_type ug::Jacobi< TAlgebra >::matrix_operator_type |
Matrix Operator type.
typedef TAlgebra::matrix_type ug::Jacobi< TAlgebra >::matrix_type |
Matrix type.
typedef TAlgebra::vector_type ug::Jacobi< TAlgebra >::vector_type |
Vector type.
|
inline |
default constructor
References ug::Jacobi< TAlgebra >::m_bBlock, and ug::ILinearIterator< TAlgebra::vector_type >::set_damp().
|
inline |
constructor setting the damping parameter
References ug::Jacobi< TAlgebra >::m_bBlock, and ug::ILinearIterator< TAlgebra::vector_type >::set_damp().
|
inline |
clone constructor
References ug::Jacobi< TAlgebra >::m_bBlock, and ug::Jacobi< TAlgebra >::set_block().
|
inlinevirtual |
Destructor.
|
inlineprotectedvirtual |
compute new correction c = B*d
This method implements the virtual method of the ILinearIterator-interface. Basically, besides some common checks the request is forwarded to the (virtual) 'step'-method.
[out] | c | correction |
[in] | d | defect |
Reimplemented from ug::IPreconditioner< TAlgebra >.
References ug::Jacobi< TAlgebra >::approx_operator(), ug::ILinearIterator< TAlgebra::vector_type >::damping(), ug::IPreconditioner< TAlgebra >::m_bInit, ug::Jacobi< TAlgebra >::name(), num_cols(), num_rows(), PROFILE_BEGIN_GROUP, ug::PST_ADDITIVE, ug::PST_CONSISTENT, ug::Jacobi< TAlgebra >::step(), THROW_IF_NOT_EQUAL_4, UG_LOG, and UG_THROW.
|
inlineprotected |
underlying matrix based operator used for the preconditioner
Referenced by ug::Jacobi< TAlgebra >::apply().
|
inlinevirtual |
|
inlineprotectedvirtual |
Name of preconditioner.
Implements ug::IPreconditioner< TAlgebra >.
Referenced by ug::Jacobi< TAlgebra >::apply().
|
inlineprotectedvirtual |
Postprocess routine.
Implements ug::IPreconditioner< TAlgebra >.
|
inlineprotectedvirtual |
Preprocess routine.
Implements ug::IPreconditioner< TAlgebra >.
References ug::ParallelVector< TVector >::change_storage_type(), ug::CheckVectorInvertible(), ug::ILinearIterator< TAlgebra::vector_type >::damping(), ug::GetDiag(), ug::GetInverse(), ug::Jacobi< TAlgebra >::m_bBlock, ug::Jacobi< TAlgebra >::m_diagInv, PROFILE_BEGIN_GROUP, ug::PST_ADDITIVE, ug::PST_CONSISTENT, ug::ParallelVector< TVector >::set_layouts(), ug::ParallelVector< TVector >::set_storage_type(), and UG_LOG.
|
inline |
sets if blocked jacobi is used (inverting block-diagonal), or plain (scalar) diagonal if false
References ug::Jacobi< TAlgebra >::m_bBlock.
Referenced by ug::Jacobi< TAlgebra >::Jacobi().
|
inlineprotectedvirtual |
computes a new correction c = B*d
This method computes a new correction c = B*d. It can only be called, when the preprocess has been done.
[in] | mat | underlying matrix (i.e. L in L*u = f) |
[out] | c | correction |
[in] | d | defect |
Implements ug::IPreconditioner< TAlgebra >.
References ug::Jacobi< TAlgebra >::m_diagInv, ug::MatMult(), PROFILE_BEGIN_GROUP, ug::PST_ADDITIVE, ug::PST_CONSISTENT, and UG_LOG.
Referenced by ug::Jacobi< TAlgebra >::apply().
|
inlinevirtual |
returns if parallel solving is supported
Implements ug::ILinearIterator< TAlgebra::vector_type >.
|
protected |
|
protected |
storage of the inverse diagonal in parallel
Referenced by ug::Jacobi< TAlgebra >::preprocess(), and ug::Jacobi< TAlgebra >::step().