Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
IterationNumberControl Class Reference

#include <deal.II/lac/solver_control.h>

Inheritance diagram for IterationNumberControl:
[legend]

Public Types

enum  State { iterate = 0 , success , failure }
 

Public Member Functions

 IterationNumberControl (const unsigned int maxiter=100, const double tolerance=1e-12, const bool log_history=false, const bool log_result=true)
 
 IterationNumberControl (const SolverControl &c)
 
IterationNumberControloperator= (const SolverControl &c)
 
virtual ~IterationNumberControl () override=default
 
virtual State check (const unsigned int step, const double check_value) override
 
void parse_parameters (ParameterHandler &param)
 
State last_check () const
 
double initial_value () const
 
double last_value () const
 
unsigned int last_step () const
 
unsigned int max_steps () const
 
unsigned int set_max_steps (const unsigned int)
 
void set_failure_criterion (const double rel_failure_residual)
 
void clear_failure_criterion ()
 
double tolerance () const
 
double set_tolerance (const double)
 
void enable_history_data ()
 
const std::vector< double > & get_history_data () const
 
double average_reduction () const
 
double final_reduction () const
 
double step_reduction (unsigned int step) const
 
void log_history (const bool)
 
bool log_history () const
 
unsigned int log_frequency (unsigned int)
 
void log_result (const bool)
 
bool log_result () const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &param)
 
static ::ExceptionBaseExcHistoryDataRequired ()
 

Protected Attributes

unsigned int maxsteps
 
double tol
 
State lcheck
 
double initial_val
 
double lvalue
 
unsigned int lstep
 
bool check_failure
 
double relative_failure_residual
 
double failure_residual
 
bool m_log_history
 
unsigned int m_log_frequency
 
bool m_log_result
 
bool history_data_enabled
 
std::vector< double > history_data
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
static std::mutex mutex
 
void check_no_subscribers () const noexcept
 

Detailed Description

Specialization of SolverControl which returns success if a given number of iteration was performed, irrespective of the actual residual. This is useful in cases where you don't want to solve exactly, but rather want to perform a fixed number of iterations, e.g. in an inner solver. The arguments given to this class are exactly the same as for the SolverControl class and the solver terminates similarly when one of the given tolerance or the maximum iteration count were reached. The only difference to SolverControl is that the solver returns success in the latter case.

Definition at line 510 of file solver_control.h.

Member Enumeration Documentation

◆ State

enum SolverControl::State
inherited

Enum denoting the different states a solver can be in. See the general documentation of this class for more information.

Enumerator
iterate 

Continue iteration.

success 

Stop iteration, goal reached.

failure 

Stop iteration, goal not reached.

Definition at line 72 of file solver_control.h.

Constructor & Destructor Documentation

◆ IterationNumberControl() [1/2]

IterationNumberControl::IterationNumberControl ( const unsigned int  maxiter = 100,
const double  tolerance = 1e-12,
const bool  log_history = false,
const bool  log_result = true 
)
explicit

Constructor. Provide exactly the same arguments as the constructor of the SolverControl class.

Definition at line 312 of file solver_control.cc.

◆ IterationNumberControl() [2/2]

IterationNumberControl::IterationNumberControl ( const SolverControl c)

Initialize with a SolverControl object. The result will emulate SolverControl by setting the reduction target to zero.

◆ ~IterationNumberControl()

virtual IterationNumberControl::~IterationNumberControl ( )
overridevirtualdefault

Virtual destructor is needed as there are virtual functions in this class.

Member Function Documentation

◆ operator=()

IterationNumberControl & IterationNumberControl::operator= ( const SolverControl c)

Assign a SolverControl object to ReductionControl. The result of the assignment will emulate SolverControl by setting the reduction target to zero.

◆ check()

SolverControl::State IterationNumberControl::check ( const unsigned int  step,
const double  check_value 
)
overridevirtual

Decide about success or failure of an iteration. This function bases success solely on the fact if a given number of iterations was reached or the check value reached exactly zero.

Reimplemented from SolverControl.

Definition at line 322 of file solver_control.cc.

◆ declare_parameters()

void SolverControl::declare_parameters ( ParameterHandler param)
staticinherited

Interface to parameter file.

Definition at line 204 of file solver_control.cc.

◆ parse_parameters()

void SolverControl::parse_parameters ( ParameterHandler param)
inherited

Read parameters from file.

Definition at line 215 of file solver_control.cc.

◆ last_check()

SolverControl::State SolverControl::last_check ( ) const
inherited

Return the result of the last check operation.

Definition at line 106 of file solver_control.cc.

◆ initial_value()

double SolverControl::initial_value ( ) const
inherited

Return the initial convergence criterion.

Definition at line 113 of file solver_control.cc.

◆ last_value()

double SolverControl::last_value ( ) const
inherited

Return the convergence value of last iteration step for which check was called by the solver.

Definition at line 120 of file solver_control.cc.

◆ last_step()

unsigned int SolverControl::last_step ( ) const
inherited

Number of last iteration step.

Definition at line 127 of file solver_control.cc.

◆ max_steps()

unsigned int SolverControl::max_steps ( ) const
inherited

Maximum number of steps.

◆ set_max_steps()

unsigned int SolverControl::set_max_steps ( const unsigned int  )
inherited

Change maximum number of steps.

◆ set_failure_criterion()

void SolverControl::set_failure_criterion ( const double  rel_failure_residual)
inherited

Enables the failure check. Solving is stopped with ReturnState failure if residual>failure_residual with failure_residual := rel_failure_residual*first_residual.

◆ clear_failure_criterion()

void SolverControl::clear_failure_criterion ( )
inherited

Disables failure check and resets relative_failure_residual and failure_residual to zero.

◆ tolerance()

double SolverControl::tolerance ( ) const
inherited

Tolerance.

◆ set_tolerance()

double SolverControl::set_tolerance ( const double  )
inherited

Change tolerance.

◆ enable_history_data()

void SolverControl::enable_history_data ( )
inherited

Enables writing residuals of each step into a vector for later analysis.

Definition at line 145 of file solver_control.cc.

◆ get_history_data()

const std::vector< double > & SolverControl::get_history_data ( ) const
inherited

Provide read access to the collected residual data.

Definition at line 153 of file solver_control.cc.

◆ average_reduction()

double SolverControl::average_reduction ( ) const
inherited

Average error reduction over all steps.

Requires enable_history_data()

Definition at line 169 of file solver_control.cc.

◆ final_reduction()

double SolverControl::final_reduction ( ) const
inherited

Error reduction of the last step; for stationary iterations, this approximates the norm of the iteration matrix.

Requires enable_history_data()

Definition at line 197 of file solver_control.cc.

◆ step_reduction()

double SolverControl::step_reduction ( unsigned int  step) const
inherited

Error reduction of any iteration step.

Requires enable_history_data()

Definition at line 185 of file solver_control.cc.

◆ log_history() [1/2]

void SolverControl::log_history ( const bool  )
inherited

Log each iteration step. Use log_frequency for skipping steps.

◆ log_history() [2/2]

bool SolverControl::log_history ( ) const
inherited

Return the log_history flag.

◆ log_frequency()

unsigned int SolverControl::log_frequency ( unsigned int  f)
inherited

Set logging frequency.

Definition at line 134 of file solver_control.cc.

◆ log_result() [1/2]

void SolverControl::log_result ( const bool  )
inherited

Log start and end step.

◆ log_result() [2/2]

bool SolverControl::log_result ( ) const
inherited

Return the log_result flag.

Member Data Documentation

◆ maxsteps

unsigned int SolverControl::maxsteps
protectedinherited

Maximum number of steps.

Definition at line 333 of file solver_control.h.

◆ tol

double SolverControl::tol
protectedinherited

Prescribed tolerance to be achieved.

Definition at line 338 of file solver_control.h.

◆ lcheck

State SolverControl::lcheck
protectedinherited

Result of last check operation.

Definition at line 343 of file solver_control.h.

◆ initial_val

double SolverControl::initial_val
protectedinherited

Initial value.

Definition at line 348 of file solver_control.h.

◆ lvalue

double SolverControl::lvalue
protectedinherited

Last value of the convergence criterion.

Definition at line 353 of file solver_control.h.

◆ lstep

unsigned int SolverControl::lstep
protectedinherited

Last step.

Definition at line 358 of file solver_control.h.

◆ check_failure

bool SolverControl::check_failure
protectedinherited

Is set to true by set_failure_criterion and enables failure checking.

Definition at line 364 of file solver_control.h.

◆ relative_failure_residual

double SolverControl::relative_failure_residual
protectedinherited

Stores the rel_failure_residual set by set_failure_criterion

Definition at line 369 of file solver_control.h.

◆ failure_residual

double SolverControl::failure_residual
protectedinherited

failure_residual equals the first residual multiplied by relative_crit set by set_failure_criterion (see there).

Until the first residual is known it is 0.

Definition at line 377 of file solver_control.h.

◆ m_log_history

bool SolverControl::m_log_history
protectedinherited

Log convergence history to deallog.

Definition at line 382 of file solver_control.h.

◆ m_log_frequency

unsigned int SolverControl::m_log_frequency
protectedinherited

Log only every nth step.

Definition at line 387 of file solver_control.h.

◆ m_log_result

bool SolverControl::m_log_result
protectedinherited

Log iteration result to deallog. If true, after finishing the iteration, a statement about failure or success together with lstep and lvalue are logged.

Definition at line 394 of file solver_control.h.

◆ history_data_enabled

bool SolverControl::history_data_enabled
protectedinherited

Control over the storage of history data. Set by enable_history_data().

Definition at line 399 of file solver_control.h.

◆ history_data

std::vector<double> SolverControl::history_data
protectedinherited

Vector storing the result after each iteration step for later statistical analysis.

Use of this vector is enabled by enable_history_data().

Definition at line 407 of file solver_control.h.


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