xref: /petsc/src/mat/matfd/fdmatrix.c (revision c688d0420b4e513ff34944d1e1ad7d4e50aafa8d)
1 
2 /*
3    This is where the abstract matrix operations are defined that are
4   used for finite difference computations of Jacobians using coloring.
5 */
6 
7 #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
8 #include <petsc/private/isimpl.h>
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "MatFDColoringSetF"
12 PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
13 {
14   PetscErrorCode ierr;
15 
16   PetscFunctionBegin;
17   if (F) {
18     ierr     = VecCopy(F,fd->w1);CHKERRQ(ierr);
19     fd->fset = PETSC_TRUE;
20   } else {
21     fd->fset = PETSC_FALSE;
22   }
23   PetscFunctionReturn(0);
24 }
25 
26 #include <petscdraw.h>
27 #undef __FUNCT__
28 #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
29 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
30 {
31   MatFDColoring  fd = (MatFDColoring)Aa;
32   PetscErrorCode ierr;
33   PetscInt       i,j,nz,row;
34   PetscReal      x,y;
35   MatEntry       *Jentry=fd->matentry;
36 
37   PetscFunctionBegin;
38   /* loop over colors  */
39   nz = 0;
40   for (i=0; i<fd->ncolors; i++) {
41     for (j=0; j<fd->nrows[i]; j++) {
42       row = Jentry[nz].row;
43       y   = fd->M - row - fd->rstart;
44       x   = (PetscReal)Jentry[nz++].col;
45       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
46     }
47   }
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "MatFDColoringView_Draw"
53 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
54 {
55   PetscErrorCode ierr;
56   PetscBool      isnull;
57   PetscDraw      draw;
58   PetscReal      xr,yr,xl,yl,h,w;
59 
60   PetscFunctionBegin;
61   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
62   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
63   if (isnull) PetscFunctionReturn(0);
64 
65   xr   = fd->N; yr  = fd->M; h = yr/10.0; w = xr/10.0;
66   xr  += w;     yr += h;    xl = -w;     yl = -h;
67   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
68   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
69   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
70   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr);
71   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "MatFDColoringView"
77 /*@C
78    MatFDColoringView - Views a finite difference coloring context.
79 
80    Collective on MatFDColoring
81 
82    Input  Parameters:
83 +  c - the coloring context
84 -  viewer - visualization context
85 
86    Level: intermediate
87 
88    Notes:
89    The available visualization contexts include
90 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
91 .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
92         output where only the first processor opens
93         the file.  All other processors send their
94         data to the first processor to print.
95 -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
96 
97    Notes:
98      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
99    involves more than 33 then some seemingly identical colors are displayed making it look
100    like an illegal coloring. This is just a graphical artifact.
101 
102 .seealso: MatFDColoringCreate()
103 
104 .keywords: Mat, finite differences, coloring, view
105 @*/
106 PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
107 {
108   PetscErrorCode    ierr;
109   PetscInt          i,j;
110   PetscBool         isdraw,iascii;
111   PetscViewerFormat format;
112 
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
115   if (!viewer) {
116     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
117   }
118   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
119   PetscCheckSameComm(c,1,viewer,2);
120 
121   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
122   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
123   if (isdraw) {
124     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
125   } else if (iascii) {
126     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
127     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr);
128     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",(double)c->umin);CHKERRQ(ierr);
129     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
130 
131     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
132     if (format != PETSC_VIEWER_ASCII_INFO) {
133       PetscInt row,col,nz;
134       nz = 0;
135       for (i=0; i<c->ncolors; i++) {
136         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
137         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
138         for (j=0; j<c->ncolumns[i]; j++) {
139           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
140         }
141         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
142         if (c->matentry) {
143           for (j=0; j<c->nrows[i]; j++) {
144             row  = c->matentry[nz].row;
145             col  = c->matentry[nz++].col;
146             ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);CHKERRQ(ierr);
147           }
148         }
149       }
150     }
151     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
152   }
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "MatFDColoringSetParameters"
158 /*@
159    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
160    a Jacobian matrix using finite differences.
161 
162    Logically Collective on MatFDColoring
163 
164    The Jacobian is estimated with the differencing approximation
165 .vb
166        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
167        h = error_rel*u[i]                 if  abs(u[i]) > umin
168          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
169        dx_{i} = (0, ... 1, .... 0)
170 .ve
171 
172    Input Parameters:
173 +  coloring - the coloring context
174 .  error_rel - relative error
175 -  umin - minimum allowable u-value magnitude
176 
177    Level: advanced
178 
179 .keywords: Mat, finite differences, coloring, set, parameters
180 
181 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
182 
183 @*/
184 PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
185 {
186   PetscFunctionBegin;
187   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
188   PetscValidLogicalCollectiveReal(matfd,error,2);
189   PetscValidLogicalCollectiveReal(matfd,umin,3);
190   if (error != PETSC_DEFAULT) matfd->error_rel = error;
191   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "MatFDColoringSetBlockSize"
197 /*@
198    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
199 
200    Logically Collective on MatFDColoring
201 
202    Input Parameters:
203 +  coloring - the coloring context
204 .  brows - number of rows in the block
205 -  bcols - number of columns in the block
206 
207    Level: intermediate
208 
209 .keywords: Mat, coloring
210 
211 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
212 
213 @*/
214 PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
215 {
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
218   PetscValidLogicalCollectiveInt(matfd,brows,2);
219   PetscValidLogicalCollectiveInt(matfd,bcols,3);
220   if (brows != PETSC_DEFAULT) matfd->brows = brows;
221   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "MatFDColoringSetUp"
227 /*@
228    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
229 
230    Collective on Mat
231 
232    Input Parameters:
233 +  mat - the matrix containing the nonzero structure of the Jacobian
234 .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
235 -  color - the matrix coloring context
236 
237    Level: beginner
238 
239 .keywords: MatFDColoring, setup
240 
241 .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
242 @*/
243 PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
244 {
245   PetscErrorCode ierr;
246 
247   PetscFunctionBegin;
248   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
249   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
250   if (color->setupcalled) PetscFunctionReturn(0);
251 
252   ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
253   if (mat->ops->fdcoloringsetup) {
254     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
255   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
256 
257   color->setupcalled = PETSC_TRUE;
258    ierr   = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "MatFDColoringGetFunction"
264 /*@C
265    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
266 
267    Not Collective
268 
269    Input Parameters:
270 .  coloring - the coloring context
271 
272    Output Parameters:
273 +  f - the function
274 -  fctx - the optional user-defined function context
275 
276    Level: intermediate
277 
278 .keywords: Mat, Jacobian, finite differences, set, function
279 
280 .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
281 
282 @*/
283 PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
284 {
285   PetscFunctionBegin;
286   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
287   if (f) *f = matfd->f;
288   if (fctx) *fctx = matfd->fctx;
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "MatFDColoringSetFunction"
294 /*@C
295    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
296 
297    Logically Collective on MatFDColoring
298 
299    Input Parameters:
300 +  coloring - the coloring context
301 .  f - the function
302 -  fctx - the optional user-defined function context
303 
304    Calling sequence of (*f) function:
305     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
306     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
307 
308    Level: advanced
309 
310    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
311      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
312      calling MatFDColoringApply()
313 
314    Fortran Notes:
315     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
316   be used without SNES or within the SNES solvers.
317 
318 .keywords: Mat, Jacobian, finite differences, set, function
319 
320 .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
321 
322 @*/
323 PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
324 {
325   PetscFunctionBegin;
326   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
327   matfd->f    = f;
328   matfd->fctx = fctx;
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "MatFDColoringSetFromOptions"
334 /*@
335    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
336    the options database.
337 
338    Collective on MatFDColoring
339 
340    The Jacobian, F'(u), is estimated with the differencing approximation
341 .vb
342        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
343        h = error_rel*u[i]                 if  abs(u[i]) > umin
344          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
345        dx_{i} = (0, ... 1, .... 0)
346 .ve
347 
348    Input Parameter:
349 .  coloring - the coloring context
350 
351    Options Database Keys:
352 +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
353 .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
354 .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
355 .  -mat_fd_coloring_view - Activates basic viewing
356 .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
357 -  -mat_fd_coloring_view draw - Activates drawing
358 
359     Level: intermediate
360 
361 .keywords: Mat, finite differences, parameters
362 
363 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
364 
365 @*/
366 PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
367 {
368   PetscErrorCode ierr;
369   PetscBool      flg;
370   char           value[3];
371 
372   PetscFunctionBegin;
373   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
374 
375   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
376   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
377   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
378   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
379   if (flg) {
380     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
381     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
382     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
383   }
384   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
385   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
386   if (flg && matfd->bcols > matfd->ncolors) {
387     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
388     matfd->bcols = matfd->ncolors;
389   }
390 
391   /* process any options handlers added with PetscObjectAddOptionsHandler() */
392   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr);
393   PetscOptionsEnd();CHKERRQ(ierr);
394   PetscFunctionReturn(0);
395 }
396 
397 #undef __FUNCT__
398 #define __FUNCT__ "MatFDColoringSetType"
399 /*@C
400    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
401 
402    Collective on MatFDColoring
403 
404    Input Parameters:
405 +  coloring - the coloring context
406 -  type - either MATMFFD_WP or MATMFFD_DS
407 
408    Options Database Keys:
409 .  -mat_fd_type - "wp" or "ds"
410 
411    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
412          but the process of computing the entries is the same as as with the MatMFFD operation so we should reuse the names instead of
413          introducing another one.
414 
415    Level: intermediate
416 
417 .keywords: Mat, finite differences, parameters
418 
419 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
420 
421 @*/
422 PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
423 {
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
426   /*
427      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
428      and this function is being provided as patch for a release so it shouldn't change the implementaton
429   */
430   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
431   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
432   else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
433   PetscFunctionReturn(0);
434 }
435 
436 #undef __FUNCT__
437 #define __FUNCT__ "MatFDColoringViewFromOptions"
438 PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
439 {
440   PetscErrorCode    ierr;
441   PetscBool         flg;
442   PetscViewer       viewer;
443   PetscViewerFormat format;
444 
445   PetscFunctionBegin;
446   if (prefix) {
447     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
448   } else {
449     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
450   }
451   if (flg) {
452     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
453     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
454     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
455     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
456   }
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "MatFDColoringCreate"
462 /*@
463    MatFDColoringCreate - Creates a matrix coloring context for finite difference
464    computation of Jacobians.
465 
466    Collective on Mat
467 
468    Input Parameters:
469 +  mat - the matrix containing the nonzero structure of the Jacobian
470 -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
471 
472     Output Parameter:
473 .   color - the new coloring context
474 
475     Level: intermediate
476 
477 .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
478           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
479           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
480 @*/
481 PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
482 {
483   MatFDColoring  c;
484   MPI_Comm       comm;
485   PetscErrorCode ierr;
486   PetscInt       M,N;
487 
488   PetscFunctionBegin;
489   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
490   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
491   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
492   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
493   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
494   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
495   ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
496 
497   c->ctype = iscoloring->ctype;
498 
499   if (mat->ops->fdcoloringcreate) {
500     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
501   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
502 
503   ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
504   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
505   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
506   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
507 
508   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
509   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
510   c->currentcolor = -1;
511   c->htype        = "wp";
512   c->fset         = PETSC_FALSE;
513   c->setupcalled  = PETSC_FALSE;
514 
515   *color = c;
516   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
517   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
518   PetscFunctionReturn(0);
519 }
520 
521 #undef __FUNCT__
522 #define __FUNCT__ "MatFDColoringDestroy"
523 /*@
524     MatFDColoringDestroy - Destroys a matrix coloring context that was created
525     via MatFDColoringCreate().
526 
527     Collective on MatFDColoring
528 
529     Input Parameter:
530 .   c - coloring context
531 
532     Level: intermediate
533 
534 .seealso: MatFDColoringCreate()
535 @*/
536 PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
537 {
538   PetscErrorCode ierr;
539   PetscInt       i;
540   MatFDColoring  color = *c;
541 
542   PetscFunctionBegin;
543   if (!*c) PetscFunctionReturn(0);
544   if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);}
545 
546   for (i=0; i<color->ncolors; i++) {
547     ierr = PetscFree(color->columns[i]);CHKERRQ(ierr);
548   }
549   ierr = PetscFree(color->ncolumns);CHKERRQ(ierr);
550   ierr = PetscFree(color->columns);CHKERRQ(ierr);
551   ierr = PetscFree(color->nrows);CHKERRQ(ierr);
552   if (color->htype[0] == 'w') {
553     ierr = PetscFree(color->matentry2);CHKERRQ(ierr);
554   } else {
555     ierr = PetscFree(color->matentry);CHKERRQ(ierr);
556   }
557   ierr = PetscFree(color->dy);CHKERRQ(ierr);
558   if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);}
559   ierr = VecDestroy(&color->w1);CHKERRQ(ierr);
560   ierr = VecDestroy(&color->w2);CHKERRQ(ierr);
561   ierr = VecDestroy(&color->w3);CHKERRQ(ierr);
562   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 
566 #undef __FUNCT__
567 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
568 /*@C
569     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
570       that are currently being perturbed.
571 
572     Not Collective
573 
574     Input Parameters:
575 .   coloring - coloring context created with MatFDColoringCreate()
576 
577     Output Parameters:
578 +   n - the number of local columns being perturbed
579 -   cols - the column indices, in global numbering
580 
581    Level: intermediate
582 
583 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
584 
585 .keywords: coloring, Jacobian, finite differences
586 @*/
587 PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
588 {
589   PetscFunctionBegin;
590   if (coloring->currentcolor >= 0) {
591     *n    = coloring->ncolumns[coloring->currentcolor];
592     *cols = coloring->columns[coloring->currentcolor];
593   } else {
594     *n = 0;
595   }
596   PetscFunctionReturn(0);
597 }
598 
599 #undef __FUNCT__
600 #define __FUNCT__ "MatFDColoringApply"
601 /*@
602     MatFDColoringApply - Given a matrix for which a MatFDColoring context
603     has been created, computes the Jacobian for a function via finite differences.
604 
605     Collective on MatFDColoring
606 
607     Input Parameters:
608 +   mat - location to store Jacobian
609 .   coloring - coloring context created with MatFDColoringCreate()
610 .   x1 - location at which Jacobian is to be computed
611 -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
612 
613     Options Database Keys:
614 +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
615 .    -mat_fd_coloring_view - Activates basic viewing or coloring
616 .    -mat_fd_coloring_view draw - Activates drawing of coloring
617 -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
618 
619     Level: intermediate
620 
621 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
622 
623 .keywords: coloring, Jacobian, finite differences
624 @*/
625 PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
626 {
627   PetscErrorCode ierr;
628   PetscBool      flg = PETSC_FALSE;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
632   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
633   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
634   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
635   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
636   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
637 
638   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
639   ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
640   if (flg) {
641     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
642   } else {
643     PetscBool assembled;
644     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
645     if (assembled) {
646       ierr = MatZeroEntries(J);CHKERRQ(ierr);
647     }
648   }
649 
650   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
651   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr);
652   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
653   if (!coloring->viewed) {
654     ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr);
655     coloring->viewed = PETSC_TRUE;
656   }
657   PetscFunctionReturn(0);
658 }
659