

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 $ N + 1 $ calls to computeF() where $ N $ is the number
of unknowns in the problem.
\\[ J_{ij} = \\frac{\\partial F_i}{\\partial x_j} =
\\frac{F_i(x+\\delta\\mathbf{e}_j) - F_i(x)}{\\delta} \\]
where $J$ is the Jacobian, $F$ is the function evaluation, $x$ is the
solution vector, and $\\delta$ is a small perturbation to the $x_j$
entry.
The perturbation, $ \\delta $, is calculated based on one of the
following equations:
\\[ \\delta = \\alpha * | x_j | + \\beta \\] \\[ \\delta
= \\alpha * | x_j | + \\beta_j \\]
where $ \\alpha $ is a scalar value (defaults to 1.0e-4) and $
\\beta $ can be either a scalar or a vector (defaults to a scalar
value of 1.0e-6). The choice is defined by the type of constructor
used. All parameters are supplied in the constructor. In addition to
the forward difference derivative approximation, backward or centered
differences can be used via the setDifferenceMethod function. Note
that centered difference provides second order spatial accuracy but at
the cost of twice as many function evaluations.
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.
C++ includes: NOX_Epetra_FiniteDifference.H
| def PyTrilinos.NOX.Epetra.FiniteDifference.__init__ | ( | self, | |
| args | |||
| ) |
__init__(NOX::Epetra::FiniteDifference self, ParameterList printingParams, Teuchos::RCP< NOX::Epetra::Interface::Required > const & i,
Vector initialGuess, double beta=1.0e-6, double alpha=1.0e-4) -> FiniteDifference
__init__(NOX::Epetra::FiniteDifference self, ParameterList printingParams, Teuchos::RCP< NOX::Epetra::Interface::Required > const & i,
Vector initialGuess, Teuchos::RCP< Epetra_Vector const > const & beta,
double alpha=1.0e-4) -> FiniteDifference
__init__(NOX::Epetra::FiniteDifference self, ParameterList printingParams, Teuchos::RCP< NOX::Epetra::Interface::Required > const & i,
Vector initialGuess, Teuchos::RCP< Epetra_CrsGraph > const & g, double beta=1.0e-6,
double alpha=1.0e-4) -> FiniteDifference
__init__(NOX::Epetra::FiniteDifference self, ParameterList printingParams, Teuchos::RCP< NOX::Epetra::Interface::Required > const & i,
Vector initialGuess, Teuchos::RCP< Epetra_CrsGraph > const & g, Teuchos::RCP< Epetra_Vector const > const & beta,
double alpha=1.0e-4) -> FiniteDifference
FiniteDifference::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.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
Reimplemented in PyTrilinos.NOX.Epetra.FiniteDifferenceColoring.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Apply | ( | self, | |
| args | |||
| ) |
Apply(FiniteDifference self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int FiniteDifference::Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const Return the result on an Epetra_Operator applied to an Epetra_MultiVector X in Y.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.ApplyInverse | ( | self, | |
| args | |||
| ) |
ApplyInverse(FiniteDifference self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int FiniteDifference::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.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Comm | ( | self, | |
| args | |||
| ) |
Comm(FiniteDifference self) -> Comm const Epetra_Comm & FiniteDifference::Comm() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.computeJacobian | ( | self, | |
| args | |||
| ) |
computeJacobian(FiniteDifference self, Epetra_Vector x, Operator Jac) -> bool computeJacobian(FiniteDifference self, Epetra_Vector x) -> bool bool FiniteDifference::computeJacobian(const Epetra_Vector &x) Compute Jacobian given the specified input vector, x. Returns true if computation was successful.
Reimplemented from PyTrilinos.NOX.Epetra.Interface.Jacobian.
Reimplemented in PyTrilinos.NOX.Epetra.FiniteDifferenceColoring.
| def PyTrilinos.NOX.Epetra.FiniteDifference.computePreconditioner | ( | self, | |
| args | |||
| ) |
computePreconditioner(FiniteDifference self, Epetra_Vector x, Operator Prec, ParameterList precParams=None) -> bool bool FiniteDifference::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.
Reimplemented from PyTrilinos.NOX.Epetra.Interface.Preconditioner.
| def PyTrilinos.NOX.Epetra.FiniteDifference.ExtractDiagonalCopy | ( | self, | |
| args | |||
| ) |
ExtractDiagonalCopy(FiniteDifference self, Epetra_Vector Diagonal) -> int int FiniteDifference::ExtractDiagonalCopy(Epetra_Vector &Diagonal) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.ExtractMyRowCopy | ( | self, | |
| args | |||
| ) |
ExtractMyRowCopy(FiniteDifference self, int MyRow, int Length, int & NumEntries, double * Values, int * Indices) -> int int FiniteDifference::ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Filled | ( | self, | |
| args | |||
| ) |
Filled(FiniteDifference self) -> bool bool FiniteDifference::Filled() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.getUnderlyingMatrix | ( | self, | |
| args | |||
| ) |
getUnderlyingMatrix(FiniteDifference self) -> CrsMatrix Epetra_CrsMatrix & FiniteDifference::getUnderlyingMatrix() const An accessor method for the underlying Epetra_CrsMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.HasNormInf | ( | self, | |
| args | |||
| ) |
HasNormInf(FiniteDifference self) -> bool bool FiniteDifference::HasNormInf() const Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.InvColSums | ( | self, | |
| args | |||
| ) |
InvColSums(FiniteDifference self, Epetra_Vector x) -> int int FiniteDifference::InvColSums(Epetra_Vector &x) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.InvRowSums | ( | self, | |
| args | |||
| ) |
InvRowSums(FiniteDifference self, Epetra_Vector x) -> int int FiniteDifference::InvRowSums(Epetra_Vector &x) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Label | ( | self, | |
| args | |||
| ) |
Label(FiniteDifference self) -> char const * const char * FiniteDifference::Label() const Returns a character string describing the name of the operator.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.LeftScale | ( | self, | |
| args | |||
| ) |
LeftScale(FiniteDifference self, Epetra_Vector x) -> int int FiniteDifference::LeftScale(const Epetra_Vector &x) See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.LowerTriangular | ( | self, | |
| args | |||
| ) |
LowerTriangular(FiniteDifference self) -> bool bool FiniteDifference::LowerTriangular() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Map | ( | self, | |
| args | |||
| ) |
Map(FiniteDifference self) -> BlockMap const Epetra_BlockMap & FiniteDifference::Map() const See Epetra_SrcDistObj documentation.
Reimplemented from PyTrilinos.NOX.Epetra.SrcDistObject.
| def PyTrilinos.NOX.Epetra.FiniteDifference.MaxNumEntries | ( | self, | |
| args | |||
| ) |
MaxNumEntries(FiniteDifference self) -> int int FiniteDifference::MaxNumEntries() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Multiply | ( | self, | |
| args | |||
| ) |
Multiply(FiniteDifference self, bool TransA, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int FiniteDifference::Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NormInf | ( | self, | |
| args | |||
| ) |
NormInf(FiniteDifference self) -> double double FiniteDifference::NormInf() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NormOne | ( | self, | |
| args | |||
| ) |
NormOne(FiniteDifference self) -> double double FiniteDifference::NormOne() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalCols | ( | self, | |
| args | |||
| ) |
NumGlobalCols(FiniteDifference self) -> int int FiniteDifference::NumGlobalCols() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalCols64 | ( | self, | |
| args | |||
| ) |
NumGlobalCols64(FiniteDifference self) -> long long long long FiniteDifference::NumGlobalCols64() const
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalDiagonals | ( | self, | |
| args | |||
| ) |
NumGlobalDiagonals(FiniteDifference self) -> int int FiniteDifference::NumGlobalDiagonals() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalDiagonals64 | ( | self, | |
| args | |||
| ) |
NumGlobalDiagonals64(FiniteDifference self) -> long long long long FiniteDifference::NumGlobalDiagonals64() const
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalNonzeros | ( | self, | |
| args | |||
| ) |
NumGlobalNonzeros(FiniteDifference self) -> int int FiniteDifference::NumGlobalNonzeros() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalNonzeros64 | ( | self, | |
| args | |||
| ) |
NumGlobalNonzeros64(FiniteDifference self) -> long long long long FiniteDifference::NumGlobalNonzeros64() const
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalRows | ( | self, | |
| args | |||
| ) |
NumGlobalRows(FiniteDifference self) -> int int FiniteDifference::NumGlobalRows() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumGlobalRows64 | ( | self, | |
| args | |||
| ) |
NumGlobalRows64(FiniteDifference self) -> long long long long FiniteDifference::NumGlobalRows64() const
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumMyCols | ( | self, | |
| args | |||
| ) |
NumMyCols(FiniteDifference self) -> int int FiniteDifference::NumMyCols() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumMyDiagonals | ( | self, | |
| args | |||
| ) |
NumMyDiagonals(FiniteDifference self) -> int int FiniteDifference::NumMyDiagonals() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumMyNonzeros | ( | self, | |
| args | |||
| ) |
NumMyNonzeros(FiniteDifference self) -> int int FiniteDifference::NumMyNonzeros() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumMyRowEntries | ( | self, | |
| args | |||
| ) |
NumMyRowEntries(FiniteDifference self, int MyRow, int & NumEntries) -> int int FiniteDifference::NumMyRowEntries(int MyRow, int &NumEntries) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.NumMyRows | ( | self, | |
| args | |||
| ) |
NumMyRows(FiniteDifference self) -> int int FiniteDifference::NumMyRows() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.OperatorDomainMap | ( | self, | |
| args | |||
| ) |
OperatorDomainMap(FiniteDifference self) -> Map const Epetra_Map & FiniteDifference::OperatorDomainMap() const Returns the Epetra_BlockMap object associated with the domain of this matrix operator.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.OperatorRangeMap | ( | self, | |
| args | |||
| ) |
OperatorRangeMap(FiniteDifference self) -> Map const Epetra_Map & FiniteDifference::OperatorRangeMap() const Returns the Epetra_BlockMap object associated with the range of this matrix operator.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.RightScale | ( | self, | |
| args | |||
| ) |
RightScale(FiniteDifference self, Epetra_Vector x) -> int int FiniteDifference::RightScale(const Epetra_Vector &x) See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.RowMatrixColMap | ( | self, | |
| args | |||
| ) |
RowMatrixColMap(FiniteDifference self) -> Map const Epetra_Map & FiniteDifference::RowMatrixColMap() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.RowMatrixImporter | ( | self, | |
| args | |||
| ) |
RowMatrixImporter(FiniteDifference self) -> Import const Epetra_Import * FiniteDifference::RowMatrixImporter() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.RowMatrixRowMap | ( | self, | |
| args | |||
| ) |
RowMatrixRowMap(FiniteDifference self) -> Map const Epetra_Map & FiniteDifference::RowMatrixRowMap() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.setDifferenceMethod | ( | self, | |
| args | |||
| ) |
setDifferenceMethod(FiniteDifference self, NOX::Epetra::FiniteDifference::DifferenceType type) void FiniteDifference::setDifferenceMethod(DifferenceType type) Set the type of perturbation method used (default is Forward)
| def PyTrilinos.NOX.Epetra.FiniteDifference.setGroupForComputeF | ( | self, | |
| args | |||
| ) |
setGroupForComputeF(FiniteDifference self, Group group) void FiniteDifference::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.
| def PyTrilinos.NOX.Epetra.FiniteDifference.SetUseTranspose | ( | self, | |
| args | |||
| ) |
SetUseTranspose(FiniteDifference self, bool UseTranspose) -> int int FiniteDifference::SetUseTranspose(bool UseTranspose) If set true, the transpose of this operator will be applied.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
| def PyTrilinos.NOX.Epetra.FiniteDifference.Solve | ( | self, | |
| args | |||
| ) |
Solve(FiniteDifference self, bool Upper, bool Trans, bool UnitDiagonal, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int FiniteDifference::Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.UpperTriangular | ( | self, | |
| args | |||
| ) |
UpperTriangular(FiniteDifference self) -> bool bool FiniteDifference::UpperTriangular() const See Epetra_RowMatrix documentation.
Reimplemented from PyTrilinos.NOX.Epetra.RowMatrix.
| def PyTrilinos.NOX.Epetra.FiniteDifference.UseTranspose | ( | self, | |
| args | |||
| ) |
UseTranspose(FiniteDifference self) -> bool bool FiniteDifference::UseTranspose() const Returns the current use transpose setting.
Reimplemented from PyTrilinos.NOX.Epetra.Operator.
1.7.6.1