Hermes2D  2.0
Hermes::Hermes2D::NewtonSolver< Scalar > Class Template Reference

#include <newton_solver.h>

+ Inheritance diagram for Hermes::Hermes2D::NewtonSolver< Scalar >:

Public Member Functions

 NewtonSolver (DiscreteProblem< Scalar > *dp)
 
 NewtonSolver (const WeakForm< Scalar > *wf, const Space< Scalar > *space)
 
 NewtonSolver (const WeakForm< Scalar > *wf, Hermes::vector< const Space< Scalar > * > spaces)
 
void init_linear_solver ()
 
virtual bool isOkay () const
 State querying helpers.
 
std::string getClassName () const
 Get class name, for the purpose of messaging.
 
void solve (Scalar *coeff_vec=NULL)
 
void solve (Solution< Scalar > *initial_guess)
 
void solve (Hermes::vector< Solution< Scalar > * > initial_guess)
 
void solve_keep_jacobian (Scalar *coeff_vec=NULL)
 
void solve_keep_jacobian (Solution< Scalar > *initial_guess)
 
void solve_keep_jacobian (Hermes::vector< Solution< Scalar > * > initial_guess)
 
void set_max_allowed_residual_norm (double max_allowed_residual_norm_to_set)
 
void set_min_allowed_damping_coeff (double min_allowed_damping_coeff_to_set)
 
virtual void set_iterative_method (const char *iterative_method_name)
 Call NonlinearSolver::set_iterative_method() and set the method to the linear solver (if applicable).
 
virtual void set_preconditioner (const char *preconditioner_name)
 Call NonlinearSolver::set_preconditioner() and set the method to the linear solver (if applicable).
 
void set_residual_as_function ()
 
void set_newton_tol (double newton_tol)
 
void set_newton_max_iter (int newton_max_iter)
 
virtual void set_time (double time)
 
virtual void set_time_step (double time_step)
 
virtual void set_spaces (Hermes::vector< const Space< Scalar > * > spaces)
 See the class Hermes::Hermes2D::Mixins::SettableSpaces.
 
virtual void set_space (const Space< Scalar > *space)
 
virtual Hermes::vector< const
Space< Scalar > * > 
get_spaces () const
 Get all spaces as a Hermes::vector.
 
void set_manual_damping_coeff (bool onOff, double coeff=1.0)
 
void set_initial_auto_damping_coeff (double coeff)
 
void set_auto_damping_ratio (double ratio)
 
void set_sufficient_improvement_factor (double ratio)
 
void set_necessary_successful_steps_to_increase (unsigned int steps)
 
void set_weak_formulation (const WeakForm< Scalar > *wf)
 Set the weak forms.
 
- Public Member Functions inherited from Hermes::Hermes2D::Mixins::SettableSpaces< Scalar >
virtual const Space< Scalar > * get_space (int n) const
 
- Public Member Functions inherited from Hermes::Hermes2D::Mixins::MatrixRhsOutput< Scalar >
 MatrixRhsOutput ()
 
void output_matrix (int firstIterations=-1)
 
void set_print_zero_matrix_entries (bool to_set)
 Sets this instance to output matrix entries even though they are zero or not.
 
void set_matrix_filename (std::string name)
 
void set_matrix_varname (std::string name)
 
void set_matrix_E_matrix_dump_format (EMatrixDumpFormat format)
 
void set_matrix_number_format (char *number_format)
 
void output_rhs (int firstIterations=-1)
 
void set_rhs_filename (std::string name)
 
void set_rhs_varname (std::string name)
 
void set_rhs_E_matrix_dump_format (EMatrixDumpFormat format)
 
void set_rhs_number_format (char *number_format)
 
- Public Member Functions inherited from Hermes::Hermes2D::Mixins::StateQueryable
void check () const
 Method to handle the state.
 

Protected Member Functions

void init_attributes ()
 Internal setting of default values (see individual set methods).
 

Protected Attributes

const bool own_dp
 This instance owns its DP.
 
SparseMatrix< Scalar > * kept_jacobian
 Used by method solve_keep_jacobian().
 
SparseMatrix< Scalar > * jacobian
 Jacobian.
 
Vector< Scalar > * residual
 Residual.
 
LinearMatrixSolver< Scalar > * linear_solver
 Linear solver.
 
double newton_tol
 
int newton_max_iter
 
bool residual_as_function
 
double max_allowed_residual_norm
 
double min_allowed_damping_coeff
 
double currentDampingCofficient
 
bool manual_damping
 Manual / auto.
 
double auto_damping_ratio
 The ratio between two damping coeffs when changing.
 
double initial_auto_damping_ratio
 The initial (and maximum) damping coefficient.
 
double sufficient_improvement_factor
 Sufficient improvement for continuing.
 
unsigned int necessary_successful_steps_to_increase
 necessary number of steps to increase back the damping coeff.
 
- Protected Attributes inherited from Hermes::Hermes2D::Mixins::MatrixRhsOutput< Scalar >
bool print_matrix_zero_values
 
bool output_matrixOn
 
int output_matrixIterations
 
std::string matrixFilename
 
std::string matrixVarname
 
EMatrixDumpFormat matrixFormat
 
char * matrix_number_format
 
bool output_rhsOn
 
int output_rhsIterations
 
std::string RhsFilename
 
std::string RhsVarname
 
EMatrixDumpFormat RhsFormat
 
char * rhs_number_format
 

Detailed Description

template<typename Scalar>
class Hermes::Hermes2D::NewtonSolver< Scalar >

Class for Newton's method.
Typical usage:
// Initialize Newton's solver.
// Here wf is Hermes2D::WeakForm<double>, space is Hermes2D::Space<double>
Hermes::Hermes2D::NewtonSolver<double> newton_solver(&wf, &space);
Set a whole bunch of parameters according to your liking.
See the class documentation for all possible parameters.
newton_solver.set_newton_tol(1e-6);
newton_solver.set_newton_max_iter(15);
newton_solver.set_max_allowed_residual_norm(1e6);
newton_solver.set_min_allowed_damping_coeff(1e-3);

// Solve the linear problem.
try
{
 // Just call solve().
 newton_solver.solve();

 // Get the solution vector from the solver.
 double* sln_vector = newton_solver.get_sln_vector();

 // Translate the solution vector into the previously initialized Solution<double> using the static method vector_to_solution.
 Hermes::Hermes2D::Solution<double>::vector_to_solution(sln_vector, &space, &sln);
}
// All kinds of Exceptions may happen (Linear algebraic solver, some bad parameters, some data not initialized...)
catch(Hermes::Exceptions::Exception& e)
{
 e.print_msg();
 return -1;
}
// For illustrative purposes, otherwise one can just catch std::exception directly, as Hermes::Exceptions::Exception derive from it.
catch(std::exception& e)
{
 std::cout << e.what();
 return -1;
}

Definition at line 71 of file newton_solver.h.

Member Function Documentation

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_auto_damping_ratio ( double  ratio)

Set the ratio to the automatic damping. When the damping coefficient is decided to be descreased or increased, this is the ratio how it will be changed (this is the bigger ( > 1.0 ) of the two possible values). I.e. when the damping coefficient is shortened 3 times if deemed too big, make the parameter not 0.333333, but 3.0. Default: 2.0

Parameters
[in]ratioThe ratio (again, it must be > 1.0, and it represents the inverse of the shortening factor).

Definition at line 213 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_initial_auto_damping_coeff ( double  coeff)

Make the automatic damping start with this coefficient. This will also be the top bound for the coefficient. Default: 1.0

Parameters
[in]coeffThe initial damping coefficient. Must be > 0 and <= 1.0.

Definition at line 202 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_manual_damping_coeff ( bool  onOff,
double  coeff = 1.0 
)

Turn on or off manual damping (default is the automatic) and optionally sets manual damping coefficient. Default: default is the automatic damping, default coefficient if manual damping used is set by this method.

Parameters
[in]onOffon(true)-manual damping, off(false)-automatic damping.
[in]coeffThe (perpetual) damping coefficient in the case of manual damping. Ignored in the case of automatic damping.

Definition at line 188 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_max_allowed_residual_norm ( double  max_allowed_residual_norm_to_set)

Sets the maximum allowed norm of the residual during the calculation. Default: 1E9

Definition at line 102 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_min_allowed_damping_coeff ( double  min_allowed_damping_coeff_to_set)

Sets minimum damping coefficient. Default: 1E-4

Definition at line 110 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_necessary_successful_steps_to_increase ( unsigned int  steps)

Set how many successful steps are necessary for the damping coefficient to be increased, by multiplication by the parameter set by set_auto_damping_ratio(). The coefficient is then increased after each 'successful' step, if the sequence of such is not interrupted by an 'unsuccessful' step. Default: 1

Parameters
[in]stepsNumber of steps.

Definition at line 233 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_newton_max_iter ( int  newton_max_iter)

Set the maximum number of Newton's iterations. Default: 15

Definition at line 94 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_newton_tol ( double  newton_tol)

Set the residual norm tolerance for ending the Newton's loop. Default: 1E-8.

Definition at line 82 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_residual_as_function ( )

Interpret the residual as a function. Translate the residual vector into a residual function (or multiple functions) in the corresponding finite element space(s) and measure their norm(s) there. This is more meaningful than just measuring the l2-norm of the residual vector, since in the FE space not all components in the residual vector have the same weight. On the other hand, this is slower as it requires global norm calculation, and thus numerical integration over the entire domain. Therefore this option is off by default.

Definition at line 118 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_sufficient_improvement_factor ( double  ratio)

Set the ratio of the current residual norm and the previous residual norm necessary to deem a step 'successful'. It can be either > 1.0, meaning that even if the norm increased, the step will be 'successful', or < 1.0, meaning that even though the residual norm goes down, we will further decrease the damping coefficient. Default: 0.95 param[in] ratio The ratio, must be positive.

Definition at line 223 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::set_time ( double  time)
virtual

Set time information for time-dependent problems. See the class Hermes::Mixins::TimeMeasurable.

Definition at line 124 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve ( Scalar *  coeff_vec = NULL)

Solve.

Parameters
[in]coeff_vecCeofficient vector to start from.

Definition at line 261 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve ( Solution< Scalar > *  initial_guess)

Solve.

Parameters
[in]initial_guessSolution to start from (which is projected to obtain the initial coefficient vector.

Definition at line 243 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve ( Hermes::vector< Solution< Scalar > * >  initial_guess)

Solve.

Parameters
[in]initial_guessSolutions to start from (which is projected to obtain the initial coefficient vector.

Definition at line 251 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve_keep_jacobian ( Scalar *  coeff_vec = NULL)

Solve which keeps jacobian. A solve() method where the jacobian is reused.

Definition at line 524 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve_keep_jacobian ( Solution< Scalar > *  initial_guess)

Solve which keeps jacobian.

Parameters
[in]initial_guessSolution to start from (which is projected to obtain the initial coefficient vector.

Definition at line 506 of file newton_solver.cpp.

template<typename Scalar >
void Hermes::Hermes2D::NewtonSolver< Scalar >::solve_keep_jacobian ( Hermes::vector< Solution< Scalar > * >  initial_guess)

Solve which keeps jacobian.

Parameters
[in]initial_guessSolutions to start from (which is projected to obtain the initial coefficient vector.

Definition at line 514 of file newton_solver.cpp.

Member Data Documentation

template<typename Scalar >
double Hermes::Hermes2D::NewtonSolver< Scalar >::max_allowed_residual_norm
protected

Maximum allowed residual norm. If this number is exceeded, the methods solve() return 'false'. By default set to 1E6. Possible to change via method set_max_allowed_residual_norm().

Definition at line 214 of file newton_solver.h.


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