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