NOX Development
|
Broyden direction More...
#include <NOX_Direction_Broyden.H>
Classes | |
class | BroydenMemory |
Utility class for NOX::Direction::Broyden method to manage the information stored in "limited" memory. More... | |
class | BroydenMemoryUnit |
Utility class for NOX::Direction::Broyden::BroydenMemory. More... | |
Public Member Functions | |
Broyden (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms) | |
Constructor. | |
virtual | ~Broyden () |
Destructor. | |
virtual bool | reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms) |
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&) | |
![]() | |
Generic () | |
Constructor. | |
virtual | ~Generic () |
Destructor. | |
virtual bool | reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms)=0 |
Reset direction based on possibly new parameters. | |
virtual bool | compute (NOX::Abstract::Vector &dir, NOX::Abstract::Group &grp, const NOX::Solver::Generic &solver)=0 |
Compute the direction vector, dir , for a specific method given the current group, grp . | |
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&) | |
Broyden direction
We will calculate a limited-memory Broyden direction of the form
Here
References
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,
"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.
|
virtual |
Not supported for this direction - only works for line search based solver.
Implements NOX::Direction::Generic.
|
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::Abstract::Group::computeJacobian(), NOX::DeepCopy, NOX::Utils::Details, NOX::Abstract::Group::getF(), NOX::Solver::LineSearchBased::getNumIterations(), NOX::Solver::LineSearchBased::getPreviousSolutionGroup(), NOX::Abstract::Vector::innerProduct(), NOX::Abstract::Group::Ok, NOX::Abstract::Vector::scale(), and NOX::Abstract::Vector::update().
|
virtual |
Reset direction based on possibly new parameters.
Implements NOX::Direction::Generic.
Referenced by Broyden().