NOX Development
Loading...
Searching...
No Matches
NOX::Direction::Broyden Class Reference

Broyden direction More...

#include <NOX_Direction_Broyden.H>

Inheritance diagram for NOX::Direction::Broyden:
Collaboration diagram for NOX::Direction::Broyden:

Classes

class  BroydenMemoryUnit
 Utility class for NOX::Direction::Broyden::BroydenMemory. More...
class  BroydenMemory
 Utility class for NOX::Direction::Broyden method to manage the information stored in "limited" memory. More...

Public Member Functions

 Broyden (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Constructor.
virtual ~Broyden ()
 Destructor.
virtual bool reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Reset direction based on possibly new parameters.
virtual bool compute (NOX::Abstract::Vector &dir, NOX::Abstract::Group &grp, const NOX::Solver::Generic &solver)
 Not supported for this direction - only works for line search based solver.
virtual bool compute (NOX::Abstract::Vector &dir, NOX::Abstract::Group &grp, const NOX::Solver::LineSearchBased &solver)
 Same as compute(NOX::Abstract::Vector&, NOX::Abstract::Group&, const NOX::Solver::Generic&).
Public Member Functions inherited from NOX::Direction::Generic
 Generic ()
 Constructor.
virtual ~Generic ()
 Destructor.

Detailed Description

Broyden direction

We will calculate a limited-memory Broyden direction of the form

   $d_k = -B_k^{-1} F_k.$

Here $B_k$ is a limited-memory Broyden approximation to the Jacobian of $F$ at $x_k$, and $F_k = F(x_k)$. It is based on apply Broyden updates to the Jacobian from some previous step.

Note
The Broyden direction can only be used with NOX::Solver::LineSearchBased. It cannot be used with any other solver, include NOX::Solver::TrustRegionBased.

References

  • C. T. Kelley, Iterative Methods for Linear and Nonlinear Equations, SIAM, 1995.

Parameters

To use this direction, specify that the "Method" is "Broyden" in the "Direction" sublist of the parameters that are passed to the solver (see NOX::Direction::Manager for more information on choosing the search direction).

In "Direction"/"Broyden":

  • "Restart Frequency" - How often the Jacobian should be refreshed. A value of 5, for example, means that the Jacobian should be updated every 5 iterations. Defaults to 10.

  • "Max Convergence Rate" - Maximum convergence rate allowed when reusing the Jacobian. The Jacobian will be refreshed if the convergence rate, $ \alpha $, is larger than this value. The convergence rate is calculated by $ \alpha = \frac{\| F_k \| }{\|
F_{k-1} \|} $ where F is the nonlinear residual and $ k $ is the nonlinear iteration. Defaults to 1.0.

  • "Memory" - The maximum number of past updates that can be saved in memory. Defaults to the value of "Restart Frequency".

  • "Linear Solver" - optional SUBLIST of linear solver parameters.

  • "Linear Solver"/"Tolerance" - Desired tolerance for linear solve. Defaults to 1.0e-4. The tolerance can be computed using adaptive forcing terms. See NOX::Direction::Utils::InexactNewton for additional options.

Member Function Documentation

◆ compute() [1/2]

bool NOX::Direction::Broyden::compute ( NOX::Abstract::Vector & dir,
NOX::Abstract::Group & grp,
const NOX::Solver::Generic & solver )
virtual

Not supported for this direction - only works for line search based solver.

Implements NOX::Direction::Generic.

◆ compute() [2/2]

bool NOX::Direction::Broyden::compute ( NOX::Abstract::Vector & dir,
NOX::Abstract::Group & grp,
const NOX::Solver::LineSearchBased & solver )
virtual

Same as compute(NOX::Abstract::Vector&, NOX::Abstract::Group&, const NOX::Solver::Generic&).

Enables direct support for line search based solvers for the purpose of efficiency since the LineSearchBased object has a getStep() function that some directions require.

If it is not redefined in the derived class, it will just call the compute with the NOX::Solver::Generic argument.

Add this direction to the memory

Reimplemented from NOX::Direction::Generic.

References NOX::Abstract::Group::clone(), NOX::Abstract::Group::computeF(), NOX::DeepCopy, NOX::Utils::Details, NOX::Abstract::Group::getF(), NOX::Solver::LineSearchBased::getNumIterations(), NOX::Solver::LineSearchBased::getPreviousSolutionGroup(), NOX::Abstract::Group::Ok, NOX::Abstract::Vector::scale(), and NOX::Abstract::Vector::update().

◆ reset()

bool NOX::Direction::Broyden::reset ( const Teuchos::RCP< NOX::GlobalData > & gd,
Teuchos::ParameterList & params )
virtual

Reset direction based on possibly new parameters.

Implements NOX::Direction::Generic.

Referenced by Broyden().


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