OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
FuelCell::ApplicationCore::Newton3pp Class Reference

Application class performing Newton's iteration. More...

#include <newton_w_3pp.h>

Inheritance diagram for FuelCell::ApplicationCore::Newton3pp:
Inheritance graph
[legend]
Collaboration diagram for FuelCell::ApplicationCore::Newton3pp:
Collaboration graph
[legend]

Public Member Functions

 Newton3pp (ApplicationBase &app)
 Constructor. More...
 
void declare_parameters (ParameterHandler &param)
 Declare any additional parameters necessary for the solver. More...
 
void initialize (ParameterHandler &param)
 Initialize NewtonBase and any additional parameters. More...
 
virtual void solve (FuelCell::ApplicationCore::FEVector &u, const FuelCell::ApplicationCore::FEVectors &in_vectors)
 The actual Newton solver. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::newtonBase
 newtonBase (ApplicationBase &app)
 Constructor, receiving the application computing the residual and solving the linear problem. More...
 
void _initialize (ParameterHandler &param)
 Read the parameters local to newtonBase. More...
 
double threshold (double new_value)
 Set the maximal residual reduction allowed without triggering assembling in the next step. More...
 
void initialize_initial_guess (BlockVector< double > &dst)
 Control object for the Newton iteration. More...
 
virtual void assemble ()
 Instead of assembling, this function only sets a flag, such that the inner application will be required to assemble a new derivative matrix next time solve() is called. More...
 
virtual double residual (FuelCell::ApplicationCore::FEVector &dst, const FuelCell::ApplicationCore::FEVectors &rhs)
 Returns the L2-norm of the residual and the residual vector in "dst" using the residual function in the ApplicationBase used to initialize the application. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::ApplicationWrapper
 ApplicationWrapper (ApplicationBase &app)
 Constructor for a derived application. More...
 
 ~ApplicationWrapper ()
 Destructor. More...
 
virtual void remesh ()
 Generate the next mesh depending on the mesh generation parameters. More...
 
virtual void init_vector (FEVector &dst) const
 Initialize vector to problem size. More...
 
virtual double residual (FEVector &dst, const FEVectors &src, bool apply_boundaries=true)
 Compute residual of src and store it into dst. More...
 
virtual void Tsolve (FEVector &dst, const FEVectors &src)
 Solve the dual system assembled with right hand side rhs and return the result in start. More...
 
virtual double estimate (const FEVectors &src)
 Estimate cell-wise errors. More...
 
virtual double evaluate (const FEVectors &src)
 Evaluate a functional. More...
 
virtual void grid_out (const std::string &filename) const
 
virtual void data_out (const std::string &filename, const FEVectors &src)
 Write data in the format specified by the ParameterHandler. More...
 
virtual std::string id () const
 Return a unique identification string for this application. More...
 
virtual void notify (const Event &reason)
 Add a reason for assembling. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::ApplicationBase
 ApplicationBase (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >())
 Constructor for an application. More...
 
 ApplicationBase (const ApplicationBase &other)
 Copy constructor. More...
 
virtual ~ApplicationBase ()
 Virtual destructor. More...
 
void print_parameters_to_file (ParameterHandler &param, const std::string &file_name, const ParameterHandler::OutputStyle &style)
 Print default parameters for the application to a file. More...
 
virtual void start_vector (FEVector &dst, std::string) const
 Initialize vector to problem size. More...
 
virtual void grid_out (const std::string &)
 Write the mesh in the format specified by the ParameterHandler. More...
 
boost::shared_ptr
< ApplicationData
get_data ()
 Get access to the protected variable data. More...
 
const boost::shared_ptr
< ApplicationData
get_data () const
 Get read-only access to the protected variable data. More...
 
virtual void clear ()
 All true in notifications. More...
 
virtual void clear_events ()
 All false in notifications. More...
 
virtual unsigned int get_solution_index ()
 Returns solution index. More...
 

Static Public Attributes

static const
FuelCell::ApplicationCore::Event 
bad_derivative
 
- Static Public Attributes inherited from FuelCell::ApplicationCore::newtonBase
static const
FuelCell::ApplicationCore::Event 
bad_derivative
 The Event set if convergence is becoming bad and a new matrix should be assembled. More...
 

Private Member Functions

void set_predefined_parameters ()
 Define default parameters for convergence. More...
 
double sum_of_squares (double residual_norm)
 Sum of squares. More...
 
double three_point_step (int line_search_iterations, double lambda_current, double lambda_previous, double mf_original, double mf_current, double mf_previous)
 Merit Function. More...
 
bool armijo (double step_length, double residual_norm, double new_residual_norm)
 Determines step size by three point parabolic function. More...
 

Private Attributes

double max_reduction_factor
 Convergance criterion. More...
 
double min_reduction_factor
 min_reduction_factor corresponds to $ \sigma_0 $ More...
 
double alpha
 This parameter is used in order to assess when the line search should be stopped. More...
 
bool with_prediction
 Boolean variable that specifies if the three-point method is used with (true) or without (false) prediction, i.e. More...
 
bool use_predefined
 Flag so that predefined parameters are used. More...
 

Additional Inherited Members

- Public Attributes inherited from FuelCell::ApplicationCore::newtonBase
ReductionControl control
 Control object for the Newton iteration. More...
 
double numIter
 Number of Iterations;. More...
 
- Protected Member Functions inherited from FuelCell::ApplicationCore::newtonBase
void debug_output (const FEVector &sol, const FEVector &update, const FEVector &residual) const
 Function used to output any necessary information at each Newton iteration for debugging purposes. More...
 
- Protected Member Functions inherited from FuelCell::ApplicationCore::ApplicationWrapper
SmartPointer< ApplicationBaseget_wrapped_application ()
 Gain access to the inner application. More...
 
- Protected Member Functions inherited from FuelCell::ApplicationCore::ApplicationBase
void print_caller_name (const std::string &caller_name) const
 Print caller name. More...
 
- Protected Attributes inherited from FuelCell::ApplicationCore::newtonBase
bool assemble_now
 This flag is set by the function assemble(), indicating that the matrix must be assembled anew upon start. More...
 
double assemble_threshold
 Threshold for re-assembling matrix. More...
 
bool debug_solution
 Print updated solution after each step into file Newton_uNNN? More...
 
bool debug_update
 Print Newton update after each step into file Newton_dNNN? More...
 
bool debug_residual
 Print Newton residual after each step into file Newton_rNNN? More...
 
unsigned int debug
 Write debug output to FcstUtilities::log; the higher the number, the more output. More...
 
unsigned int step
 The number of a basic Newton iteration. More...
 
std::vector< unsigned int > blocks
 This vector specifies the blocks of the global solution which are supposed to be treated specially. More...
 
unsigned int n_blocks
 The total number of blocks. More...
 
- Protected Attributes inherited from FuelCell::ApplicationCore::ApplicationWrapper
SmartPointer< ApplicationBaseapp
 Pointer to the application this one depends upon. More...
 
- Protected Attributes inherited from FuelCell::ApplicationCore::ApplicationBase
boost::shared_ptr
< ApplicationData
data
 Object for auxiliary data. More...
 
Event notifications
 Accumulate reasons for assembling here. More...
 

Detailed Description

Application class performing Newton's iteration.

At each iteration a line search is performed using a three point parabolic algorithm to enhance robustness of the algorithm and to improve convergence to a global solution. Either a standard three point parabolic algorithm or a three point parabolic algorithm with prediction can be used depending on the input parameter "set Include step size prediction" in the Newton section of the input file. If this option is set to true, the algorithm with prediction is used. Othewise, the algorithm without prediction is used.

The three-point parabolic model, described in section 8.3.1 of [1], attempts to find a step size length, $ h $, such that it minimizes the scalar function

\[ f (h) = \| \mathbf{F}(\mathbf{x}_{i−1} + h \delta \mathbf{x} ) \|_2^2. \]

where the meaning of $ \delta \mathbf{x} $ and how it is obtain is described in the documentation for the base class, i.e. newtonBase.

We estimate the norm above using the following three-point interpolating polynomial

\[ p(h) = f (h) = f(0) + \frac{h}{h_{j-1} − h_{j−2}} \left( \frac{(h − h_{j−2}) ( f(h_{j−1}) − f(0) )}{h_{j−1}} + \frac{(h_{j−1} − h)(f(h_{j−2})-f(0)}{h_{j−2}} \right ) \]

where 0, $ h_{j−1} $, $ h_{j−2} $ are the initial and the two lastly rejected step sizes in the interpolation.

Using the equation above, and assuming that $ p > 0 $, the step size $ h_j $ is given by

\[ h_j = − \frac{p'(0)}{p''(0)} \]

If $ p'' ≤ 0 $, $ h_j $ is reduced by a fixed amount.

The three-point parabolic method is continued until a sufficient decrease in the residual $ \| \mathbf{F}(x_i) \| $ occurs based on the condition

\[ \| \mathbf{F}(\mathbf{x_{i−1}} + h_j \delta \mathbf{x}) \| < (1 − \alpha h_j ) \| \mathbf{F}(x_{i−1}) \|. \]

where $ \alpha ∈ [0, 1] $. In practice, $ \alpha $ is typically $ 10^{−4}$ [2].

When implementing this method, there are a few considerations to take into account.

To prevent a reduction in the step size that is either too large or too small to be effective, we allow the user to specify a maximum and minimum reduction rate, i.e. $ \sigma_0 , \sigma_1 $. Therefore, $ h_j $ must equal $ \sigma h_{j−1} $, where $ \sigma ∈ [ \sigma_0 , \sigma_1] $. This is guaranteed by making the following check after $ h_j $ is determined

\[ h_j = \left\{ \begin{array}{ll} \sigma_0 h_{j−1} & \mbox{ if } h_j < \sigma_0 h_{j−1} \\ \sigma_1 h_{j−1} & \mbox{ if } h_j > \sigma_1 h_{j−1} \\ h_j & \mbox{otherwise}. \end{array} \right. \]

The three-point interpolating polynomial also requires some initial values for $ h_0 $ and $ h_1 $. If a maximum reduction factor is set, then it is convenient to choose

\[ h_0 = 1, h_1 = \sigma_1 \]

Three point parabolic with prediction

The value of $ h_0 $ can be seen as an initial guess for the three-point parabolic method. In the previous section, we simply set $ h_0 = 1 $. Of course, it is very likely that 1 is a poor choice for the initial guess. This results in additional line-search iterations.

In an effort to reduce the number of line-search iterations, we can use the previous step size to more accurately choose an initial step size $ h_{0,i} $ for the current Newton iteration $ i $. However, there is also a possibility that the previous step size is a poor choice for an initial guess, therefore we use the condition

\[ h_{0,i} = \left\{ \begin{array}{ll} h_{i−1} & \mbox{ if } h_{i-1} < h_{i−2}(1-\alpha) \\ min(1, 2 h_{i-1}) & \mbox{otherwise}. \end{array} \right. \]

The step sizes $ h_{i−1}$ and $ h_{i−2} $ are from previous two Newton iterations. Initially, $ h_{0,1} = h_{0,2} = 1 $.

Using this strategy results in a predictor-corrector approach to the three-point parabolic method. The condition above predicts the value of the step size. After, the three-point parabolic method simply corrects for the differences between $ h_{i−1} $ and $ h_i $.

Note
This is the method by default (set Include step size prediction = true)

Parameters

The parameters that control all Newton solvers are defined in the subsection Newton in the data input file. The parameters are the following:

* subsection Newton
* set Include step size prediction = true # Use prediction with the three-point parabolic method
* set Use predefined parameters = false # Use optimized reduction factor and line search parameters.
* set Max. reduction factor = 0.9
* set Min. reduction factor = 0.5
* set Line search convergence, alpha [0-1] = 0.01
* // All the parameters below are from NetwonBase:
* set Max steps = 100 # Maximum number of iterations
* set Tolerance = 1.e-8 # Absolute tolerance
* set Reduction = 1.e-20 # Maximum allowed reduction during line search. If the residual is not reduced by at least this number the iteration is discarded.
* set Assemble threshold = 0.0 # Used as a threshold to state if the Jacobian needs to be re-evaluated
* set Debug level = 0
* set Debug residual = false # Would you like the code to output the residual at every Newton iteration?
* set Debug solution = false # Would you like the code to output the solution at every Newton iteration?
* set Debug update = false
* end
*

References

[1] C. T. Kelley. Iterative methods for linear and nonlinear equations, volume 16 of Frontiers in Applied Mathematics. Society for Industrial and Applied Mathematics (SIAM), Philadelphia, PA, 1995. With separately available software.

[2] C. T. Kelley. Solving nonlinear equations with Newton’s method. Fundamentals of Algorithms. Society for Industrial and Applied Mathematics (SIAM), Philadelphia, PA, 2003.

Author
Jason Boisvert, Raymond Spiteri and Marc Secanell, 2009-13

Constructor & Destructor Documentation

FuelCell::ApplicationCore::Newton3pp::Newton3pp ( ApplicationBase app)

Constructor.

Member Function Documentation

bool FuelCell::ApplicationCore::Newton3pp::armijo ( double  step_length,
double  residual_norm,
double  new_residual_norm 
)
private

Determines step size by three point parabolic function.

void FuelCell::ApplicationCore::Newton3pp::declare_parameters ( ParameterHandler &  param)
virtual

Declare any additional parameters necessary for the solver.

In this case the only additional parameter is a bool parameter estating if we want to use a Newton step prediction, i.e.

In subsection Newton you have the following option:

  • set Include step size prediction = true (false)

If true then a three-point parabolic method with prediction is used. Otherwise, the standard three-point parabolic is used.

Reimplemented from FuelCell::ApplicationCore::newtonBase.

void FuelCell::ApplicationCore::Newton3pp::initialize ( ParameterHandler &  param)
virtual

Initialize NewtonBase and any additional parameters.

Reimplemented from FuelCell::ApplicationCore::newtonBase.

void FuelCell::ApplicationCore::Newton3pp::set_predefined_parameters ( )
inlineprivate

Define default parameters for convergence.

These depend on refinement level.

References alpha, FuelCell::ApplicationCore::ApplicationBase::get_data(), max_reduction_factor, and min_reduction_factor.

Here is the call graph for this function:

virtual void FuelCell::ApplicationCore::Newton3pp::solve ( FuelCell::ApplicationCore::FEVector u,
const FuelCell::ApplicationCore::FEVectors in_vectors 
)
virtual

The actual Newton solver.

In this section, the Newton solver performs the iterations necessary to solve the nonlinear system. During the process it will call assemble() and residual() as needed. Once it is done, it will return the solution as parameter {u}.

Implements FuelCell::ApplicationCore::newtonBase.

double FuelCell::ApplicationCore::Newton3pp::sum_of_squares ( double  residual_norm)
private

Sum of squares.

double FuelCell::ApplicationCore::Newton3pp::three_point_step ( int  line_search_iterations,
double  lambda_current,
double  lambda_previous,
double  mf_original,
double  mf_current,
double  mf_previous 
)
private

Merit Function.

Member Data Documentation

double FuelCell::ApplicationCore::Newton3pp::alpha
private

This parameter is used in order to assess when the line search should be stopped.

This value should be between zero and one. In practice, $ \alpha $ is typically $ 10^{−4}$.

The larger the value of $ \alpha $, the strictest the convergence criteria, i.e.

\[ \| \mathbf{F}(\mathbf{x_{i−1}} + h_j \delta \mathbf{x}) \| < (1 − \alpha h_j ) \| \mathbf{F}(x_{i−1}) \|. \]

will be.

This value can be modified in the input file in subsection Newton using "set alpha = 1e-4".

Referenced by set_predefined_parameters().

const FuelCell::ApplicationCore::Event FuelCell::ApplicationCore::Newton3pp::bad_derivative
static
double FuelCell::ApplicationCore::Newton3pp::max_reduction_factor
private

Convergance criterion.

To prevent a reduction in the step size that is either too large or too small to be effective, we allow the user to specify a maximum and minimum reduction rate, i.e. $ \sigma_0 , \sigma_1 $. Therefore, $ h_j $ must equal $ \sigma h_{j−1} $, where $ \sigma ∈ [ \sigma_0 , \sigma_1] $. This is guaranteed by making the following check after $ h_j $ is determined

\[ h_j = \left\{ \begin{array}{ll} \sigma_0 h_{j−1} & \mbox{ if } h_j < \sigma_0 h_{j−1} \\ \sigma_1 h_{j−1} & \mbox{ if } h_j > \sigma_1 h_{j−1} \\ h_j & \mbox{otherwise}. \end{array} \right. \]

max_reduction_factor corresponds to $ \sigma_1 $

Referenced by set_predefined_parameters().

double FuelCell::ApplicationCore::Newton3pp::min_reduction_factor
private

min_reduction_factor corresponds to $ \sigma_0 $

Referenced by set_predefined_parameters().

bool FuelCell::ApplicationCore::Newton3pp::use_predefined
private

Flag so that predefined parameters are used.

bool FuelCell::ApplicationCore::Newton3pp::with_prediction
private

Boolean variable that specifies if the three-point method is used with (true) or without (false) prediction, i.e.

the first iteration step is based on previous iteration steps.


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