NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
NOX::Solver::TensorBased Class Reference

Nonlinear solver based on a rank-1 tensor method. More...

#include <NOX_Solver_TensorBased.H>

Inheritance diagram for NOX::Solver::TensorBased:
Inheritance graph
[legend]
Collaboration diagram for NOX::Solver::TensorBased:
Collaboration graph
[legend]

Public Member Functions

 TensorBased (const Teuchos::RCP< NOX::Abstract::Group > &grp, const Teuchos::RCP< NOX::StatusTest::Generic > &tests, const Teuchos::RCP< Teuchos::ParameterList > &params)
 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 > &params)
 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::GlobalDataglobalDataPtr
 Pointer to the global data object.
 
Teuchos::RCP< NOX::UtilsutilsPtr
 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::ParameterListlinearParamsPtr
 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
print
 Common line search printing utilities.
 
NOX::LineSearchCounterscounter
 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::Observerobserver
 Pointer to a user defined NOX::Observer object.
 

Detailed Description

Nonlinear solver based on a rank-1 tensor method.

Solves $ F(x)=0$ using a rank-1 tensor method and a linesearch globalization.

At the kth nonlinear iteration, the solver does the following:

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.

Output Parameters

Every time solve() is called, a sublist for output parameters called "Output" will be created and will contain the following parameters:

References

Author
Brett Bader (SNL 9233)

Member Enumeration Documentation

Algorithms used to determine convergence of the line search.

Enumerator
ArmijoGoldstein 

Sufficient decrease condition.

AredPred 

Ared/Pred condition.

None 

Just accept the first step.

Flag for the direction to be computed this iteration.

Enumerated list for each type of line search

Type of recovery step to use.

Enumerator
Constant 

Use a constant value.

LastComputedStep 

Use the last value computed in the line search algorithm.

Constructor & Destructor Documentation

NOX::Solver::TensorBased::TensorBased ( const Teuchos::RCP< NOX::Abstract::Group > &  grp,
const Teuchos::RCP< NOX::StatusTest::Generic > &  tests,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)

Constructor.

See reset() for description.

References NOX::DeepCopy, reset(), and NOX::ShapeCopy.

Member Function Documentation

void NOX::Solver::TensorBased::reset ( )
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().

NOX::StatusTest::StatusType NOX::Solver::TensorBased::solve ( )
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.

Member Data Documentation

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::aVecPtr
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.

Teuchos::ParameterList* NOX::Solver::TensorBased::linearParamsPtr
protected

Line Search parameters.

Direction parameters. Parameters for the Linear Solver of the local model.

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::newtonVecPtr
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.

Teuchos::RCP<NOX::Abstract::Group> NOX::Solver::TensorBased::oldSolnPtr
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().

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::residualVecPtr
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.

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::sVecPtr
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.

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::tensorVecPtr
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.

Teuchos::RCP<NOX::Abstract::Vector> NOX::Solver::TensorBased::tmpVecPtr
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.


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