NOX
Development
|
Concrete implementation of NOX::Abstract::Group for Trilinos/Epetra. More...
#include <NOX_Epetra_Group.H>
Public Member Functions | |
Group (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess) | |
Constructor with NO linear system (VERY LIMITED). More... | |
Group (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< NOX::Epetra::LinearSystem > &linSys) | |
Standard Constructor. | |
Group (const NOX::Epetra::Group &source, NOX::CopyType type=NOX::DeepCopy) | |
Copy constructor. If type is DeepCopy, takes ownership of valid shared linear system. | |
virtual | ~Group () |
Destructor. | |
virtual NOX::Abstract::Group & | operator= (const NOX::Abstract::Group &source) |
Copies the source group into this group. More... | |
virtual NOX::Abstract::Group & | operator= (const NOX::Epetra::Group &source) |
See operator=(const NOX::Abstract::Group&);. | |
virtual Teuchos::RCP < NOX::Abstract::Group > | clone (CopyType type=DeepCopy) const |
Create a new Group of the same derived type as this one by cloning this one, and return a ref count pointer to the new group. More... | |
virtual Teuchos::RCP < NOX::Epetra::Interface::Required > | getRequiredInterface () |
Return the userInterface. | |
virtual Teuchos::RCP< const NOX::Epetra::LinearSystem > | getLinearSystem () const |
Return the Linear System. | |
virtual Teuchos::RCP < NOX::Epetra::LinearSystem > | getLinearSystem () |
Return the Linear System. | |
virtual NOX::Abstract::Group::ReturnType | computeJacobianConditionNumber (int maxIters, double tolerance, int krylovSubspaceSize=100, bool printOutput=false) |
virtual double | getJacobianConditionNumber () const |
Returns the condition number of the Jacobian matrix. | |
virtual void | disableLinearResidualComputation (const bool disableChoice) |
Sets option to disable linear resid computation. If disabled, this saves on a MatVec per Newton but disallows inexact Newton methods. | |
"Compute" functions. | |
virtual void | setX (const NOX::Epetra::Vector &y) |
virtual void | setX (const NOX::Abstract::Vector &y) |
Set the solution vector x to y. More... | |
virtual void | computeX (const Group &grp, const NOX::Epetra::Vector &d, double step) |
virtual void | computeX (const NOX::Abstract::Group &grp, const NOX::Abstract::Vector &d, double step) |
Compute x = grp.x + step * d. More... | |
virtual NOX::Abstract::Group::ReturnType | computeF () |
Compute and store F(x). More... | |
virtual NOX::Abstract::Group::ReturnType | computeJacobian () |
Compute and store Jacobian. More... | |
virtual NOX::Abstract::Group::ReturnType | computeGradient () |
Compute and store gradient. More... | |
virtual NOX::Abstract::Group::ReturnType | computeNewton (Teuchos::ParameterList ¶ms) |
Compute the Newton direction, using parameters for the linear solve. More... | |
"Is" functions | |
Checks to see if various objects have been computed. Returns true if the corresponding "compute" function has been called since the last update to the solution vector (via instantiation or computeX). | |
virtual bool | isF () const |
Return true if F is valid. | |
virtual bool | isJacobian () const |
Return true if the Jacobian is valid. More... | |
virtual bool | isGradient () const |
Return true if the gradient is valid. More... | |
virtual bool | isNewton () const |
Return true if the Newton direction is valid. More... | |
virtual bool | isNormNewtonSolveResidual () const |
Returns true if the value of the Norm of the linear model for a full Newton step ||Js + f|| is valid with respect to the current solution vector. | |
virtual bool | isPreconditioner () const |
Returns true if an explicitly constructed preconditioner exists (i.e. one that is computed and saved for further use in multiple calls to applyRightPreconditioner). | |
virtual bool | isConditionNumber () const |
Returns true if the condition number has been computed. | |
"Get" functions | |
Note that these function do not check whether or not the vectors are valid. Must use the "Is" functions for that purpose. | |
virtual const NOX::Abstract::Vector & | getX () const |
Return solution vector. | |
virtual const NOX::Abstract::Vector & | getF () const |
Return F(x) | |
virtual double | getNormF () const |
Return 2-norm of F(x). More... | |
virtual const NOX::Abstract::Vector & | getGradient () const |
Return gradient. | |
virtual const NOX::Abstract::Vector & | getNewton () const |
Return Newton direction. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getXPtr () const |
Return RCP to solution vector. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getFPtr () const |
Return RCP to F(x) | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getGradientPtr () const |
Return RCP to gradient. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getNewtonPtr () const |
Return RCP to Newton direction. | |
virtual NOX::Abstract::Group::ReturnType | getNormLastLinearSolveResidual (double &residual) const |
Returns the 2-norm of the residual of the linear model used in the Newton solve computation, ||Js+f||. This does not account for line search adjustments to the step length! | |
Public Member Functions inherited from NOX::Abstract::Group | |
Group () | |
Constructor. More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianInverse (Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result. More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
applyJacobian for multiple right-hand sides More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTransposeMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
applyJacobianTranspose for multiple right-hand sides More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianInverseMultiVector (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
applyJacobianInverse for multiple right-hand sides More... | |
virtual NOX::Abstract::Group::ReturnType | applyRightPreconditioningMultiVector (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
applyRightPreconditioning for multiple right-hand sides More... | |
virtual const NOX::Abstract::Vector & | getScaledX () const |
Protected Member Functions | |
virtual void | resetIsValid () |
resets the isValid flags to false | |
virtual void | logLastLinearSolveStats (NOX::SolverStats &stats) const |
Adds statistics from last linear solve to the SovlerStats object. | |
virtual bool | computeNormNewtonSolveResidual () |
Computes the 2-norm of the residual of the linear model used in the Newton solve computation, ||Js+f||. | |
Protected Attributes | |
const NOX::Utils | utils |
Printing Utilities object. | |
double | normNewtonSolveResidual |
2-Norm of the Newton solve residual: ||Js+f|| | |
double | conditionNumber |
condition number of Jacobian | |
Teuchos::RCP < AztecOOConditionNumber > | azConditionNumberPtr |
Pointer to the condition number object. | |
bool | linearResidCompDisabled |
Teuchos::RCP < NOX::Epetra::Interface::Required > | userInterfacePtr |
Reference to the user supplied interface functions. | |
bool | linearSolveConverged |
int | numIterations |
double | achievedTol |
Vectors | |
Teuchos::RCP< NOX::Epetra::Vector > | xVectorPtr |
Solution vector pointer. | |
NOX::Epetra::Vector & | xVector |
Solution vector. | |
Teuchos::RCP< NOX::Epetra::Vector > | RHSVectorPtr |
Right-hand-side vector (function evaluation). | |
NOX::Epetra::Vector & | RHSVector |
Right-hand-side vector pointer (function evaluation). | |
Teuchos::RCP< NOX::Epetra::Vector > | gradVectorPtr |
Gradient vector pointer(steepest descent vector). | |
NOX::Epetra::Vector & | gradVector |
Gradient vector (steepest descent vector). | |
Teuchos::RCP< NOX::Epetra::Vector > | NewtonVectorPtr |
Newton direction vector pointer. | |
NOX::Epetra::Vector & | NewtonVector |
Newton direction vector. | |
Teuchos::RCP< Epetra_Vector > | tmpVectorPtr |
An extra temporary vector, only allocated if needed. | |
IsValid flags | |
True if the current solution is up-to-date with respect to the currect xVector. | |
bool | isValidRHS |
bool | isValidJacobian |
bool | isValidGrad |
bool | isValidNewton |
bool | isValidNormNewtonSolveResidual |
bool | isValidPreconditioner |
bool | isValidSolverJacOp |
bool | isValidConditionNumber |
Shared Operators | |
Teuchos::RCP < NOX::SharedObject < NOX::Epetra::LinearSystem, NOX::Epetra::Group > > | sharedLinearSystemPtr |
Pointer to shared Jacobian matrix. | |
NOX::SharedObject < NOX::Epetra::LinearSystem, NOX::Epetra::Group > & | sharedLinearSystem |
Reference to shared Jacobian matrix. | |
Jacobian operations. | |
Operations using the Jacobian matrix. These may not be defined in matrix-free scenarios. | |
virtual NOX::Abstract::Group::ReturnType | const |
Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result. More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobian (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const |
virtual NOX::Abstract::Group::ReturnType | applyJacobian (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Applies Jacobian to the given input vector and puts the answer in the result. More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTranspose (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Applies Jacobian-Transpose to the given input vector and puts the answer in the result. More... | |
virtual NOX::Abstract::Group::ReturnType | applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const |
virtual NOX::Abstract::Group::ReturnType | applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Apply right preconditiong to the given input vector. More... | |
Additional Inherited Members | |
Public Types inherited from NOX::Abstract::Group | |
enum | ReturnType { Ok, NotDefined, BadDependency, NotConverged, Failed } |
The computation of, say, the Newton direction in computeNewton() may fail in many different ways, so we have included a variety of return codes to describe the failures. Of course, we also have a code for success. More... | |
Concrete implementation of NOX::Abstract::Group for Trilinos/Epetra.
This group is set up to use the linear algebra services provided
through the Trilinos/Epetra package with AztecOO for the linear solver.
Group::Group | ( | Teuchos::ParameterList & | printingParams, |
const Teuchos::RCP< NOX::Epetra::Interface::Required > & | i, | ||
const NOX::Epetra::Vector & | initialGuess | ||
) |
Constructor with NO linear system (VERY LIMITED).
WARNING: If this constructor is used, then methods that require a Jacobian or preconditioning will not be available. You will be limited to simple algorithms like nonlinear-CG with no preconditioning.
References resetIsValid().
|
virtual |
Applies Jacobian to the given input vector and puts the answer in the result.
Computes
where is the Jacobian, is the input vector, and is the result vector.
Reimplemented from NOX::Abstract::Group.
|
virtual |
Applies Jacobian-Transpose to the given input vector and puts the answer in the result.
Computes
where is the Jacobian, is the input vector, and is the result vector.
Reimplemented from NOX::Abstract::Group.
|
virtual |
Apply right preconditiong to the given input vector.
Let be a right preconditioner for the Jacobian ; in other words, is a matrix such that
Compute
where is the input vector and is the result vector.
If useTranspose is true, then the transpose of the preconditioner is applied:
The transpose preconditioner is currently only required for Tensor methods.
The "Tolerance" parameter specifies that the solution should be such that
params
The parameters are from the "Linear %Solver" sublist of the "Direction" sublist that is passed to solver during construction.
Reimplemented from NOX::Abstract::Group.
|
virtual |
Create a new Group of the same derived type as this one by cloning this one, and return a ref count pointer to the new group.
If type is NOX::DeepCopy, then we need to create an exact replica of "this". Otherwise, if type is NOX::ShapeCopy, we need only replicate the shape of "this" (only the memory is allocated, the values are not copied into the vectors and Jacobian). Returns NULL if clone is not supported.
Implements NOX::Abstract::Group.
Reimplemented in LOCA::Epetra::Group.
References Teuchos::rcp().
|
virtual |
Compute and store F(x).
Implements NOX::Abstract::Group.
Reimplemented in LOCA::Epetra::Group.
References NOX::Epetra::Interface::Required::computeF(), NOX::Epetra::Vector::getEpetraVector(), isF(), NOX::Abstract::Group::Ok, NOX::Epetra::Interface::Required::Residual, RHSVector, userInterfacePtr, and xVector.
Referenced by LOCA::Epetra::Group::computeF().
|
virtual |
Compute and store gradient.
We can pose the nonlinear equation problem as an optimization problem as follows:
In that case, the gradient (of ) is defined as
Reimplemented from NOX::Abstract::Group.
References gradVector, isF(), isGradient(), isJacobian(), NOX::Abstract::Group::Ok, RHSVector, and sharedLinearSystem.
|
virtual |
Compute and store Jacobian.
Recall that
The Jacobian is denoted by and defined by
Reimplemented from NOX::Abstract::Group.
Reimplemented in LOCA::Epetra::Group.
References isJacobian(), NOX::Abstract::Group::Ok, sharedLinearSystem, and xVector.
Referenced by LOCA::Epetra::Group::computeJacobian().
|
virtual |
Compute the Newton direction, using parameters for the linear solve.
The Newton direction is the solution, s, of
The parameters are from the "Linear %Solver" sublist of the "Direction" sublist that is passed to solver during construction.
The "Tolerance" parameter may be added/modified in the sublist of "Linear Solver" parameters that is passed into this function. The solution should be such that
params
Reimplemented from NOX::Abstract::Group.
References NOX::Abstract::Group::applyJacobianInverse(), computeNormNewtonSolveResidual(), NOX::Epetra::Vector::init(), isF(), isJacobian(), isNewton(), NewtonVector, NOX::Abstract::Group::Ok, RHSVector, and NOX::Epetra::Vector::scale().
|
virtual |
Compute x = grp.x + step * d.
Let denote this group's solution vector. Let denote the result of grp.getX(). Then set
Implements NOX::Abstract::Group.
|
virtual |
Return 2-norm of F(x).
In other words,
Implements NOX::Abstract::Group.
References isF(), NOX::Epetra::Vector::norm(), and RHSVectorPtr.
|
virtual |
Return true if the gradient is valid.
Reimplemented from NOX::Abstract::Group.
Referenced by computeGradient(), and getGradient().
|
virtual |
Return true if the Jacobian is valid.
Reimplemented from NOX::Abstract::Group.
References sharedLinearSystem.
Referenced by computeGradient(), computeJacobian(), and computeNewton().
|
virtual |
Return true if the Newton direction is valid.
Reimplemented from NOX::Abstract::Group.
Referenced by computeNewton(), and getNewton().
|
virtual |
Copies the source group into this group.
Implements NOX::Abstract::Group.
Reimplemented in LOCA::Epetra::Group.
Referenced by LOCA::Epetra::Group::operator=().
|
virtual |
Set the solution vector x to y.
Implements NOX::Abstract::Group.
NOX::Abstract::Group::ReturnType NOX::Epetra::Group::const |
Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result.
Computes
where is the Jacobian, is the input vector, and is the result vector.
The "Tolerance" parameter specifies that the solution should be such that
\return <ul> <li> NOX::Abstract::Group::NotDefined - Returned by default implementation
in NOX::Abstract::Group NOX::Abstract::Group::BadDependency - If has not been computed NOX::Abstract::Group::NotConverged - If the linear solve fails to satisfy the "Tolerance" specified in params
NOX::Abstract::Group::Failed - If the computation fails NOX::Abstract::Group::Ok - Otherwise
The parameter "Tolerance" may be added/modified in the list of parameters - this is the ideal solution tolerance for an iterative linear solve.
The parameter "Reuse Preconditioner" is a boolean that tells the group to turn off control of preconditioner recalculation. This is a dangerous flag but can really speed the computations if the user knows what they are doing. Toggling this flag is left to the user (ideally it should be done through a status test). Defaults to false.