PETSc provides some tools to help approximate the Jacobian matrices efficiently via finite differences. These tools are intended for use in certain situations where one is unable to compute Jacobian matrices analytically, and matrix-free methods do not work well without a preconditioner, due to very poor conditioning. The approximation requires several steps:
ISColoring iscoloring; MatFDColoring fdcoloring; MatStructure str;/* This initializes the nonzero structure of the Jacobian. This is artificial because clearly if we had a routine to compute the Jacobian we wouldn't need to use finite differences. */ FormJacobian(snes,x,&J,&J,&str,&user);
/* Color the matrix, i.e. determine groups of columns that share no common rows. These columns in the Jacobian can all be computed simulataneously. */ MatGetColoring(J,MATCOLORING_SL,&iscoloring);
/* Create the data structure that SNESDefaultComputeJacobianColor() uses to compute the actual Jacobians via finite differences. */ MatFDColoringCreate(J,iscoloring,&fdcoloring); ISColoringDestroy(iscoloring); MatFDColoringSetFromOptions(fdcoloring);
/* Tell SNES to use the routine SNESDefaultComputeJacobianColor() to compute Jacobians. */ SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fdcoloring);
Of course, we are cheating a bit. If we do not have an analytic formula for computing the Jacobian, then how do we know what its nonzero structure is so that it may be colored? Determining the structure is problem dependent, but fortunately, for most grid-based problems (the class of problems for which PETSc is designed) if one knows the stencil used for the nonlinear function one can usually fairly easily obtain an estimate of the location of nonzeros in the matrix.
One need not necessarily use the routine MatGetColoring() to determine a coloring. For example, if a grid can be colored directly (without using the associated matrix), then that coloring can be provided to MatFDColoringCreate(). Note that the user must always preset the nonzero structure in the matrix regardless of which coloring routine is used.
For sequential matrices PETSc provides three matrix coloring routines from the
MINPACK package [14]: smallest-last ( sl), largest-first ( lf),
and incidence-degree ( id). These colorings, as well as the ``natural'' coloring
for which each column has its own unique color, may be accessed with the command line options
-mat_coloring [sl,id,lf,natural]Alternatively, one can set a coloring type of COLORING_SL, COLORING_ID, COLORING_LF, or COLORING_NATURAL when calling MatGetColoring().
As for the matrix-free computation of Jacobians (see Section
Matrix-Free Methods
), two parameters affect the accuracy of the
finite difference Jacobian approximation. These are set with the command
ierr = MatFDColoringSetParameters(MatFDColoring fdcoloring,double rerror,double umin);The parameter rerror is the square root of the relative error in the function evaluations, erel; the default is 10-8 , which assumes that the functions are evaluated to full double-precision accuracy. The second parameter, umin, is a bit more involved; its default is 10e-8 . Column i of the Jacobian matrix (denoted by F_:i) is approximated by the formula
where h is computed via
These parameters may be set from the options database with
-mat_fd_coloring_err <err> -mat_fd_coloring_umin <umin>
Note that the MatGetColoring() routine currently works only on sequential routines. Extensions may be forthcoming. However, if one can compute the coloring iscoloring some other way, the routine MatFDColoringCreate() does work in parallel. An example of this for 2D distributed arrays is given below that uses the utility routine DAGetColoring().
ierr = DAGetColoring(da,&iscoloring,&J); ierr = MatFDColoringCreate(J,iscoloring,&fdcoloring); ierr = MatFDColoringSetFromOptions(fdcoloring); ierr = ISColoringDestroy(iscoloring);Note that the routine MatFDColoringCreate() currently is only supported for the AIJ matrix format.