NOX
Development
|
Nonlinear solver based on a rank-1 tensor method. More...
#include <NOX_Solver_TensorBased.H>
Public Member Functions | |
TensorBased (const Teuchos::RCP< NOX::Abstract::Group > &grp, const Teuchos::RCP< NOX::StatusTest::Generic > &tests, const Teuchos::RCP< Teuchos::ParameterList > ¶ms) | |
Constructor. More... | |
virtual | ~TensorBased () |
Destructor. | |
virtual void | reset (const NOX::Abstract::Vector &initialGuess, const Teuchos::RCP< NOX::StatusTest::Generic > &tests) |
Resets the solver, sets a new status test, and sets a new initial guess. | |
virtual void | reset (const NOX::Abstract::Vector &initialGuess) |
Resets the solver and sets a new initial guess. | |
virtual void | reset () |
Resets the solver for another solve. This resets the counters and status only. Uses the final solution from the last solve as the initial guess for the next solve. More... | |
virtual NOX::StatusTest::StatusType | step () |
Do one nonlinear step in the iteration sequence and return status. | |
virtual NOX::StatusTest::StatusType | solve () |
Solve the nonlinear problem and return final status. More... | |
virtual const NOX::Abstract::Group & | getSolutionGroup () const |
Return a reference to the current solution group. | |
virtual const NOX::Abstract::Group & | getPreviousSolutionGroup () const |
Return a reference to the previous solution group. | |
virtual NOX::StatusTest::StatusType | getStatus () const |
Returns the current status of the solver. | |
virtual int | getNumIterations () const |
Get number of iterations. | |
virtual const Teuchos::ParameterList & | getList () const |
Return a reference to the solver parameters. | |
virtual Teuchos::RCP< const NOX::Abstract::Group > | getSolutionGroupPtr () const |
Return a RCP to the solution group. | |
virtual Teuchos::RCP< const NOX::Abstract::Group > | getPreviousSolutionGroupPtr () const |
Return a RCP to the previous solution group. | |
virtual Teuchos::RCP< const Teuchos::ParameterList > | getListPtr () const |
Return a RCP to the solver parameters. | |
virtual Teuchos::RCP< const NOX::SolverStats > | getSolverStatistics () const |
Return a RCP to the solver statistics. | |
Public Member Functions inherited from NOX::Solver::Generic | |
Generic () | |
Constructor (does nothing) | |
virtual | ~Generic () |
Destructor (does nothing) | |
Protected Types | |
enum | StepType { TensorStep, NewtonStep } |
Types of steps. | |
enum | LineSearchType { Curvilinear, Standard, Dual, FullStep, Newton } |
Flag for the direction to be computed this iteration. More... | |
enum | ConvergenceCriteriaType { ArmijoGoldstein, AredPred, None } |
Algorithms used to determine convergence of the line search. More... | |
enum | LambdaSelectionType { Halving, Quadratic } |
Types of lambda selection. | |
enum | RecoveryStepType { Constant, LastComputedStep } |
Type of recovery step to use. More... | |
Protected Member Functions | |
virtual void | init () |
Print out initialization information. | |
virtual void | printUpdate () |
Prints the current iteration information. | |
virtual bool | reset (const Teuchos::RCP< NOX::Abstract::Group > &grp, const Teuchos::RCP< NOX::StatusTest::Generic > &tests, const Teuchos::RCP< Teuchos::ParameterList > ¶ms) |
Constructor initialization routine. | |
bool | computeTensorDirection (NOX::Abstract::Group &soln, const NOX::Solver::Generic &solver) |
Subroutine for computing the tensor and Newton directions. | |
double | calculateBeta (double qa, double qb, double qc, double &qval, double &lambdaBar, double lambda=1.0) const |
Subroutine for calculating beta. | |
bool | computeCurvilinearStep (NOX::Abstract::Vector &dir, const NOX::Abstract::Group &soln, const NOX::Solver::Generic &s, double &lambda) |
Subroutine for computing the curvilinear step. | |
bool | implementGlobalStrategy (NOX::Abstract::Group &newGrp, double &step, const NOX::Solver::Generic &s) |
Subroutine for executing the tensor linesearch. | |
bool | performLinesearch (NOX::Abstract::Group &newsoln, double &step, const NOX::Abstract::Vector &lsDir, const NOX::Solver::Generic &s) |
Performs a standard tensor linesearch (tensor or Newton direction) | |
double | getNormModelResidual (const NOX::Abstract::Vector &dir, const NOX::Abstract::Group &soln, bool isTensorModel) const |
Compute the residual norm of the local model. | |
void | printDirectionInfo (std::string dirName, const NOX::Abstract::Vector &dir, const NOX::Abstract::Group &soln, bool isTensorModel) const |
Print pertinent information about the direction. | |
double | getDirectionalDerivative (const NOX::Abstract::Vector &dir, const NOX::Abstract::Group &soln) const |
Calculate the directional derivative. | |
double | selectLambda (double newf, double oldf, double oldfprime, double lambda) |
Select lambda for linesearch (quadratic or halving) | |
void | throwError (const std::string &functionName, const std::string &errorMsg) const |
Throw an error with a method's name and error message. | |
Protected Attributes | |
Teuchos::RCP< NOX::GlobalData > | globalDataPtr |
Pointer to the global data object. | |
Teuchos::RCP< NOX::Utils > | utilsPtr |
Printing Utils. | |
Teuchos::RCP < NOX::Abstract::Group > | solnPtr |
Current solution. | |
Teuchos::RCP < NOX::Abstract::Group > | oldSolnPtr |
Previous solution pointer. More... | |
Teuchos::RCP < NOX::Abstract::Vector > | newtonVecPtr |
Current Newton direction (pointer). More... | |
Teuchos::RCP < NOX::Abstract::Vector > | tensorVecPtr |
Current tensor direction (pointer). More... | |
Teuchos::RCP < NOX::Abstract::Vector > | aVecPtr |
Current tensor term vector (pointer). More... | |
Teuchos::RCP < NOX::Abstract::Vector > | sVecPtr |
Vector to previous point (pointer). More... | |
Teuchos::RCP < NOX::Abstract::Vector > | tmpVecPtr |
Working vector (pointer). More... | |
Teuchos::RCP < NOX::Abstract::Vector > | residualVecPtr |
Residual vector (pointer). More... | |
Teuchos::RCP < NOX::StatusTest::Generic > | testPtr |
Stopping test. | |
NOX::StatusTest::CheckType | checkType |
Type of check to use for status tests. See NOX::StatusTest for more details. | |
Teuchos::RCP < Teuchos::ParameterList > | paramsPtr |
Input parameters. | |
Teuchos::ParameterList * | linearParamsPtr |
Line Search parameters. More... | |
double | stepSize |
Current step. | |
double | beta |
Value of sc'*dt. | |
int | nIter |
Number of nonlinear iterations. | |
NOX::StatusTest::StatusType | status |
Status of nonlinear solver. | |
StepType | requestedBaseStep |
Flag for the base direction to compute after the first iteration. | |
LineSearchType | lsType |
Choice of line search. | |
ConvergenceCriteriaType | convCriteria |
Choice of convergence criteria (currently unused) | |
LambdaSelectionType | lambdaSelection |
Flag for lambda selection (Halving/Quadratic) | |
RecoveryStepType | recoveryStepType |
Choice of the recovery step type; uses "Recovery Step Type" parameter. | |
bool | useModifiedMethod |
Flag for using modifications that force quadratic to have real root. | |
bool | isNewtonDirection |
Flag for Newton direction. | |
bool | doRescue |
Flag for rescuing Linear Solver from a bad solve. | |
double | minStep |
Minimum step length (i.e., when we give up) | |
double | defaultStep |
Default step. | |
double | recoveryStep |
Default step for linesearch failure. | |
int | maxIters |
Maximum iterations. | |
double | alpha |
Scaling factor for the Armijo-Goldstein condition. | |
double | sTinvJF |
Value of s'*inv(J)*F. | |
double | sTinvJa |
Value of s'*inv(J)*a. | |
Teuchos::RCP < NOX::LineSearch::Utils::Printing > | |
Common line search printing utilities. | |
NOX::LineSearchCounters * | counter |
Common common counters for line searches. | |
NOX::LineSearch::Utils::Slope | slopeObj |
Common slope calculations for line searches. | |
int | numJvMults |
Counter for number of Jacobian-vector products. | |
int | numJ2vMults |
Counter for number of "double" Jacobian-vector products. | |
Teuchos::RCP< NOX::Observer > | observer |
Pointer to a user defined NOX::Observer object. | |
Nonlinear solver based on a rank-1 tensor method.
Solves using a rank-1 tensor method and a linesearch globalization.
At the kth nonlinear iteration, the solver does the following:
Computes the tensor direction by finding the root or smallest magnitude minimizer of the local model
where
and
Modifies the step according to a global strategy and updates as via a linesearch method, where is some function of . For instance, the curvilinear step is a function of the linesearch parameter and is a parametric step that spans the directions of the tensor step and the Newton step. At , the curvilinear step equals the full tensor step, and as nears 0, the curvilinear step approaches the Newton direction. This step provides a monotonic decrease in the norm of the local tensor model as varies from 0 to 1.
The solver iterates until the status tests (see NOX::StatusTest) determine either failure or convergence.
Input Parameters
To use this solver, set the "Nonlinear Solver" parameter to be "Tensor Based". Then, specify the following sublists with the appropriate parameters as indicated below.
"Direction" - Sublist of the direction parameters, passed to the NOX::Direction::Factory constructor. Defaults to an empty list.
"Method" - Name of the direction to be computed in this solver. "Tensor" and "Newton" are the only two valid choices. A sublist by this name specifies all of the parameters to be passed to the linear solver. See below under "Linear Solver".
"Rescue Bad Newton Solve" (Boolean) - If the linear solve does not meet the tolerance specified by the forcing term, then use the step anyway. Defaults to true.
"Linear Solver" - Sublist for the specific linear solver parameters that are passed to NOX::Abstract::Group::computeNewton() and NOX::Abstract::Group::applyJacobianInverse(). "Linear Solver" is itself a sublist of the list specified in "Method" above (i.e., "Tensor" or "Newton"). Below is a partial list of standard parameters usually available in common linear solvers. Check with the specific linear solver being used for other parameters.
"Max Iterations" - Maximum number of Arnoldi iterations (also max Krylov space dimension)
"Tolerance" - Relative tolerance for solving local model [default = 1e-4]
"Output Frequency" - Print output at every number of iterations [default = 20]
"Line Search" - Sublist of the line search parameters. Because the tensor step is not guaranteed to be a descent direction on the function, not all "basic" line search approaches would be appropriate. Thus, the LineSearch classes available to Newton's method (e.g., Polynomial, More-Thuente) are not used here. Instead, this solver class approriately handles technical considerations for tensor methods with its own set of global strategies. The following parameters specify the specific options for this line search:
"Method" - Name of the line search available to tensor methods Valid choices are:
<ul> <li> "Curvilinear" - Backtrack along the "curvilinear"
path that spans the tensor direction and the Newton direction and that maintains monotonicity on the tensor model. Recommended because it tends to be more robust and efficient than the other choices. [Default]
"Standard" - Backtrack along tensor direction unless it is not a descent direction, in which case backtrack along Newton direction.
"Dual" - Backtrack along both the Newton and tensor directions and choose the better of the two.
"Full Step" - Only use the full step and do not backtrack along both the Newton and tensor directions and choose the better of the two.
</ul>
"Lambda selection" - Flag for how to calculate the next linesearch parameter lambda. Valid choices are "Quadratic" and "Halving" (default). Quadratic constructs a quadratic interpolating polynomial from the last trial point and uses the minimum of this function as the next trial lambda (bounded by 0.1). Halving divides the linesearch parameter by 2 before each trial, which is simpler but tends to generate longer steps than quadratic.
"Default Step" - Starting value of the linesearch parameter (defaults to 1.0)
"Minimum Step" - Minimum acceptable linesearch parameter before the linesearch terminates (defaults to 1.0e-12). If there are many linesearch failures, then lowering this value is one thing to try.
"Recovery Step Type" - Determines the step size to take when the line search fails. Choices are:
"Constant" [default] - Uses a constant value set in "Recovery Step".
"Last Computed Step" - Uses the last value computed by the line search algorithm.
"Recovery Step" - Step parameter to take when the line search fails (defaults to value for "Default Step")
"Max Iters" - Maximum number of iterations (i.e., backtracks)
"Solver Options" - Sublist of general solver options.
"User Defined Pre/Post Operator" is supported. See NOX::Parameter::PrePostOperator for more details.
Output Parameters
Every time solve() is called, a sublist for output parameters called "Output" will be created and will contain the following parameters:
"Nonlinear Iterations" - Number of nonlinear iterations
"2-Norm of Residual" - L-2 norm of the final residual .
References
B. W. Bader, Tensor-Krylov methods for solving large-scale systems of nonlinear equations, Ph.D. Thesis, 2003, University of Colorado, Boulder, Colorado.
B. W. Bader, Tensor-Krylov methods for solving large-scale systems of nonlinear equations, submitted to SIAM J. Numer. Anal.
B. W. Bader and R. B. Schnabel, Curvilinear linesearch for tensor methods, SISC, 25(2):604-622.
R. B. Schnabel and P. D. Frank, Tensor methods for nonlinear equations, SIAM J. Numer. Anal., 21(5):815-843.
|
protected |
|
protected |
Flag for the direction to be computed this iteration.
Enumerated list for each type of line search
|
protected |
NOX::Solver::TensorBased::TensorBased | ( | const Teuchos::RCP< NOX::Abstract::Group > & | grp, |
const Teuchos::RCP< NOX::StatusTest::Generic > & | tests, | ||
const Teuchos::RCP< Teuchos::ParameterList > & | params | ||
) |
|
virtual |
Resets the solver for another solve. This resets the counters and status only. Uses the final solution from the last solve as the initial guess for the next solve.
NOTE: All NOX solvers will call reset() automatically at teh beginning of the solve() method. We add the reset() method to the solver interface for the application to call in case the application needs to reset counters and status manually before the next call to solve() is made.
Implements NOX::Solver::Generic.
References NOX::StatusTest::Unconverged.
Referenced by TensorBased().
|
virtual |
Solve the nonlinear problem and return final status.
By "solve", we call iterate() until the NOX::StatusTest value is either NOX::StatusTest::Converged or NOX::StatusTest::Failed.
Implements NOX::Solver::Generic.
References Teuchos::ParameterList::set(), Teuchos::ParameterList::sublist(), and NOX::StatusTest::Unconverged.
|
protected |
Current tensor term vector (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
|
protected |
|
protected |
Current Newton direction (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
|
protected |
Previous solution pointer.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by getPreviousSolutionGroupPtr().
|
protected |
Residual vector (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
|
protected |
Vector to previous point (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
|
protected |
Current tensor direction (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
|
protected |
Working vector (pointer).
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.