|
NOX Development
|
Concrete implementation for creating an Epetra_RowMatrix Jacobian via finite differencing of the residual. More...
#include <NOX_Epetra_FiniteDifference.H>


Public Types | |
| enum | DifferenceType { Forward , Backward , Centered } |
| Define types for use of the perturbation parameter | |
Public Member Functions | |
| FiniteDifference (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, double beta=1.0e-6, double alpha=1.0e-4) | |
| Constructor with scalar beta. | |
| FiniteDifference (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< const Epetra_Vector > &beta, double alpha=1.0e-4) | |
| Constructor with vector beta. | |
| FiniteDifference (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< Epetra_CrsGraph > &g, double beta=1.0e-6, double alpha=1.0e-4) | |
| Constructor that takes a pre-constructed Epetra_CrsGraph so it does not have to determine the non-zero entries in the matrix. | |
| FiniteDifference (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< Epetra_CrsGraph > &g, const Teuchos::RCP< const Epetra_Vector > &beta, double alpha=1.0e-4) | |
| Constructor with output control that takes a pre-constructed Epetra_CrsGraph so it does not have to determine the non-zero entries in the matrix. | |
| virtual | ~FiniteDifference () |
| Pure virtual destructor. | |
| virtual const char * | Label () const |
| Returns a character std::string describing the name of the operator. | |
| virtual int | SetUseTranspose (bool UseTranspose) |
| If set true, the transpose of this operator will be applied. | |
| virtual int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Return the result on an Epetra_Operator applied to an Epetra_MultiVector X in Y. | |
| virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Return the result on an Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. | |
| virtual bool | UseTranspose () const |
| Returns the current use transpose setting. | |
| virtual bool | HasNormInf () const |
| Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
| virtual const Epetra_Map & | OperatorDomainMap () const |
| Returns the Epetra_BlockMap object associated with the domain of this matrix operator. | |
| virtual const Epetra_Map & | OperatorRangeMap () const |
| Returns the Epetra_BlockMap object associated with the range of this matrix operator. | |
| virtual bool | Filled () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | NumMyRowEntries (int MyRow, int &NumEntries) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | MaxNumEntries () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | InvRowSums (Epetra_Vector &x) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | LeftScale (const Epetra_Vector &x) |
| See Epetra_RowMatrix documentation. | |
| virtual int | InvColSums (Epetra_Vector &x) const |
| See Epetra_RowMatrix documentation. | |
| virtual int | RightScale (const Epetra_Vector &x) |
| See Epetra_RowMatrix documentation. | |
| virtual double | NormInf () const |
| See Epetra_RowMatrix documentation. | |
| virtual double | NormOne () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | NumGlobalNonzeros () const |
| See Epetra_RowMatrix documentation. | |
| virtual long long | NumGlobalNonzeros64 () const |
| virtual int | NumGlobalRows () const |
| See Epetra_RowMatrix documentation. | |
| virtual long long | NumGlobalRows64 () const |
| virtual int | NumGlobalCols () const |
| See Epetra_RowMatrix documentation. | |
| virtual long long | NumGlobalCols64 () const |
| virtual int | NumGlobalDiagonals () const |
| See Epetra_RowMatrix documentation. | |
| virtual long long | NumGlobalDiagonals64 () const |
| virtual int | NumMyNonzeros () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | NumMyRows () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | NumMyCols () const |
| See Epetra_RowMatrix documentation. | |
| virtual int | NumMyDiagonals () const |
| See Epetra_RowMatrix documentation. | |
| virtual bool | LowerTriangular () const |
| See Epetra_RowMatrix documentation. | |
| virtual bool | UpperTriangular () const |
| See Epetra_RowMatrix documentation. | |
| virtual const Epetra_Comm & | Comm () const |
| See Epetra_RowMatrix documentation. | |
| virtual const Epetra_Map & | RowMatrixRowMap () const |
| See Epetra_RowMatrix documentation. | |
| virtual const Epetra_Map & | RowMatrixColMap () const |
| See Epetra_RowMatrix documentation. | |
| virtual const Epetra_Import * | RowMatrixImporter () const |
| See Epetra_RowMatrix documentation. | |
| virtual const Epetra_BlockMap & | Map () const |
| See Epetra_SrcDistObj documentation. | |
| 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 bool | computeJacobian (const Epetra_Vector &x) |
| Compute Jacobian given the specified input vector, x. Returns true if computation was successful. | |
| virtual bool | computePreconditioner (const Epetra_Vector &x, Epetra_Operator &Prec, Teuchos::ParameterList *precParams=0) |
| Compute an Epetra_RowMatrix to be used by Aztec preconditioners 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) | |
| virtual Epetra_CrsMatrix & | getUnderlyingMatrix () const |
| An accessor method for the underlying Epetra_CrsMatrix. | |
| virtual void | Print (std::ostream &) const |
| Output the underlying matrix. | |
| void | setGroupForComputeF (NOX::Abstract::Group &group) |
| Register 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. | |
Public Member Functions inherited from NOX::Epetra::Interface::Jacobian | |
| Jacobian () | |
| Constructor. | |
| virtual | ~Jacobian () |
| Destructor. | |
Public Member Functions inherited from NOX::Epetra::Interface::Preconditioner | |
| Preconditioner () | |
| Constructor. | |
| virtual | ~Preconditioner () |
| Destructor. | |
Protected Types | |
| enum | BetaType { Scalar , Vector } |
Protected Member Functions | |
| Teuchos::RCP< Epetra_CrsMatrix > | createGraphAndJacobian (Interface::Required &i, const Epetra_Vector &x) |
| Constructs an Epetra_CrsGraph and Epetra_RowMatrix for the Jacobian. This is only called if the user does not supply an Epetra_CrsGraph. | |
| bool | computeF (const Epetra_Vector &input, Epetra_Vector &result, NOX::Epetra::Interface::Required::FillType) |
Protected Attributes | |
| const NOX::Utils | utils |
| Printing Utilities object. | |
| Teuchos::RCP< Epetra_CrsGraph > | graph |
| Pointer to the Jacobian graph. | |
| Teuchos::RCP< Epetra_CrsMatrix > | jacobian |
| Pointer to the Jacobian. | |
| Teuchos::RCP< NOX::Epetra::Interface::Required > | interface |
| User provided interface function. | |
| Epetra_Vector | x_perturb |
| Perturbed solution vector - a work array that needs to be mutable. | |
| Epetra_Vector | fo |
| Function evaluation at currentX - a work array that needs to be mutable. | |
| Epetra_Vector | fp |
| Function evaluation at perturbX - a work array that needs to be mutable. | |
| Teuchos::RCP< Epetra_Vector > | fmPtr |
| Optional pointer to function evaluation at -perturbX - needed only for centered finite differencing. | |
| Epetra_Vector | Jc |
| Column vector of the jacobian - a work array that needs to be mutable. | |
| double | alpha |
| Constant for the perturbation calculation. | |
| double | beta |
| Constant for the perturbation calculation. | |
| Teuchos::RCP< const Epetra_Vector > | betaVector |
| Vector for the perturbation calculation. | |
| BetaType | betaType |
| Flag that sets whether | |
| DifferenceType | diffType |
| Define types for use of the perturbation parameter | |
| std::string | label |
| label for the Epetra_RowMatrix | |
| bool | useGroupForComputeF |
| Flag to enables the use of a group instead of the interface for the computeF() calls in the directional difference calculation. | |
| Teuchos::RCP< NOX::Abstract::Group > | groupPtr |
| Pointer to the group for possible use in computeF() calls. | |
Concrete implementation for creating an Epetra_RowMatrix Jacobian via finite differencing of the residual.
The Jacobian entries are calculated via 1st order finite differencing. This requires 

![\[ J_{ij} = \frac{\partial F_i}{\partial x_j} = \frac{F_i(x+\delta\mathbf{e}_j) - F_i(x)}{\delta} \]](form_668.png)
where 




The perturbation, 
![\[ \delta = \alpha * | x_j | + \beta \]](form_672.png)
![\[ \delta = \alpha * | x_j | + \beta_j \]](form_673.png)
where 

Since this inherits from the Epetra_RowMatrix class, it can be used as the preconditioning matrix for AztecOO preconditioners. This method is very inefficient when computing the Jacobian and is not recommended for large-scale systems but only for debugging purposes.
|
protected |
Define types for the 

|
virtual |
Compute Jacobian given the specified input vector, x. Returns true if computation was successful.
Reimplemented in NOX::Epetra::FiniteDifferenceColoring, and NOX::Epetra::FiniteDifferenceColoringWithUpdate.
References computeJacobian().
|
virtual |
Compute Jacobian given the specified input vector, x. Returns true if computation was successful.
Implements NOX::Epetra::Interface::Jacobian.
Reimplemented in NOX::Epetra::FiniteDifferenceColoring, and NOX::Epetra::FiniteDifferenceColoringWithUpdate.
References alpha, beta, betaType, diffType, NOX::Epetra::Interface::Required::FD_Res, fmPtr, fo, fp, jacobian, Jc, NOX::Utils::out(), utils, and x_perturb.
Referenced by computeJacobian(), and computePreconditioner().
|
virtual |
Compute an Epetra_RowMatrix to be used by Aztec preconditioners given the specified input vector, x. Returns true if computation was successful.
Implements NOX::Epetra::Interface::Preconditioner.
References computeJacobian().