Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/solver_control.h>
Public Member Functions | |
ConsecutiveControl (const unsigned int maxiter=100, const double tolerance=1.e-10, const unsigned int n_consecutive_iterations=2, const bool log_history=false, const bool log_result=false) | |
ConsecutiveControl (const SolverControl &c) | |
ConsecutiveControl & | operator= (const SolverControl &c) |
virtual | ~ConsecutiveControl ()=default |
virtual State | check (const unsigned int step, const double check_value) |
Public Member Functions inherited from SolverControl | |
SolverControl (const unsigned int n=100, const double tol=1.e-10, const bool log_history=false, const bool log_result=true) | |
virtual | ~SolverControl ()=default |
void | parse_parameters (ParameterHandler ¶m) |
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 |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Attributes | |
unsigned int | n_consecutive_iterations |
unsigned int | n_converged_iterations |
Protected Attributes inherited from SolverControl | |
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 |
Additional Inherited Members | |
Public Types inherited from SolverControl | |
enum | State { iterate = 0, success, failure } |
Static Public Member Functions inherited from SolverControl | |
static void | declare_parameters (ParameterHandler ¶m) |
static ::ExceptionBase & | ExcHistoryDataRequired () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Specialization of SolverControl
which returns SolverControl::State::success if and only if a certain positive number of consecutive iterations satisfy the specified tolerance. This is useful in cases when solving nonlinear problems using inexact Hessian.
For example: The requested number of consecutively converged iterations is 2, the tolerance is 0.2. The ConsecutiveControl will return SolverControl::State::success only at the last step in the sequence 0.5, 0.0005, 1.0, 0.05, 0.01.
Definition at line 540 of file solver_control.h.
ConsecutiveControl::ConsecutiveControl | ( | const unsigned int | maxiter = 100 , |
const double | tolerance = 1.e-10 , |
||
const unsigned int | n_consecutive_iterations = 2 , |
||
const bool | log_history = false , |
||
const bool | log_result = false |
||
) |
Constructor. n_consecutive_iterations
is the number of consecutive iterations which should satisfy the prescribed tolerance for convergence. Other arguments have the same meaning as those of the constructor of the SolverControl.
Definition at line 348 of file solver_control.cc.
ConsecutiveControl::ConsecutiveControl | ( | const SolverControl & | c | ) |
Initialize with a SolverControl object. The result will emulate SolverControl by setting n_consecutive_iterations
to one.
Definition at line 364 of file solver_control.cc.
|
virtualdefault |
Virtual destructor is needed as there are virtual functions in this class.
ConsecutiveControl & ConsecutiveControl::operator= | ( | const SolverControl & | c | ) |
Assign a SolverControl object to ConsecutiveControl. The result of the assignment will emulate SolverControl by setting n_consecutive_iterations
to one.
Definition at line 374 of file solver_control.cc.
|
virtual |
Decide about success or failure of an iteration, see the class description above.
Reimplemented from SolverControl.
Definition at line 385 of file solver_control.cc.
|
protected |
The number of consecutive iterations which should satisfy the prescribed tolerance for convergence.
Definition at line 586 of file solver_control.h.
|
protected |
Counter for the number of consecutively converged iterations.
Definition at line 591 of file solver_control.h.