NOX::Epetra::FiniteDifferenceColoring Class Reference

Concrete implementation for creating an Epetra_RowMatrix Jacobian via finite differencing of the residual using coloring. More...

#include <NOX_Epetra_FiniteDifferenceColoring.H>

Inheritance diagram for NOX::Epetra::FiniteDifferenceColoring:

[legend]
Collaboration diagram for NOX::Epetra::FiniteDifferenceColoring:
[legend]
List of all members.

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_MapColoringcolorMap
 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_MapcMap
 Coloring Map created by external algorithm.
Epetra_ImportImporter
 Importer needed for mapping Color Map to unColored Map.
Epetra_VectorcolorVect
 Color vector based on Color Map containing perturbations.
Epetra_VectorbetaColorVect
 Color vector based on Color Map containing beta value(s).
Epetra_Vector mappedColorVect
 Color vector based on unColorred Map containing perturbations.
Epetra_VectorxCol_perturb
 Perturbed solution vector based on column map.
const Epetra_BlockMapcolumnMap
 Overlap Map (Column Map of Matrix Graph) needed for parallel.
Epetra_ImportrowColImporter
 An Import object needed in parallel to map from row-space to column-space.

Detailed Description

Concrete implementation for creating an Epetra_RowMatrix Jacobian via finite differencing of the residual using coloring.

The Jacobian entries are calculated via 1st or 2nd order finite differencing. This requires $ N + 1 $ or $ 2N + 1 $ calls to computeF(), respectively, where $ N $ is the number of colors.

\[ 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.

Instead of perturbing each $ N_{dof} $ problem degrees of freedom sequentially and then evaluating all $ N_{dof} $ 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 $\mathbf{J}$ from $ N_{dof}^2 $ as is required using FiniteDifference to $ N\cdot N_{dof} $, 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, $ \delta $, is calculated using the following equation:

\[ \delta = \alpha * | x_j | + \beta \]

where $ \alpha $ is a scalar value (defaults to 1.0e-4) and $ \beta $ 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 $ J $ 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.


The documentation for this class was generated from the following files:
Generated on Thu Sep 18 12:42:27 2008 for NOX by doxygen 1.3.9.1