FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
FEAT::Solver::IterativeSolver< Vector_ > Class Template Referenceabstract

Abstract base-class for iterative solvers. More...

#include <iterative.hpp>

Inheritance diagram for FEAT::Solver::IterativeSolver< Vector_ >:
FEAT::Solver::SolverBase< Vector_ > FEAT::Solver::PreconditionedIterativeSolver< Matrix_::VectorTypeR > FEAT::Solver::PreconditionedIterativeSolver< Functional_::VectorTypeR > FEAT::Solver::PreconditionedIterativeSolver< Vector_ > FEAT::Solver::BiCGStab< Matrix_, Filter_ > FEAT::Solver::BiCGStabL< Matrix_, Filter_ > FEAT::Solver::Descent< Matrix_, Filter_ > FEAT::Solver::FGMRES< Matrix_, Filter_ > FEAT::Solver::GMRES< Matrix_, Filter_ > FEAT::Solver::GroppPCG< Matrix_, Filter_ > FEAT::Solver::IDRS< Matrix_, Filter_ > FEAT::Solver::PCG< Matrix_, Filter_ > FEAT::Solver::PCR< Matrix_, Filter_ > FEAT::Solver::PMR< Matrix_, Filter_ > FEAT::Solver::PSD< Matrix_, Filter_ > FEAT::Solver::PipePCG< Matrix_, Filter_ > FEAT::Solver::RBiCGStab< Matrix_, Filter_ > FEAT::Solver::RGCR< Matrix_, Filter_ > FEAT::Solver::Richardson< Matrix_, Filter_ > FEAT::Solver::NLOptLS< Functional_, Filter_ >

Public Types

typedef SolverBase< VectorTypeBaseClass
 The base class. More...
 
typedef VectorType::DataType DataType
 Floating point type. More...
 
typedef Vector_ VectorType
 The vector type this solver can be applied to. More...
 

Public Member Functions

virtual Status apply (Vector_ &vec_cor, const Vector_ &vec_def)=0
 Solver application method. More...
 
DataType calc_convergence_rate () const
 Computes the overall convergence rate: (defect_final / defect_initial) ^ (1 / number of iterations) More...
 
DataType calc_defect_reduction () const
 Computes the overall defect reduction factor: (defect_final / defect_inital) More...
 
virtual Status correct (VectorType &vec_sol, const VectorType &vec_rhs)=0
 Solver correction method. More...
 
virtual void done ()
 Finalization method. More...
 
virtual void done_numeric ()
 Numeric finalization method. More...
 
virtual void done_symbolic ()
 Symbolic finalization method. More...
 
virtual void force_defect_norm_calc (bool force)
 Forces the calculation of defect norms in every iteration. More...
 
DataType get_def_final () const
 Returns the final defect. More...
 
DataType get_def_initial () const
 Returns the initial defect. More...
 
DataType get_div_abs () const
 Returns the absolute divergence. More...
 
DataType get_div_rel () const
 Returns the relative divergence. More...
 
Index get_max_iter () const
 Returns the maximum number of iterations. More...
 
Index get_min_iter () const
 Returns the minimal number of iterations. More...
 
Index get_min_stag_iter () const
 Returns the minimum stagnation iteration count. More...
 
Index get_num_iter () const
 Returns number of performed iterations. More...
 
String get_plot_name () const
 Returns the plot name of the solver. More...
 
DataType get_stag_rate () const
 Returns the stagnation rate. More...
 
Status get_status () const
 Returns the status. More...
 
virtual String get_summary () const
 Returns a summary string. More...
 
DataType get_tol_abs () const
 Returns the absolute tolerance. More...
 
DataType get_tol_abs_low () const
 Returns the lower absolute tolerance. More...
 
DataType get_tol_rel () const
 Returns the relative tolerance. More...
 
virtual void init ()
 Initialization method. More...
 
virtual void init_numeric ()
 Numeric initialization method. More...
 
virtual void init_symbolic ()
 Symbolic initialization method. More...
 
bool is_converged () const
 checks for convergence More...
 
bool is_converged (const DataType def_cur) const
 checks for convergence More...
 
bool is_diverged () const
 checks for divergence More...
 
bool is_diverged (const DataType def_cur) const
 checks for divergence More...
 
virtual String name () const =0
 Returns a descriptive string. More...
 
virtual void plot_summary () const
 Plot a summary of the last solver run. More...
 
void set_div_abs (DataType div_abs)
 Sets the absolute divergence for the solver. More...
 
void set_div_rel (DataType div_rel)
 Sets the relative divergence for the solver. More...
 
void set_max_iter (Index max_iter)
 Sets the maximum iteration count for the solver. More...
 
void set_min_iter (Index min_iter)
 Sets the minimum iteration count for the solver. More...
 
void set_min_stag_iter (Index min_iter)
 Sets the minimum stagnate iteration count for the solver. More...
 
void set_plot_interval (const Index plot_interval)
 Sets the interval between two plot outputs, if any. More...
 
void set_plot_mode (const PlotMode plot_mode)
 Sets the plot mode of the solver. More...
 
void set_plot_name (const String &plot_name)
 Sets the plot name of the solver. More...
 
void set_stag_rate (DataType rate)
 Sets the stagnation rate for the solver. More...
 
void set_tol_abs (DataType tol_abs)
 Sets the absolute tolerance for the solver. More...
 
void set_tol_abs_low (DataType tol_abs_low)
 Sets the lower absolute tolerance for the solver. More...
 
void set_tol_rel (DataType tol_rel)
 Sets the relative tolerance for the solver. More...
 

Protected Member Functions

 IterativeSolver (const String &plot_name)
 Protected constructor. More...
 
 IterativeSolver (const String &plot_name, const String &section_name, const PropertyMap *section)
 Constructor using a PropertyMap. More...
 
virtual Status _analyze_defect (Index num_iter, DataType def_cur, DataType def_prev, bool check_stag)
 Internal function: analyze the current defect. More...
 
virtual DataType _calc_def_norm (const VectorType &vec_def, const VectorType &vec_sol)
 Computes the defect norm. More...
 
bool _plot_iter (Status st=Status::progress) const
 Plot the current iteration? More...
 
virtual void _plot_iter_line (Index num_iter, DataType def_cur, DataType def_prev)
 Plots an iteration line. More...
 
bool _plot_summary () const
 Plot summary? More...
 
void _print_line (const String &line) const
 Prints a line. More...
 
bool _progress () const
 Progress iteration? More...
 
void _set_comm (const Dist::Comm *comm)
 Sets the communicator for the solver directly. More...
 
template<typename Matrix_ >
void _set_comm_by_matrix (const Matrix_ &matrix)
 Sets the communicator for the solver from a matrix. More...
 
void _set_comm_by_vector (const Vector_ &vector)
 Sets the communicator for the solver from a vector. More...
 
virtual Status _set_initial_defect (const VectorType &vec_def, const VectorType &vec_sol)
 Internal function: sets the initial defect vector. More...
 
virtual Status _set_new_defect (const VectorType &vec_def, const VectorType &vec_sol)
 Internal function: sets the new (next) defect vector. More...
 
virtual Status _update_defect (const DataType def_cur_norm)
 Internal function: sets the new (next) defect norm. More...
 

Protected Attributes

const Dist::Comm_comm
 Communicator of the solver. More...
 
DataType _def_cur
 current defect More...
 
DataType _def_init
 initial defect More...
 
DataType _def_prev
 previous iteration defect More...
 
DataType _div_abs
 absolute divergence parameter More...
 
DataType _div_rel
 relative divergence parameter More...
 
bool _force_def_norm_calc
 whether skipping of defect computation is allowed or not More...
 
Index _iter_digits
 iteration count digits for plotting More...
 
Index _max_iter
 maximum number of iterations More...
 
Index _min_iter
 minimum number of iterations More...
 
Index _min_stag_iter
 minimum number of stagnation iterations More...
 
Index _num_iter
 number of performed iterations More...
 
Index _num_stag_iter
 number of consecutive stagnated iterations More...
 
Index _plot_interval
 plot output interval More...
 
PlotMode _plot_mode
 whether to plot something More...
 
String _plot_name
 name of the solver in plots More...
 
DataType _stag_rate
 stagnation rate More...
 
Status _status
 current status of the solver More...
 
DataType _tol_abs
 absolute tolerance parameter More...
 
DataType _tol_abs_low
 absolute low tolerance parameter More...
 
DataType _tol_rel
 relative tolerance parameter More...
 

Detailed Description

template<typename Vector_>
class FEAT::Solver::IterativeSolver< Vector_ >

Abstract base-class for iterative solvers.

This class template acts as an abstract base class for iterative solvers and it implements various auxiliary features for convergence control, output formatting and a posteriori convergence analysis.

Template Parameters
Vector_The class of the vector that is passed to the solver in the solve() method.

Details about the convergence control mechanisms offered by this class:
This class offers a wide variety of convergence control mechanisms and stopping criterions, which are summarized in this section.

In the following, let \(k\geq 0\) denote the current iteration and let \(r_0 := b - Ax_0\) denote the initial defect/residual vector and let \(r_k := b - Ax_k\) denote the defect/residual vector of the current iteration.

1. Divergence Criterion: The first stopping criterion checks whether the residual vector norm has risen above a relative or absolute divergence tolerance and, if so, the solver stops and returns the status code Status::diverged.
The divergence criterion is defined as follows:

\[ (\vert \vert r_k\vert\vert \geq \text{div\_abs}) \lor (\vert \vert r_k\vert\vert \geq \text{div\_rel}\cdot\vert \vert r_0\vert\vert) \]

By default, the divergence tolerances are chosen as

  • div_abs = 1/eps^2
  • dil_rel = 1/eps

where eps is the machine precision constant for the solvers underlying data type, which evaluates to approximately div_abs = 1E+32 and div_rel = 1E+16 for double precision.

2. Maximum Iteration Number Criterion: The second stopping criterion checks whether the solver has reached the maximum number of allowed iterations, i.e. whether k >= max_iter, and, if so, the solver returns the status code Status::max_iter.

3. Convergence Criterion: The third stopping criterion checks whether the residual vector has fallen below a relative and/or absolute convergence tolerance and, if so, the solver returns the status code Status::converged. The convergence criterion consists of a relative tolerance and two absolute tolerances and is defined as follows:

\[ (\vert \vert r_k\vert\vert \leq \text{tol\_abs}) \land ((\vert \vert r_k\vert\vert \leq \text{tol\_rel}\cdot\vert \vert r_0\vert\vert) \lor (\vert \vert r_k\vert\vert \leq \text{tol\_abs\_low})) \]

By default, the convergence tolerances are chosen as

  • tol_rel = sqrt(eps)
  • tol_abs = 1/eps^2
  • tol_abs_low = 0

where eps is the machine precision constant for the solvers underlying data type, which evaluates to approximately tol_abs = 1E-8 and tol_rel = 1E+32 for double precision, which effectively "disables" both absolute tolerances therefore resulting in a convergence criterion that only depends on the relative tolerance.

The tol_abs_low parameter can be used to "disarm" an overly strict relative tolerance tol_rel in the case where the initial defect is already quite small and therefore fulfilling the relative tolerance may become impossible due to finite precision arithmetic combined with a large condition number. A typical use case would be to set tol_abs to 1E-8 and tol_rel to 1E-5 (thus asking that the solver should gain 5 digits) and tol_abs_low to 1E-12, which will "disable" the tol_rel criterion if (and only if) the initial defect norm is less than 1E-7. If the initial defect norm is larger than 1E-7, then the tol_abs_low parameter will become irrelevant and the tol_abs and tol_rel criterions have to be fulfilled for convergence.

  • If you want the solver to fulfill only a relative tolerance, then set tol_rel to the corresponding tolerance and leave tol_abs_low set to 0 and tol_abs to set a large value (1/eps^2 by default).
  • If you want the solver to fulfill only an absolute tolerance, then set both tol_abs and tol_abs_low to the corresponding tolerance and set tol_rel to any value you like, e.g. 0.
  • If you want the solver to fulfill either a relative or an absolute tolerance, then set tol_rel and tol_abs_low to the corresponding tolerances and leave tol_abs to set a large value (1/eps^2 by default).
  • If you want the solver to fulfill both a relative and an absolute tolerance, then set tol_rel and tol_abs to the corresponding tolerances and leave tol_abs_low set to 0.

5. Stagnation Criterion: The fifth stopping criterion checks whether the residual vector norms do not decrease by a minimum factor over the last few iterations, thus indicating that the solver seems to be stagnating and that further iterations will probably not improve the solution any further and, if so, the solver stops and returns the status code Status::stagnated.
The stagnation criterion is fulfilled if min_stag_iter is set to a value > 0 and if

\[ \vert \vert r_k\vert\vert \geq \text{stag\_rate}\cdot\vert \vert r_{k-1}\vert\vert \]

for at least min_stag_iter consecutive iterations.

By default, the stagnation criterion is disabled and its parameters are chosen as

  • stag_rate = 0.95
  • min_stag_iter = 0

It is important to keep in mind that the stagnation criterion can usually only be helpful for simple fix-point type iterative solvers, which have more or less monotonous convergence behavior. Other types of solvers, especially Krylov subspace solvers, often show a somewhat "ballistic" convergence curve, where the defect norms seem to indicate divergence for the first few dozen/hundred iterations before the defect norms start falling down very quickly – in this case, the stagnation criterion will trigger a false positive during the initial climb of the defect norms. Also, some solvers can run into a stagnation scenario where the defect norms oscillate, which will not be detected by the stagnation criterion because the defect norm drops every other iteration before rising back up again in the next iteration.

Author
Peter Zajac

Definition at line 196 of file iterative.hpp.

Member Typedef Documentation

◆ BaseClass

template<typename Vector_ >
typedef SolverBase<VectorType> FEAT::Solver::IterativeSolver< Vector_ >::BaseClass

The base class.

Definition at line 205 of file iterative.hpp.

◆ DataType

template<typename Vector_ >
typedef VectorType::DataType FEAT::Solver::IterativeSolver< Vector_ >::DataType

Floating point type.

Definition at line 203 of file iterative.hpp.

◆ VectorType

template<typename Vector_ >
typedef Vector_ FEAT::Solver::IterativeSolver< Vector_ >::VectorType

The vector type this solver can be applied to.

Definition at line 201 of file iterative.hpp.

Constructor & Destructor Documentation

◆ IterativeSolver() [1/2]

template<typename Vector_ >
FEAT::Solver::IterativeSolver< Vector_ >::IterativeSolver ( const String plot_name)
inlineexplicitprotected

Protected constructor.

This constructor initializes the following values:

  • relative tolerance: sqrt(eps) (~1E-8 for double)
  • absolute tolerance: 1/eps^2 (~1E+32 for double)
  • lower absolute tolerance: 0
  • relative divergence: 1/eps (~1E+16 for double)
  • absolute divergence: 1/eps^2 (~1E+32 for double)
  • stagnation rate: 0.95
  • minimum iterations: 0
  • maximum iterations: 100
  • minimum stagnation iterations: 0
  • convergence plot: false
Parameters
[in]plot_nameSpecifies the name of the iterative solver. This is used as a prefix for the convergence plot.

Definition at line 270 of file iterative.hpp.

◆ IterativeSolver() [2/2]

template<typename Vector_ >
FEAT::Solver::IterativeSolver< Vector_ >::IterativeSolver ( const String plot_name,
const String section_name,
const PropertyMap section 
)
inlineexplicitprotected

Constructor using a PropertyMap.

See also
Solver configuration via PropertyMaps
Parameters
[in]section_nameThe name of the config section, which it does not know by itself
[in]sectionA pointer to the PropertyMap section configuring this solver
[in]plot_nameThe name of the solver.

Definition at line 311 of file iterative.hpp.

References FEAT::PropertyMap::get_entry().

Member Function Documentation

◆ _analyze_defect()

template<typename Vector_ >
virtual Status FEAT::Solver::IterativeSolver< Vector_ >::_analyze_defect ( Index  num_iter,
DataType  def_cur,
DataType  def_prev,
bool  check_stag 
)
inlineprotectedvirtual

Internal function: analyze the current defect.

Note
This function is called by _set_new_defect() and _update_defect().
Parameters
[in]num_iterCurrent number of iterations; usually = this->_num_iter.
[in]def_curCurrent defect norm; usually = this->_def_cur.
[in]def_prevPrevious defect norm; usually = this->_def_prev.
[in]check_stagSpecifies whether to check (and update) the stagnation criterion. This is typically set to false if one wants to check anything else than the 'true' next defect norm.
Returns
The updated status code.

Definition at line 866 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_min_iter, FEAT::Solver::aborted, FEAT::Solver::diverged, FEAT::Solver::IterativeSolver< Vector_ >::is_converged(), FEAT::Solver::IterativeSolver< Vector_ >::is_diverged(), FEAT::Math::isfinite(), FEAT::Solver::max_iter, FEAT::Solver::progress, FEAT::Solver::stagnated, and FEAT::Solver::success.

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::_set_new_defect(), and FEAT::Solver::IterativeSolver< Vector_ >::_update_defect().

◆ _calc_def_norm()

template<typename Vector_ >
virtual DataType FEAT::Solver::IterativeSolver< Vector_ >::_calc_def_norm ( const VectorType vec_def,
const VectorType vec_sol 
)
inlineprotectedvirtual

Computes the defect norm.

Parameters
[in]vec_defThe current defect vector.
[in]vec_solThe current solution vector approximation.

Definition at line 706 of file iterative.hpp.

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::_set_initial_defect(), FEAT::Solver::NLOptLS< Functional_, Filter_ >::_set_initial_defect(), FEAT::Solver::IterativeSolver< Vector_ >::_set_new_defect(), and FEAT::Solver::NLOptLS< Functional_, Filter_ >::_set_new_defect().

◆ _plot_iter()

◆ _plot_iter_line()

template<typename Vector_ >
virtual void FEAT::Solver::IterativeSolver< Vector_ >::_plot_iter_line ( Index  num_iter,
DataType  def_cur,
DataType  def_prev 
)
inlineprotectedvirtual

Plots an iteration line.

Parameters
[in]num_iterCurrent number of iterations; usually = this->_num_iter.
[in]def_curCurrent defect norm; usually = this->_def_cur.
[in]def_prevPrevious defect norm; usually = this->_def_prev.

Definition at line 773 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_print_line(), FEAT::String::pad_front(), FEAT::stringify(), FEAT::stringify_fp_fix(), and FEAT::stringify_fp_sci().

Referenced by FEAT::Solver::RBiCGStab< Matrix_, Filter_ >::_apply_intern(), FEAT::Solver::BiCGStab< Matrix_, Filter_ >::_apply_intern(), FEAT::Solver::IterativeSolver< Vector_ >::_set_initial_defect(), FEAT::Solver::IterativeSolver< Vector_ >::_set_new_defect(), and FEAT::Solver::IterativeSolver< Vector_ >::_update_defect().

◆ _plot_summary()

template<typename Vector_ >
bool FEAT::Solver::IterativeSolver< Vector_ >::_plot_summary ( ) const
inlineprotected

◆ _print_line()

◆ _progress()

template<typename Vector_ >
bool FEAT::Solver::IterativeSolver< Vector_ >::_progress ( ) const
inlineprotected

Progress iteration?

Returns
true if the solver should process, otherwise false.

Definition at line 741 of file iterative.hpp.

References FEAT::Solver::progress.

◆ _set_comm()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::_set_comm ( const Dist::Comm comm)
inlineprotected

Sets the communicator for the solver directly.

Parameters
[in]commA pointer to the communicator that is to be used by the solver.

Definition at line 667 of file iterative.hpp.

◆ _set_comm_by_matrix()

◆ _set_comm_by_vector()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::_set_comm_by_vector ( const Vector_ &  vector)
inlineprotected

Sets the communicator for the solver from a vector.

Parameters
[in]vectorA reference to a vector. If 'Vector_' is a 'Global::Vector', then the communicator of the vector gate is taken.

Definition at line 692 of file iterative.hpp.

◆ _set_initial_defect()

template<typename Vector_ >
virtual Status FEAT::Solver::IterativeSolver< Vector_ >::_set_initial_defect ( const VectorType vec_def,
const VectorType vec_sol 
)
inlineprotectedvirtual

◆ _set_new_defect()

◆ _update_defect()

template<typename Vector_ >
virtual Status FEAT::Solver::IterativeSolver< Vector_ >::_update_defect ( const DataType  def_cur_norm)
inlineprotectedvirtual

Internal function: sets the new (next) defect norm.

This function takes a pre-calculated defect vector norm, increments the iteration count, plots an output line to std::cout and checks whether any of the stopping criterions is fulfilled.

Parameters
[in]def_cur_normThe new defect norm.
Returns
The updated Status code.
Note
This function is preferred over _set_new_defect when using asynchronous mpi operations.

Definition at line 970 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_analyze_defect(), FEAT::Solver::IterativeSolver< Vector_ >::_def_cur, FEAT::Solver::IterativeSolver< Vector_ >::_num_iter, FEAT::Solver::IterativeSolver< Vector_ >::_plot_iter(), FEAT::Solver::IterativeSolver< Vector_ >::_plot_iter_line(), FEAT::Solver::IterativeSolver< Vector_ >::get_num_iter(), and FEAT::Solver::SolverBase< Vector_ >::name().

Referenced by FEAT::Solver::PipePCG< Matrix_, Filter_ >::_apply_intern(), FEAT::Solver::RBiCGStab< Matrix_, Filter_ >::_apply_intern(), and FEAT::Solver::GroppPCG< Matrix_, Filter_ >::_apply_intern().

◆ apply()

template<typename Vector_ >
virtual Status FEAT::Solver::SolverBase< Vector_ >::apply ( Vector_ &  vec_cor,
const Vector_ &  vec_def 
)
pure virtualinherited

Solver application method.

This method applies the solver represented by this object onto a given defect vector and returns the corresponding correction vector.

Note
Solvers which derive from the IterativeSolver base class also provide a correct() method which corrects an initial solution instead of starting with the null vector.
Parameters
[out]vec_corThe vector that shall receive the solution of the linear system. It is assumed to be allocated, but its numerical contents may be undefined upon calling this method.
[in]vec_defThe vector that represents the right-hand-side of the linear system to be solved.
Attention
vec_cor and vec_def must not refer to the same vector object!
Returns
A Status code that represents the status of the solution step.

Implemented in FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::CuSolverLU, FEAT::Solver::CuSolverQR, FEAT::Solver::ScalePrecond< Vector_, Filter_ >, FEAT::Solver::DiagonalPrecond< Vector_, Filter_ >, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::NLOptPrecond< VectorType_, FilterType_ >, FEAT::Solver::AmaVanka< Matrix_, Filter_ >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, and FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >.

Referenced by FEAT::Solver::solve().

◆ calc_convergence_rate()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::calc_convergence_rate ( ) const
inline

Computes the overall convergence rate: (defect_final / defect_initial) ^ (1 / number of iterations)

Definition at line 582 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_def_cur, FEAT::Solver::IterativeSolver< Vector_ >::_def_init, FEAT::Solver::IterativeSolver< Vector_ >::_num_iter, and FEAT::Math::pow().

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::get_summary().

◆ calc_defect_reduction()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::calc_defect_reduction ( ) const
inline

Computes the overall defect reduction factor: (defect_final / defect_inital)

Definition at line 597 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_def_init, and FEAT::Math::abs().

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::get_summary().

◆ correct()

template<typename Vector_ >
virtual Status FEAT::Solver::IterativeSolver< Vector_ >::correct ( VectorType vec_sol,
const VectorType vec_rhs 
)
pure virtual

Solver correction method.

This method applies the solver represented by this object onto a given right-hand-side vector and updates the corresponding solution vector.

In contrast to the apply() method of the SolverBase base class, this method uses the vector vec_sol as the initial solution vector for the iterative solution process instead of ignoring its contents upon entry and starting with the null vector.

Parameters
[in,out]vec_solThe vector that contains the initial solution upon entry and receives the solution of the linear system upon exit.
[in]vec_rhsThe vector that represents the right-hand-side of the linear system to be solved.
Attention
vec_sol and vec_rhs must not refer to the same vector object!
Returns
A Status code that represents the status of the solution step.

Implemented in FEAT::Solver::FixedStepLinesearch< Functional_, Filter_ >, FEAT::Solver::MQCLinesearch< Functional_, Filter_ >, FEAT::Solver::NewtonRaphsonLinesearch< Functional_, Filter_ >, FEAT::Solver::SecantLinesearch< Functional_, Filter_ >, and FEAT::Solver::QPenalty< Functional_ >.

Referenced by FEAT::Solver::solve().

◆ done()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::done ( )
inlinevirtualinherited

Finalization method.

This function performs both the symbolic and numeric finalization, i.e. it simply performs

this->done_numeric();
this->done_symbolic();

Definition at line 283 of file base.hpp.

References FEAT::Solver::SolverBase< Vector_ >::done_numeric(), and FEAT::Solver::SolverBase< Vector_ >::done_symbolic().

◆ done_numeric()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::done_numeric ( )
inlinevirtualinherited

Numeric finalization method.

This method is called to release any data allocated in the numeric initialization step.

Reimplemented in FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, FEAT::Solver::CUDSS, FEAT::Solver::DirectStokesSolver< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, LAFEM::TupleFilter< FilterV_, FilterP_ > >, FEAT::Solver::DirectStokesSolver< Global::Matrix< LocalMatrix_, Mirror_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::ParaSailsPrecond< Matrix_, Filter_ >, FEAT::Solver::EuclidPrecond< Matrix_, Filter_ >, FEAT::Solver::BoomerAMG< Matrix_, Filter_ >, FEAT::Solver::PreconditionedIterativeSolver< Vector_ >, FEAT::Solver::PreconditionedIterativeSolver< Matrix_::VectorTypeR >, FEAT::Solver::PreconditionedIterativeSolver< Functional_::VectorTypeR >, FEAT::Solver::MKLDSS, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::PCGNR< Matrix_, Filter_ >, FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::Umfpack, FEAT::Solver::UmfpackMean, FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >, FEAT::Solver::GenericUmfpack< Matrix_ >, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >.

Definition at line 246 of file base.hpp.

Referenced by FEAT::Solver::SolverBase< Vector_ >::done(), FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_numeric(), FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_numeric(), FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >::done_numeric(), FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >::done_numeric(), FEAT::Solver::ParaSailsPrecond< Matrix_, Filter_ >::done_numeric(), FEAT::Solver::PreconditionedIterativeSolver< Vector_ >::done_numeric(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_numeric(), FEAT::Solver::PCGNR< Matrix_, Filter_ >::done_numeric(), FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >::done_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_numeric(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_numeric().

◆ done_symbolic()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::done_symbolic ( )
inlinevirtualinherited

Symbolic finalization method.

This method is called to release any data allocated in the symbolic initialization step.

Reimplemented in FEAT::Solver::ALGLIBMinLBFGS< Functional_, Filter_ >, FEAT::Solver::ALGLIBMinCG< Functional_, Filter_ >, FEAT::Solver::AmaVanka< Matrix_, Filter_ >, FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::BiCGStab< Matrix_, Filter_ >, FEAT::Solver::BiCGStabL< Matrix_, Filter_ >, FEAT::Solver::Chebyshev< Matrix_, Filter_ >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, FEAT::Solver::CUDSS, FEAT::Solver::Descent< Matrix_, Filter_ >, FEAT::Solver::DirectStokesSolver< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, LAFEM::TupleFilter< FilterV_, FilterP_ > >, FEAT::Solver::DirectStokesSolver< Global::Matrix< LocalMatrix_, Mirror_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::FGMRES< Matrix_, Filter_ >, FEAT::Solver::GMRES< Matrix_, Filter_ >, FEAT::Solver::GroppPCG< Matrix_, Filter_ >, FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >, FEAT::Solver::IDRS< Matrix_, Filter_ >, FEAT::Solver::ILUPrecond< Matrix_, Filter_ >, FEAT::Solver::PreconditionedIterativeSolver< Vector_ >, FEAT::Solver::PreconditionedIterativeSolver< Matrix_::VectorTypeR >, FEAT::Solver::PreconditionedIterativeSolver< Functional_::VectorTypeR >, FEAT::Solver::JacobiPrecond< Matrix_, Filter_ >, FEAT::Solver::Linesearch< Functional_, Filter_ >, FEAT::Solver::MatrixPrecond< Matrix_, Filter_ >, FEAT::Solver::MKLDSS, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::NLCG< Functional_, Filter_ >, FEAT::Solver::NLSD< Functional_, Filter_ >, FEAT::Solver::PCG< Matrix_, Filter_ >, FEAT::Solver::PCGNR< Matrix_, Filter_ >, FEAT::Solver::PCGNRILU< Matrix_, Filter_ >, FEAT::Solver::PCR< Matrix_, Filter_ >, FEAT::Solver::PipePCG< Matrix_, Filter_ >, FEAT::Solver::PMR< Matrix_, Filter_ >, FEAT::Solver::PolynomialPrecond< Matrix_, Filter_ >, FEAT::Solver::PSD< Matrix_, Filter_ >, FEAT::Solver::QPenalty< Functional_ >, FEAT::Solver::RBiCGStab< Matrix_, Filter_ >, FEAT::Solver::RGCR< Matrix_, Filter_ >, FEAT::Solver::Richardson< Matrix_, Filter_ >, FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::SORPrecond< Matrix_, Filter_ >, FEAT::Solver::SSORPrecond< Matrix_, Filter_ >, FEAT::Solver::SuperLU< Matrix_, Filter_ >, FEAT::Solver::Umfpack, FEAT::Solver::UmfpackMean, FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >, FEAT::Solver::GenericUmfpack< Matrix_ >, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::Vanka< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, Filter_ >, and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.

Definition at line 255 of file base.hpp.

Referenced by FEAT::Solver::SolverBase< Vector_ >::done(), FEAT::Solver::ALGLIBMinLBFGS< Functional_, Filter_ >::done_symbolic(), FEAT::Solver::ALGLIBMinCG< Functional_, Filter_ >::done_symbolic(), FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic(), FEAT::Solver::Chebyshev< Matrix_, Filter_ >::done_symbolic(), FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >::done_symbolic(), FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >::done_symbolic(), FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >::done_symbolic(), FEAT::Solver::PreconditionedIterativeSolver< Vector_ >::done_symbolic(), FEAT::Solver::Linesearch< Functional_, Filter_ >::done_symbolic(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::done_symbolic(), FEAT::Solver::NLCG< Functional_, Filter_ >::done_symbolic(), FEAT::Solver::NLSD< Functional_, Filter_ >::done_symbolic(), FEAT::Solver::PCGNR< Matrix_, Filter_ >::done_symbolic(), FEAT::Solver::PCGNRILU< Matrix_, Filter_ >::done_symbolic(), FEAT::Solver::QPenalty< Functional_ >::done_symbolic(), FEAT::Solver::SuperLU< Matrix_, Filter_ >::done_symbolic(), FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >::done_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::done_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::done_symbolic().

◆ force_defect_norm_calc()

template<typename Vector_ >
virtual void FEAT::Solver::IterativeSolver< Vector_ >::force_defect_norm_calc ( bool  force)
inlinevirtual

Forces the calculation of defect norms in every iteration.

Note
Please note that allowing the skipping of defect norm calculations is merely a hint and a derived class may override this function therefore always forcing the calculation of defect norms. One possible reason is that the derived solver class may require defect norms as part of its solver algorithm and therefore skipping the computation of the defect norms may break the solver; two examples are the BiCGStab and GMRES solvers.
Parameters
[in]forceSpecifies whether defect norm computation is to be enforced in every iteration

Reimplemented in FEAT::Solver::BiCGStab< Matrix_, Filter_ >, FEAT::Solver::BiCGStabL< Matrix_, Filter_ >, FEAT::Solver::FGMRES< Matrix_, Filter_ >, FEAT::Solver::GMRES< Matrix_, Filter_ >, and FEAT::Solver::RBiCGStab< Matrix_, Filter_ >.

Definition at line 490 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_force_def_norm_calc.

◆ get_def_final()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_def_final ( ) const
inline

Returns the final defect.

Definition at line 570 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_def_cur.

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::get_summary().

◆ get_def_initial()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_def_initial ( ) const
inline

Returns the initial defect.

Definition at line 564 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_def_init.

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::get_summary().

◆ get_div_abs()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_div_abs ( ) const
inline

Returns the absolute divergence.

Definition at line 419 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_div_abs.

◆ get_div_rel()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_div_rel ( ) const
inline

Returns the relative divergence.

Definition at line 413 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_div_rel.

◆ get_max_iter()

template<typename Vector_ >
Index FEAT::Solver::IterativeSolver< Vector_ >::get_max_iter ( ) const
inline

Returns the maximum number of iterations.

Definition at line 474 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_max_iter.

◆ get_min_iter()

template<typename Vector_ >
Index FEAT::Solver::IterativeSolver< Vector_ >::get_min_iter ( ) const
inline

Returns the minimal number of iterations.

Definition at line 468 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_min_iter.

◆ get_min_stag_iter()

template<typename Vector_ >
Index FEAT::Solver::IterativeSolver< Vector_ >::get_min_stag_iter ( ) const
inline

Returns the minimum stagnation iteration count.

Definition at line 443 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_min_stag_iter.

◆ get_num_iter()

◆ get_plot_name()

template<typename Vector_ >
String FEAT::Solver::IterativeSolver< Vector_ >::get_plot_name ( ) const
inline

◆ get_stag_rate()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_stag_rate ( ) const
inline

Returns the stagnation rate.

Definition at line 431 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_stag_rate.

◆ get_status()

template<typename Vector_ >
Status FEAT::Solver::IterativeSolver< Vector_ >::get_status ( ) const
inline

◆ get_summary()

◆ get_tol_abs()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_tol_abs ( ) const
inline

Returns the absolute tolerance.

Definition at line 389 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_tol_abs.

◆ get_tol_abs_low()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_tol_abs_low ( ) const
inline

Returns the lower absolute tolerance.

Definition at line 395 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_tol_abs_low.

◆ get_tol_rel()

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::get_tol_rel ( ) const
inline

Returns the relative tolerance.

Definition at line 383 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_tol_rel.

◆ init()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::init ( )
inlinevirtualinherited

Initialization method.

This function performs both the symbolic and numeric initialization, i.e. it simply performs

this->init_symbolic();
this->init_numeric();

Definition at line 268 of file base.hpp.

References FEAT::Solver::SolverBase< Vector_ >::init_numeric(), and FEAT::Solver::SolverBase< Vector_ >::init_symbolic().

◆ init_numeric()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::init_numeric ( )
inlinevirtualinherited

Numeric initialization method.

This method is called to perform numeric initialization of the solver.
Before this function can be called, the symbolic initialization must be performed.

Reimplemented in FEAT::Solver::AmaVanka< Matrix_, Filter_ >, FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::Chebyshev< Matrix_, Filter_ >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, FEAT::Solver::CUDSS, FEAT::Solver::DirectStokesSolver< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, LAFEM::TupleFilter< FilterV_, FilterP_ > >, FEAT::Solver::DirectStokesSolver< Global::Matrix< LocalMatrix_, Mirror_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >, FEAT::Solver::ParaSailsPrecond< Matrix_, Filter_ >, FEAT::Solver::EuclidPrecond< Matrix_, Filter_ >, FEAT::Solver::BoomerAMG< Matrix_, Filter_ >, FEAT::Solver::ILUPrecond< Matrix_, Filter_ >, FEAT::Solver::PreconditionedIterativeSolver< Vector_ >, FEAT::Solver::PreconditionedIterativeSolver< Matrix_::VectorTypeR >, FEAT::Solver::PreconditionedIterativeSolver< Functional_::VectorTypeR >, FEAT::Solver::JacobiPrecond< Matrix_, Filter_ >, FEAT::Solver::MatrixPrecond< Matrix_, Filter_ >, FEAT::Solver::MKLDSS, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::NonlinearOperatorPrecondWrapper< NonlinearOperator_ >, FEAT::Solver::PCGNR< Matrix_, Filter_ >, FEAT::Solver::PCGNRILU< Matrix_, Filter_ >, FEAT::Solver::PolynomialPrecond< Matrix_, Filter_ >, FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::SuperLU< Matrix_, Filter_ >, FEAT::Solver::Umfpack, FEAT::Solver::UmfpackMean, FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >, FEAT::Solver::GenericUmfpack< Matrix_ >, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::Vanka< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, Filter_ >, and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.

Definition at line 237 of file base.hpp.

Referenced by FEAT::Solver::SolverBase< Vector_ >::init(), FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_numeric(), FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_numeric(), FEAT::Solver::Chebyshev< Matrix_, Filter_ >::init_numeric(), FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >::init_numeric(), FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >::init_numeric(), FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >::init_numeric(), FEAT::Solver::PreconditionedIterativeSolver< Vector_ >::init_numeric(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_numeric(), FEAT::Solver::PCGNR< Matrix_, Filter_ >::init_numeric(), FEAT::Solver::PCGNRILU< Matrix_, Filter_ >::init_numeric(), FEAT::Solver::SuperLU< Matrix_, Filter_ >::init_numeric(), FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >::init_numeric(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_numeric(), FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_numeric(), and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >::init_numeric().

◆ init_symbolic()

template<typename Vector_ >
virtual void FEAT::Solver::SolverBase< Vector_ >::init_symbolic ( )
inlinevirtualinherited

Symbolic initialization method.

This method is called to perform symbolic initialization of the solver.

Reimplemented in FEAT::Solver::ALGLIBMinLBFGS< Functional_, Filter_ >, FEAT::Solver::ALGLIBMinCG< Functional_, Filter_ >, FEAT::Solver::AmaVanka< Matrix_, Filter_ >, FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::BiCGStab< Matrix_, Filter_ >, FEAT::Solver::BiCGStabL< Matrix_, Filter_ >, FEAT::Solver::Chebyshev< Matrix_, Filter_ >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, FEAT::Solver::CUDSS, FEAT::Solver::Descent< Matrix_, Filter_ >, FEAT::Solver::DirectStokesSolver< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, LAFEM::TupleFilter< FilterV_, FilterP_ > >, FEAT::Solver::DirectStokesSolver< Global::Matrix< LocalMatrix_, Mirror_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::FGMRES< Matrix_, Filter_ >, FEAT::Solver::GMRES< Matrix_, Filter_ >, FEAT::Solver::GroppPCG< Matrix_, Filter_ >, FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >, FEAT::Solver::IDRS< Matrix_, Filter_ >, FEAT::Solver::ILUPrecond< Matrix_, Filter_ >, FEAT::Solver::PreconditionedIterativeSolver< Vector_ >, FEAT::Solver::PreconditionedIterativeSolver< Matrix_::VectorTypeR >, FEAT::Solver::PreconditionedIterativeSolver< Functional_::VectorTypeR >, FEAT::Solver::JacobiPrecond< Matrix_, Filter_ >, FEAT::Solver::Linesearch< Functional_, Filter_ >, FEAT::Solver::MatrixPrecond< Matrix_, Filter_ >, FEAT::Solver::MKLDSS, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::NLCG< Functional_, Filter_ >, FEAT::Solver::NLSD< Functional_, Filter_ >, FEAT::Solver::PCG< Matrix_, Filter_ >, FEAT::Solver::PCGNR< Matrix_, Filter_ >, FEAT::Solver::PCGNRILU< Matrix_, Filter_ >, FEAT::Solver::PCR< Matrix_, Filter_ >, FEAT::Solver::PipePCG< Matrix_, Filter_ >, FEAT::Solver::PMR< Matrix_, Filter_ >, FEAT::Solver::PolynomialPrecond< Matrix_, Filter_ >, FEAT::Solver::PSD< Matrix_, Filter_ >, FEAT::Solver::QPenalty< Functional_ >, FEAT::Solver::RBiCGStab< Matrix_, Filter_ >, FEAT::Solver::RGCR< Matrix_, Filter_ >, FEAT::Solver::Richardson< Matrix_, Filter_ >, FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::SORPrecond< Matrix_, Filter_ >, FEAT::Solver::SSORPrecond< Matrix_, Filter_ >, FEAT::Solver::SuperLU< Matrix_, Filter_ >, FEAT::Solver::Umfpack, FEAT::Solver::UmfpackMean, FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >, FEAT::Solver::GenericUmfpack< Matrix_ >, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::Vanka< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, Filter_ >, and FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >.

Definition at line 227 of file base.hpp.

Referenced by FEAT::Solver::SolverBase< Vector_ >::init(), FEAT::Solver::ALGLIBMinLBFGS< Functional_, Filter_ >::init_symbolic(), FEAT::Solver::ALGLIBMinCG< Functional_, Filter_ >::init_symbolic(), FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic(), FEAT::Solver::Chebyshev< Matrix_, Filter_ >::init_symbolic(), FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >::init_symbolic(), FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >::init_symbolic(), FEAT::Solver::HypreSolverBase< Matrix_, Filter_, SolverBase_ >::init_symbolic(), FEAT::Solver::PreconditionedIterativeSolver< Vector_ >::init_symbolic(), FEAT::Solver::Linesearch< Functional_, Filter_ >::init_symbolic(), FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >::init_symbolic(), FEAT::Solver::NLCG< Functional_, Filter_ >::init_symbolic(), FEAT::Solver::NLSD< Functional_, Filter_ >::init_symbolic(), FEAT::Solver::PCGNR< Matrix_, Filter_ >::init_symbolic(), FEAT::Solver::PCGNRILU< Matrix_, Filter_ >::init_symbolic(), FEAT::Solver::QPenalty< Functional_ >::init_symbolic(), FEAT::Solver::SuperLU< Matrix_, Filter_ >::init_symbolic(), FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >::init_symbolic(), FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >::init_symbolic(), and FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >::init_symbolic().

◆ is_converged() [1/2]

◆ is_converged() [2/2]

template<typename Vector_ >
bool FEAT::Solver::IterativeSolver< Vector_ >::is_converged ( const DataType  def_cur) const
inline

checks for convergence

Parameters
[in]def_curThe defect norm that is to be analyzed
Returns
true, if converged, else false

Definition at line 546 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_def_init, FEAT::Solver::IterativeSolver< Vector_ >::_tol_abs, FEAT::Solver::IterativeSolver< Vector_ >::_tol_abs_low, and FEAT::Solver::IterativeSolver< Vector_ >::_tol_rel.

◆ is_diverged() [1/2]

◆ is_diverged() [2/2]

template<typename Vector_ >
bool FEAT::Solver::IterativeSolver< Vector_ >::is_diverged ( const DataType  def_cur) const
inline

◆ name()

template<typename Vector_ >
virtual String FEAT::Solver::SolverBase< Vector_ >::name ( ) const
pure virtualinherited

Returns a descriptive string.

Returns
A string describing the solver.

Implemented in FEAT::Solver::ALGLIBMinLBFGS< Functional_, Filter_ >, FEAT::Solver::ALGLIBMinCG< Functional_, Filter_ >, FEAT::Solver::AmaVanka< Matrix_, Filter_ >, FEAT::Solver::BFBT< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::BFBT< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::BiCGStab< Matrix_, Filter_ >, FEAT::Solver::BiCGStabL< Matrix_, Filter_ >, FEAT::Solver::Chebyshev< Matrix_, Filter_ >, FEAT::Solver::ConvertPrecond< VectorOuter_, VectorInner_ >, FEAT::Solver::ConvertPrecond< Global::Vector< LocalVectorOuter_, MirrorOuter_ >, Global::Vector< LocalVectorInner_, MirrorInner_ > >, FEAT::Solver::CUDSS, FEAT::Solver::CuSolverLU, FEAT::Solver::CuSolverQR, FEAT::Solver::Descent< Matrix_, Filter_ >, FEAT::Solver::DiagonalPrecond< Vector_, Filter_ >, FEAT::Solver::DirectStokesSolver< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, LAFEM::TupleFilter< FilterV_, FilterP_ > >, FEAT::Solver::DirectStokesSolver< Global::Matrix< LocalMatrix_, Mirror_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::FGMRES< Matrix_, Filter_ >, FEAT::Solver::FixedStepLinesearch< Functional_, Filter_ >, FEAT::Solver::GMRES< Matrix_, Filter_ >, FEAT::Solver::GroppPCG< Matrix_, Filter_ >, FEAT::Solver::ParaSailsPrecond< Matrix_, Filter_ >, FEAT::Solver::EuclidPrecond< Matrix_, Filter_ >, FEAT::Solver::BoomerAMG< Matrix_, Filter_ >, FEAT::Solver::IDRS< Matrix_, Filter_ >, FEAT::Solver::ILUPrecond< Matrix_, Filter_ >, FEAT::Solver::JacobiPrecond< Matrix_, Filter_ >, FEAT::Solver::MatrixPrecond< Matrix_, Filter_ >, FEAT::Solver::MKLDSS, FEAT::Solver::MQCLinesearch< Functional_, Filter_ >, FEAT::Solver::MultiGrid< SystemMatrix_, SystemFilter_, TransferOperator_ >, FEAT::Solver::NewtonRaphsonLinesearch< Functional_, Filter_ >, FEAT::Solver::NLCG< Functional_, Filter_ >, FEAT::Solver::NonlinearOperatorPrecondWrapper< NonlinearOperator_ >, FEAT::Solver::NLSD< Functional_, Filter_ >, FEAT::Solver::PCG< Matrix_, Filter_ >, FEAT::Solver::PCGNR< Matrix_, Filter_ >, FEAT::Solver::PCGNRILU< Matrix_, Filter_ >, FEAT::Solver::PCR< Matrix_, Filter_ >, FEAT::Solver::PipePCG< Matrix_, Filter_ >, FEAT::Solver::PMR< Matrix_, Filter_ >, FEAT::Solver::PolynomialPrecond< Matrix_, Filter_ >, FEAT::Solver::PreconWrapper< Matrix_, Filter_, Precon_ >, FEAT::Solver::PSD< Matrix_, Filter_ >, FEAT::Solver::QPenalty< Functional_ >, FEAT::Solver::RBiCGStab< Matrix_, Filter_ >, FEAT::Solver::RGCR< Matrix_, Filter_ >, FEAT::Solver::Richardson< Matrix_, Filter_ >, FEAT::Solver::ScalePrecond< Vector_, Filter_ >, FEAT::Solver::SchwarzPrecond< Global::Vector< LocalVector_, Mirror_ >, Global::Filter< LocalFilter_, Mirror_ > >, FEAT::Solver::SecantLinesearch< Functional_, Filter_ >, FEAT::Solver::SORPrecond< Matrix_, Filter_ >, FEAT::Solver::SSORPrecond< Matrix_, Filter_ >, FEAT::Solver::SuperLU< Matrix_, Filter_ >, FEAT::Solver::Umfpack, FEAT::Solver::UmfpackMean, FEAT::Solver::SaddleUmfpackMean< DT_, IT_, dim_ >, FEAT::Solver::GenericUmfpack< Matrix_ >, FEAT::Solver::UzawaPrecond< MatrixA_, MatrixB_, MatrixD_, FilterV_, FilterP_ >, FEAT::Solver::UzawaPrecond< Global::Matrix< MatrixA_, MirrorV_, MirrorV_ >, Global::Matrix< MatrixB_, MirrorV_, MirrorP_ >, Global::Matrix< MatrixD_, MirrorP_, MirrorV_ >, Global::Filter< FilterV_, MirrorV_ >, Global::Filter< FilterP_, MirrorP_ > >, FEAT::Solver::Vanka< LAFEM::SaddlePointMatrix< MatrixA_, MatrixB_, MatrixD_ >, Filter_ >, FEAT::Solver::VoxelAmaVanka< Matrix_, Filter_, pol_threading_, macro_type_ >, FEAT::Solver::NLOptPrecond< VectorType_, FilterType_ >, FEAT::Solver::NLOptPrecond< OperatorType_::VectorType, Filter >, FEAT::Solver::NLOptPrecond< Operator_::VectorType, Filter_ >, and FEAT::Solver::NLOptPrecond< NonlinearOperator_::SystemLevelType::GlobalSystemVectorR, NonlinearOperator_::SystemLevelType::GlobalSystemFilter >.

Referenced by FEAT::Solver::PreconditionedIterativeSolver< Vector_ >::_apply_precond(), FEAT::Solver::Linesearch< Functional_, Filter_ >::_check_convergence(), FEAT::Solver::IterativeSolver< Vector_ >::_set_initial_defect(), FEAT::Solver::NLOptLS< Functional_, Filter_ >::_set_initial_defect(), FEAT::Solver::IterativeSolver< Vector_ >::_set_new_defect(), FEAT::Solver::NLOptLS< Functional_, Filter_ >::_set_new_defect(), and FEAT::Solver::IterativeSolver< Vector_ >::_update_defect().

◆ plot_summary()

template<typename Vector_ >
virtual void FEAT::Solver::IterativeSolver< Vector_ >::plot_summary ( ) const
inlinevirtual

◆ set_div_abs()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_div_abs ( DataType  div_abs)
inline

Sets the absolute divergence for the solver.

Definition at line 407 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_div_abs.

◆ set_div_rel()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_div_rel ( DataType  div_rel)
inline

Sets the relative divergence for the solver.

Definition at line 401 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_div_rel.

◆ set_max_iter()

◆ set_min_iter()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_min_iter ( Index  min_iter)
inline

Sets the minimum iteration count for the solver.

Definition at line 449 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_min_iter.

◆ set_min_stag_iter()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_min_stag_iter ( Index  min_iter)
inline

Sets the minimum stagnate iteration count for the solver.

Definition at line 437 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_min_stag_iter.

◆ set_plot_interval()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_plot_interval ( const Index  plot_interval)
inline

Sets the interval between two plot outputs, if any.

Parameters
[in]plot_intervalThe desired interval of iteration plots

Definition at line 511 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_plot_interval.

◆ set_plot_mode()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_plot_mode ( const PlotMode  plot_mode)
inline

Sets the plot mode of the solver.

Parameters
[in]plot_modeThe desired plot mode

Definition at line 501 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_plot_mode.

◆ set_plot_name()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_plot_name ( const String plot_name)
inline

◆ set_stag_rate()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_stag_rate ( DataType  rate)
inline

Sets the stagnation rate for the solver.

Definition at line 425 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_stag_rate.

◆ set_tol_abs()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_tol_abs ( DataType  tol_abs)
inline

◆ set_tol_abs_low()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_tol_abs_low ( DataType  tol_abs_low)
inline

Sets the lower absolute tolerance for the solver.

Definition at line 377 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_tol_abs_low.

◆ set_tol_rel()

template<typename Vector_ >
void FEAT::Solver::IterativeSolver< Vector_ >::set_tol_rel ( DataType  tol_rel)
inline

Sets the relative tolerance for the solver.

Definition at line 365 of file iterative.hpp.

References FEAT::Solver::IterativeSolver< Vector_ >::_tol_rel.

Member Data Documentation

◆ _comm

template<typename Vector_ >
const Dist::Comm* FEAT::Solver::IterativeSolver< Vector_ >::_comm
protected

Communicator of the solver.

Definition at line 209 of file iterative.hpp.

Referenced by FEAT::Solver::IterativeSolver< Vector_ >::_print_line().

◆ _def_cur

◆ _def_init

◆ _def_prev

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::_def_prev
protected

◆ _div_abs

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::_div_abs
protected

◆ _div_rel

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::_div_rel
protected

◆ _force_def_norm_calc

◆ _iter_digits

◆ _max_iter

◆ _min_iter

◆ _min_stag_iter

◆ _num_iter

◆ _num_stag_iter

template<typename Vector_ >
Index FEAT::Solver::IterativeSolver< Vector_ >::_num_stag_iter
protected

◆ _plot_interval

template<typename Vector_ >
Index FEAT::Solver::IterativeSolver< Vector_ >::_plot_interval
protected

◆ _plot_mode

◆ _plot_name

◆ _stag_rate

template<typename Vector_ >
DataType FEAT::Solver::IterativeSolver< Vector_ >::_stag_rate
protected

◆ _status

◆ _tol_abs

◆ _tol_abs_low

◆ _tol_rel


The documentation for this class was generated from the following file: