Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/solver_control.h>
Public Member Functions | |
ReductionControl (const unsigned int maxiter=100, const double tolerance=1.e-10, const double reduce=1.e-2, const bool log_history=false, const bool log_result=true) | |
ReductionControl (const SolverControl &c) | |
ReductionControl & | operator= (const SolverControl &c) |
virtual | ~ReductionControl ()=default |
void | parse_parameters (ParameterHandler ¶m) |
virtual State | check (const unsigned int step, const double check_value) |
double | reduction () const |
double | set_reduction (const double) |
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) |
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler ¶m) |
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) |
Protected Attributes | |
double | reduce |
double | reduced_tol |
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 } |
Specialization of SolverControl
which returns success
if either the specified tolerance is achieved or if the initial residual (or whatever criterion was chosen by the solver class) is reduced by a given factor. This is useful in cases where you don't want to solve exactly, but rather want to gain two digits or if the maximal number of iterations is achieved. For example: The maximal number of iterations is 20, the reduction factor is 1% and the tolerance is 0.1%. The initial residual is 2.5. The process will break if 20 iteration are completed or the new residual is less then 2.5*1% or if it is less then 0.1%.
Definition at line 402 of file solver_control.h.
ReductionControl::ReductionControl | ( | const unsigned int | maxiter = 100 , |
const double | tolerance = 1.e-10 , |
||
const double | reduce = 1.e-2 , |
||
const bool | log_history = false , |
||
const bool | log_result = true |
||
) |
Constructor. Provide the reduction factor in addition to arguments that have the same meaning as those of the constructor of the SolverControl constructor.
Definition at line 226 of file solver_control.cc.
ReductionControl::ReductionControl | ( | const SolverControl & | c | ) |
Initialize with a SolverControl object. The result will emulate SolverControl by setting reduce
to zero.
Definition at line 238 of file solver_control.cc.
|
virtualdefault |
Virtual destructor is needed as there are virtual functions in this class.
ReductionControl & ReductionControl::operator= | ( | const SolverControl & | c | ) |
Assign a SolverControl object to ReductionControl. The result of the assignment will emulate SolverControl by setting reduce
to zero.
Definition at line 249 of file solver_control.cc.
|
static |
Interface to parameter file.
Definition at line 296 of file solver_control.cc.
void ReductionControl::parse_parameters | ( | ParameterHandler & | param | ) |
Read parameters from file.
Definition at line 304 of file solver_control.cc.
|
virtual |
Decide about success or failure of an iteration. This function calls the one in the base class, but sets the tolerance to reduction * initial value
upon the first iteration.
Reimplemented from SolverControl.
Definition at line 259 of file solver_control.cc.
double ReductionControl::reduction | ( | ) | const |
Reduction factor.
double ReductionControl::set_reduction | ( | const double | ) |
Change reduction factor.
|
protected |
Desired reduction factor.
Definition at line 466 of file solver_control.h.
|
protected |
Reduced tolerance. Stop iterations if either this value is achieved or if the base class indicates success.
Definition at line 472 of file solver_control.h.