#include <NOX_Epetra_FiniteDifferenceColoring.H>
Inheritance diagram for NOX::Epetra::FiniteDifferenceColoring:
Public Member Functions | |
| FiniteDifferenceColoring (Interface &i, const Epetra_Vector &initialGuess, Epetra_MapColoring &colorMap, vector< Epetra_IntVector > &columns, double beta=1.0e-6, double alpha=1.0e-4) | |
| Constructor. | |
| FiniteDifferenceColoring (Interface &i, const Epetra_Vector &initialGuess, Epetra_CrsGraph &rawGraph, Epetra_MapColoring &colorMap, vector< Epetra_IntVector > &columns, double beta=1.0e-6, double alpha=1.0e-4) | |
| Constructor. | |
| virtual | ~FiniteDifferenceColoring () |
| Pure virtual destructor. | |
| 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. | |
Protected Attributes | |
| const Epetra_MapColoring * | colorMap |
| Color Map created by external algorithm. | |
| vector< Epetra_IntVector > * | columns |
| vector of Epetra_IntVectors containing columns corresponding to a given row and color | |
| int | numColors |
| Number of colors in Color Map. | |
| int * | colorList |
| List of colors in Color Map. | |
| Epetra_Map * | cMap |
| Coloring Map created by external algorithm. | |
| Epetra_Import * | Importer |
| Importer needed for mapping Color Map to unColored Map. | |
| Epetra_Vector * | colorVect |
| Color vector based on Color Map containing perturbations. | |
| Epetra_Vector * | betaColorVect |
| Color vector based on Color Map containing beta value(s). | |
| Epetra_Vector | mappedColorVect |
| Color vector based on unColorred Map containing perturbations. | |
| Epetra_Vector * | xCol_perturb |
| Perturbed solution vector based on column map. | |
| const Epetra_BlockMap * | columnMap |
| Overlap Map (Column Map of Matrix Graph) needed for parallel. | |
| Epetra_Import * | rowColImporter |
| An Import object needed in parallel to map from row-space to column-space. | |
The Jacobian entries are calculated via 1st or 2nd order finite differencing. This requires
or
calls to computeF(), respectively, where
is the number of colors.
where
is the Jacobian,
is the function evaluation,
is the solution vector, and
is a small perturbation to the
entry.
Instead of perturbing each
problem degrees of freedom sequentially and then evaluating all
functions for each perturbation, coloring allows several degrees of freedom (all belonging to the same color) to be perturbed at the same time. This reduces the total number of function evaluations needed to compute
from
as is required using FiniteDifference to
, often representing substantial computational savings.
Coloring is based on a user-supplied color map generated using an appropriate algorithm, eg greedy-algorithm - Y. Saad, "Iterative Methods for Sparse Linear Systems, 2nd ed.," chp. 3, SIAM, 2003.. Use can be made of the coloring algorithm provided by the EpetraExt package in Trilinos. The 1Dfem_nonlinearColoring and Brusselator example problems located in the nox/epetra-examples subdirectory demonstrate use of the EpetraExt package, and the 1Dfem_nonlinearColoring directory also contains a stand-alone coloring algorithm very similar to that in EpetraExt.
The perturbation,
, is calculated using the following equation:
where
is a scalar value (defaults to 1.0e-4) and
is another scalar (defaults to 1.0e-6).
Since both FiniteDifferenceColoring and FiniteDifference inherit from the Epetra_RowMatrix class, they can be used as preconditioning matrices for AztecOO preconditioners.
As for FiniteDifference, 1st order accurate Forward and Backward differences as well as 2nd order accurate Centered difference can be specified using setDifferenceMethod with the appropriate enumerated type passed as the argument.
Using FiniteDifferenceColoring in Parallel
This class can be used to numerically approximate
for distributed parallel problems provided a very minimal of conditions are met. Namely, an Epetra_CrsGraph reflecting a row-based distribution must be used in the constructor. Next, the application must be aware that the perturbed vector used to compute residuals has already been scattered to a form consistent with the column space of the Epetra_CrsGraph. In practice, this means that the perturbed vector used by computeF() has already been scattered to a ghosted or overlapped state. The application should then not perform this step but rather simply use the vector provided with the possible exception of requiring a local index reordering to bring the column-space based vector in sync with a potentially different ghosted index ordering. See the Brusselator and 1Dfem_nonlinearColoring example problems for details.
1.3.9.1