ug4
|
describes an inverse linear mapping X->Y More...
#include <linear_operator_inverse.h>
Public Types | |
typedef Y | codomain_function_type |
Range space. | |
typedef X | domain_function_type |
Domain space. | |
![]() | |
typedef Y | codomain_function_type |
Range space. | |
typedef X | domain_function_type |
Domain space. | |
Public Member Functions | |
virtual bool | apply (Y &u, const X &f)=0 |
applies inverse operator, i.e. returns u = A^{-1} f | |
virtual bool | apply_return_defect (Y &u, X &f)=0 |
applies inverse operator, i.e. returns u = A^{-1} f and returns defect d := f - A*u | |
virtual bool | apply_update_defect (Y &u, X &f) |
compute new correction c = B*d and update defect d := d - A*c | |
virtual SmartPtr< ILinearIterator< X, Y > > | clone () |
clone | |
virtual std::string | config_string () const |
returns information about configuration parameters | |
SmartPtr< IConvergenceCheck< X > > | convergence_check () |
returns the convergence check | |
ConstSmartPtr< IConvergenceCheck< X > > | convergence_check () const |
returns the convergence check | |
number | defect () const |
returns the current defect | |
ILinearOperatorInverse () | |
constructor setting convergence check to (100, 1e-12, 1e-12, true) | |
ILinearOperatorInverse (SmartPtr< IConvergenceCheck< X > > spConvCheck) | |
Default constructor. | |
virtual bool | init (SmartPtr< ILinearOperator< Y, X > > J, const Y &u) |
initializes for the inverse for a linearized operator at linearization point u | |
virtual bool | init (SmartPtr< ILinearOperator< Y, X > > L) |
initializes for the inverse for a linear operator | |
SmartPtr< ILinearOperator< Y, X > > | linear_operator () |
returns the current Operator this Inverse Operator is initialized for | |
virtual const char * | name () const =0 |
returns the name of the operator inverse | |
number | reduction () const |
returns the current relative reduction | |
void | set_convergence_check (SmartPtr< IConvergenceCheck< X > > spConvCheck) |
set the convergence check | |
virtual int | standard_offset () const |
returns the standard offset for output | |
int | step () const |
returns the current number of steps | |
virtual bool | supports_parallel () const =0 |
returns if parallel solving is supported | |
virtual | ~ILinearOperatorInverse () |
virtual destructor | |
![]() | |
SmartPtr< IDamping< X, Y > > | damping () |
returns the scaling | |
ILinearIterator () | |
constructor | |
ILinearIterator (const ILinearIterator< X, Y > &parent) | |
copy constructor | |
void | set_damp (number factor) |
sets the damping to a constant factor | |
void | set_damp (SmartPtr< IDamping< X, Y > > spScaling) |
sets a scaling for the correction | |
virtual | ~ILinearIterator () |
virtual destructor | |
Protected Attributes | |
SmartPtr< IConvergenceCheck< X > > | m_spConvCheck |
smart pointer holding the convergence check | |
SmartPtr< ILinearOperator< Y, X > > | m_spLinearOperator |
Operator that is inverted by this Inverse Operator. | |
![]() | |
SmartPtr< IDamping< X, Y > > | m_spDamping |
the scaling | |
describes an inverse linear mapping X->Y
This class is the base class for the inversion of linear operator given in form of the ILinearOperator interface class. Given a operator L, the basic usage of this class is to invert this operator, i.e. to compute the solution u of
L*u = f i.e. u := L^{-1} f
This application has been split up into three steps:
This splitting has been made, since initialization and preparation may be computationally expansive. Thus, the user of this class has the choice when to call this initialization/preparation. E.g. when the operator is applied several times on the same vectors, those have only to be prepared once and the init of the operator is only needed once.
X | domain space |
Y | range space |
typedef Y ug::ILinearOperatorInverse< X, Y >::codomain_function_type |
Range space.
typedef X ug::ILinearOperatorInverse< X, Y >::domain_function_type |
Domain space.
|
inline |
constructor setting convergence check to (100, 1e-12, 1e-12, true)
References new().
|
inline |
Default constructor.
|
inlinevirtual |
virtual destructor
|
pure virtual |
applies inverse operator, i.e. returns u = A^{-1} f
This method applies the inverse operator, i.e. u = A^{-1} f. The domain function f remains unchanged. Note, that this method can always be implemented by creating a copy of f and calling apply_return_defect with this copy.
[in] | f | right-hand side |
[out] | u | solution |
Implements ug::ILinearIterator< X, Y >.
Implemented in ug::PrimalSubassembledMatrixInverse< TAlgebra >, ug::PrimalSubassembledMatrixInverse< algebra_type >, ug::AnalyzingSolver< M, X, Y >, and ug::IMatrixOperatorInverse< M, X, Y >.
Referenced by ug::OperatorInverseIterator< TAlgebra >::apply(), ug::SchurInverseWithOperator< TAlgebra >::apply(), ug::SchurInverseWithFullMatrix< TAlgebra >::apply(), ug::SchurInverseWithAutoFullMatrix< TAlgebra >::apply(), ug::PowerMethod< TAlgebra >::calculate_max_eigenvalue(), and ug::PowerMethod< TAlgebra >::calculate_min_eigenvalue().
|
pure virtual |
applies inverse operator, i.e. returns u = A^{-1} f and returns defect d := f - A*u
This method applies the inverse operator, i.e. u = A^{-1} f. The domain function f is changed in the way, that the defect d := f - A*u is returned in the function. This is always useful, when the inverting algorithm can (or must) update the defect during computation (this is e.g. the case for the geometric multigrid method). Note, that this method can always be implemented by calling apply and then computing d := f - A*u.
[in,out] | f | right-hand side |
[out] | u | solution |
Implemented in ug::PrimalSubassembledMatrixInverse< TAlgebra >, ug::PrimalSubassembledMatrixInverse< algebra_type >, ug::AnalyzingSolver< M, X, Y >, and ug::IMatrixOperatorInverse< M, X, Y >.
Referenced by ug::IPreconditionedLinearOperatorInverse< X >::apply(), ug::SchurInverseWithOperator< TAlgebra >::apply_return_defect(), ug::SchurInverseWithFullMatrix< TAlgebra >::apply_return_defect(), ug::SchurInverseWithAutoFullMatrix< TAlgebra >::apply_return_defect(), ug::OperatorInverseIterator< TAlgebra >::apply_update_defect(), and ug::ILinearOperatorInverse< X, Y >::apply_update_defect().
|
inlinevirtual |
compute new correction c = B*d and update defect d := d - A*c
This method applies the inverse operator, i.e. c = B*d. The domain function d is changed in the way, that the defect d := d - A*c is returned in the function. This is always useful, when the iterating algorithm can (or must) update the defect during computation (this is e.g. the case for the geometric multigrid method). Note, that this method can always be implemented by calling apply and then computing d := d - A*c.
[in,out] | d | defect |
[out] | u | correction |
Implements ug::ILinearIterator< X, Y >.
References ug::ILinearOperatorInverse< X, Y >::apply_return_defect().
Referenced by ug::LinearSolver< TVector >::compute_correction().
|
inlinevirtual |
|
inlinevirtual |
returns information about configuration parameters
this should return necessary information about parameters and possibly calling config_string of subcomponents.
Reimplemented from ug::ILinearIterator< X, Y >.
Reimplemented in ug::IPreconditionedLinearOperatorInverse< X >, ug::IPreconditionedLinearOperatorInverse< TVector >, ug::IPreconditionedLinearOperatorInverse< ug::ParallelVector >, ug::IPreconditionedLinearOperatorInverse< vector_type >, ug::AgglomeratingSolver< TAlgebra >, ug::AnalyzingSolver< M, X, Y >, ug::BiCGStab< TVector >, ug::GMRES< TVector >, and ug::LU< TAlgebra >.
References ug::ILinearOperatorInverse< X, Y >::name().
Referenced by ug::SchurInverseWithOperator< TAlgebra >::config_string(), ug::SchurInverseWithFullMatrix< TAlgebra >::config_string(), ug::SchurInverseWithAutoFullMatrix< TAlgebra >::config_string(), and ug::SchurPrecond< TAlgebra >::config_string().
|
inline |
returns the convergence check
References ug::ILinearOperatorInverse< X, Y >::m_spConvCheck.
|
inline |
returns the convergence check
References ug::ILinearOperatorInverse< X, Y >::m_spConvCheck.
Referenced by ug::LU< TAlgebra >::apply(), ug::AutoLinearSolver< TVector >::apply_return_defect(), ug::BiCGStab< TVector >::apply_return_defect(), ug::CG< TVector >::apply_return_defect(), ug::GMRES< TVector >::apply_return_defect(), ug::LinearSolver< TVector >::apply_return_defect(), ug::ILinearOperatorInverse< X, Y >::defect(), ug::AutoLinearSolver< TVector >::prepare_conv_check(), ug::BiCGStab< TVector >::prepare_conv_check(), ug::CG< TVector >::prepare_conv_check(), ug::FETISolver< TAlgebra >::prepare_conv_check(), ug::GMRES< TVector >::prepare_conv_check(), ug::LinearSolver< TVector >::prepare_conv_check(), ug::ILinearOperatorInverse< X, Y >::reduction(), and ug::ILinearOperatorInverse< X, Y >::step().
|
inline |
returns the current defect
References ug::ILinearOperatorInverse< X, Y >::convergence_check().
Referenced by ug::AutoLinearSolver< TVector >::apply_return_defect().
|
inlinevirtual |
initializes for the inverse for a linearized operator at linearization point u
This method passes the linear operator J(u) that should be inverted by this operator. As second argument the linearization point is passed. This is needed e.g. for the geometric multigrid method, that inverts a linearized operator based on coarser grid operators, that have to be initialized based on the linearization point.
[in] | J | linearized operator to invert |
[in] | u | linearization point |
Implements ug::ILinearIterator< X, Y >.
Reimplemented in ug::IMatrixOperatorInverse< M, X, Y >, ug::AnalyzingSolver< M, X, Y >, and ug::LU< TAlgebra >.
References ug::ILinearOperatorInverse< X, Y >::m_spLinearOperator.
|
inlinevirtual |
initializes for the inverse for a linear operator
This method passes the operator L that is inverted by this operator. In addition some preparation step can be made.
[in] | L | linear operator to invert |
Implements ug::ILinearIterator< X, Y >.
Reimplemented in ug::IMatrixOperatorInverse< M, X, Y >, ug::AnalyzingSolver< M, X, Y >, and ug::LU< TAlgebra >.
References ug::ILinearOperatorInverse< X, Y >::m_spLinearOperator.
Referenced by ug::PowerMethod< TAlgebra >::calculate_max_eigenvalue(), ug::PowerMethod< TAlgebra >::calculate_min_eigenvalue(), ug::OperatorInverseIterator< TAlgebra >::init(), ug::OperatorInverseIterator< TAlgebra >::init(), ug::SchurInverseWithOperator< TAlgebra >::init(), ug::SchurInverseWithFullMatrix< TAlgebra >::init(), ug::SchurInverseWithAutoFullMatrix< TAlgebra >::init(), and ug::AutoLinearSolver< TVector >::init_op().
|
inline |
returns the current Operator this Inverse Operator is initialized for
References ug::ILinearOperatorInverse< X, Y >::m_spLinearOperator, ug::ILinearOperatorInverse< X, Y >::name(), and UG_THROW.
Referenced by ug::IPreconditionedLinearOperatorInverse< X >::apply(), ug::AutoLinearSolver< TVector >::apply_return_defect(), ug::BiCGStab< TVector >::apply_return_defect(), ug::CG< TVector >::apply_return_defect(), ug::GMRES< TVector >::apply_return_defect(), ug::LinearSolver< TVector >::apply_return_defect(), and ug::AutoLinearSolver< TVector >::compute_correction().
|
pure virtual |
returns the name of the operator inverse
This method returns the name of the inverse operator. This function is typically needed, when the inverse operator is used inside of another and some debug output should be printed
Implements ug::ILinearIterator< X, Y >.
Implemented in ug::AgglomeratingSolver< TAlgebra >, ug::AnalyzingSolver< M, X, Y >, ug::AutoLinearSolver< TVector >, ug::BiCGStab< TVector >, ug::CG< TVector >, ug::IExternalSolver< TAlgebra >, ug::PrimalSubassembledMatrixInverse< TAlgebra >, ug::PrimalSubassembledMatrixInverse< algebra_type >, ug::FETISolver< TAlgebra >, ug::GMRES< TVector >, ug::LinearSolver< TVector >, ug::LinearSolver< ug::ParallelVector >, and ug::LU< TAlgebra >.
Referenced by ug::IPreconditionedLinearOperatorInverse< X >::apply(), ug::ILinearOperatorInverse< X, Y >::config_string(), ug::IPreconditionedLinearOperatorInverse< X >::config_string(), ug::IPreconditionedLinearOperatorInverse< X >::init(), ug::IPreconditionedLinearOperatorInverse< X >::init(), ug::ILinearOperatorInverse< X, Y >::linear_operator(), and ug::OperatorInverseIterator< TAlgebra >::OperatorInverseIterator().
|
inline |
returns the current relative reduction
References ug::ILinearOperatorInverse< X, Y >::convergence_check().
Referenced by ug::AutoLinearSolver< TVector >::apply_return_defect().
|
inline |
set the convergence check
References ug::ILinearOperatorInverse< X, Y >::m_spConvCheck, and ug::ILinearOperatorInverse< X, Y >::standard_offset().
|
inlinevirtual |
returns the standard offset for output
Referenced by ug::ILinearOperatorInverse< X, Y >::set_convergence_check().
|
inline |
returns the current number of steps
References ug::ILinearOperatorInverse< X, Y >::convergence_check().
Referenced by ug::BiCGStab< TVector >::apply_return_defect(), and ug::CG< TVector >::apply_return_defect().
|
pure virtual |
returns if parallel solving is supported
Implements ug::ILinearIterator< X, Y >.
Implemented in ug::AgglomeratingSolver< TAlgebra >, ug::AnalyzingSolver< M, X, Y >, ug::AutoLinearSolver< TVector >, ug::BiCGStab< TVector >, ug::CG< TVector >, ug::IExternalSolver< TAlgebra >, ug::PrimalSubassembledMatrixInverse< TAlgebra >, ug::PrimalSubassembledMatrixInverse< algebra_type >, ug::FETISolver< TAlgebra >, ug::GMRES< TVector >, ug::LinearSolver< TVector >, ug::LinearSolver< ug::ParallelVector >, and ug::LU< TAlgebra >.
Referenced by ug::PrimalSubassembledMatrixInverse< TAlgebra >::supports_parallel(), ug::FETISolver< TAlgebra >::supports_parallel(), ug::OperatorInverseIterator< TAlgebra >::supports_parallel(), ug::SchurInverseWithOperator< TAlgebra >::supports_parallel(), ug::SchurInverseWithFullMatrix< TAlgebra >::supports_parallel(), ug::SchurInverseWithAutoFullMatrix< TAlgebra >::supports_parallel(), and ug::SchurPrecond< TAlgebra >::supports_parallel().
|
protected |
smart pointer holding the convergence check
Referenced by ug::IPreconditionedLinearOperatorInverse< X >::config_string_preconditioner_convergence_check(), ug::ILinearOperatorInverse< X, Y >::convergence_check(), ug::ILinearOperatorInverse< X, Y >::convergence_check(), and ug::ILinearOperatorInverse< X, Y >::set_convergence_check().
|
protected |
Operator that is inverted by this Inverse Operator.
Referenced by ug::ILinearOperatorInverse< X, Y >::init(), ug::ILinearOperatorInverse< X, Y >::init(), and ug::ILinearOperatorInverse< X, Y >::linear_operator().