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

Concrete implementation for creating an Epetra_Operator Jacobian based on the Matrix-Free Newton-Krylov method. More...

#include <NOX_Epetra_MatrixFree.H>

Inheritance diagram for NOX::Epetra::MatrixFree:
Inheritance graph
[legend]
Collaboration diagram for NOX::Epetra::MatrixFree:
Collaboration graph
[legend]

Public Types

enum  DifferenceType { Forward, Backward, Centered }
 Define types for use of the perturbation parameter $ \delta$.
 

Public Member Functions

 MatrixFree (Teuchos::ParameterList &printParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &cloneVector, bool useNewPerturbation=false)
 Constructor. More...
 
virtual ~MatrixFree ()
 Pure virtual destructor.
 
virtual int SetUseTranspose (bool UseTranspose)
 If set true, transpose of this operator will be applied. More...
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix.
 
virtual const char * Label () const
 Returns a character std::string describing the operator.
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting.
 
virtual bool HasNormInf () const
 Returns true if the this object can provide an approximate Inf-norm, false otherwise.
 
virtual const Epetra_CommComm () const
 Returns a reference to the Epetra_Comm communicator associated with this operator.
 
virtual const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_BlockMap object associated with the domain of this matrix operator.
 
virtual const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_BlockMap object associated with the range of this matrix operator.
 
virtual bool computeJacobian (const Epetra_Vector &x, Epetra_Operator &Jac)
 Compute Jacobian given the specified input vector, x. Returns true if computation was successful.
 
virtual void setDifferenceMethod (DifferenceType type)
 Set the type of perturbation method used (default is Forward)
 
void setLambda (double lambda_)
 Allows the user to change the value of $ \lambda $ in the perturbation calculation.
 
void setComputePerturbation (bool bVal)
 
void setPerturbation (double eta_)
 Set the perturbation parameter $ \eta $.
 
double getPerturbation () const
 Returns the most recently used value of the perturbation parameter $ \eta $.
 
void setGroupForComputeF (const NOX::Abstract::Group &group)
 Clone a NOX::Abstract::Group derived object and use the computeF() method of that group for the perturbation instead of the NOX::Epetra::Interface::Required::computeF() method. This is required for LOCA to get the operators correct during homotopy.
 
void setSolverForComputeJacobian (const Teuchos::RCP< NOX::Solver::Generic > &slvr)
 Save a RCP to a solver, and use the Solver's current Group's computeF() in the computeJacobian call, which can save a function call by respecting the isValid flag.
 
- Public Member Functions inherited from NOX::Epetra::Interface::Jacobian
 Jacobian ()
 Constructor.
 
virtual ~Jacobian ()
 Destructor.
 

Protected Attributes

std::string label
 Label for matrix.
 
Teuchos::RCP
< NOX::Epetra::Interface::Required
interface
 User provided interface function.
 
NOX::Epetra::Vector currentX
 The current solution vector.
 
NOX::Epetra::Vector perturbX
 Perturbed solution vector.
 
NOX::Epetra::Vector fo
 Function evaluation at currentX.
 
NOX::Epetra::Vector fp
 Function evaluation at perturbX.
 
Teuchos::RCP< NOX::Epetra::VectorfmPtr
 Optional pointer to function evaluation at -perturbX - needed only for centered finite differencing.
 
Teuchos::RCP< const Epetra_MapepetraMap
 
DifferenceType diffType
 Define types for use of the perturbation parameter $ \delta$.
 
double lambda
 Scale factor for eta calculation.
 
double eta
 Perturbation value to use in the directional derivative.
 
double userEta
 User specified perturbation value to use in the directional derivative. Set by setPerturbation().
 
bool computeEta
 Flag that determines if we should calculate eta or use a value set by the user.
 
bool useGroupForComputeF
 Flag to enables the use of a group instead of the interface for the computeF() calls in the directional difference calculation.
 
bool useSolverForComputeJacobian
 Flag to enables the use of Solver's Group instead of the interface for the computeF() calls in the directional difference calculation, to make use of isValid flags.
 
bool useNewPerturbation
 A new perturbation formulation developed by C. T. Kelley and A. G. Salinger can be used by the constructor flag useNewPerturbation = true.
 
Teuchos::RCP
< NOX::Abstract::Group
groupPtr
 Pointer to the group for possible use in computeF() calls.
 
Teuchos::RCP
< NOX::Solver::Generic
slvrPtr
 Pointer to the Solver for possible use in computeF() calls.
 
NOX::Utils utils
 Printing utilities.
 

Detailed Description

Concrete implementation for creating an Epetra_Operator Jacobian based on the Matrix-Free Newton-Krylov method.

Matrix-Free Newton-Krylov is a method that takes advantage of the fact the Newton Krylov solvers do not require an explicit Jacobian matrix. Newton-Krylov solvers only require the matrix-vector product $ Jy$ in the iteration sequence. This product can approximated by the following:

\[ Jy = \frac{F(x + \delta y) - F(x)}{\delta} \]

where $ J$ is the Jacobian, $ F$ is the function evaluation, $ x$ is the solution vector, $ y$ is the vector to be operated on, and $\delta$ is a scalar perturbation calculated by:

\[ \delta = \lambda * (\lambda + \frac{\| x\|}{\| y\|} ) \]

where $ \lambda = 1.0e-6 $.

Constructor & Destructor Documentation

MatrixFree::MatrixFree ( Teuchos::ParameterList printParams,
const Teuchos::RCP< NOX::Epetra::Interface::Required > &  i,
const NOX::Epetra::Vector cloneVector,
bool  useNewPerturbation = false 
)

Constructor.

The vector x is used to clone the solution vector.

References currentX, epetraMap, fo, fp, NOX::Epetra::Vector::getEpetraVector(), NOX::Epetra::Vector::init(), perturbX, and Teuchos::rcp().

Member Function Documentation

int MatrixFree::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual
int MatrixFree::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

Parameters
X- A Epetra_MultiVector of dimension NumVectors to solve for.
Y-A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.
Warning
In order to work with AztecOO, any implementation of this method must support the case where X and Y are the same object.

Implements Epetra_Operator.

References NOX::Utils::out(), TEUCHOS_UNREACHABLE_RETURN, and utils.

void MatrixFree::setComputePerturbation ( bool  bVal)

Flag that toggles whether MatrixFree should compute the perturbation parameter $ \eta $ or use a value supplied by the user through setPerturbation().

References computeEta.

int MatrixFree::SetUseTranspose ( bool  UseTranspose)
virtual

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used implicitly. Setting this flag affects only the Apply() and ApplyInverse() methods. If the implementation of this interface does not support transpose use, this method should return a value of -1.

Parameters
UseTranspose-If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.

Implements Epetra_Operator.

References NOX::Utils::out(), and utils.

Member Data Documentation

Teuchos::RCP<const Epetra_Map> NOX::Epetra::MatrixFree::epetraMap
protected

Epetra_Map object used in the returns of the Epetra_Operator derived methods.

If the user is using Epetra_BlockMaps, then NOX::Epetra::MatrixFree must create an equivalent Epetra_Map from the Epetra_BlockMap that can be used as the return object of the OperatorDomainMap() and OperatorRangeMap() methods.

Referenced by MatrixFree(), OperatorDomainMap(), and OperatorRangeMap().


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