xref: /petsc/src/mat/matfd/fdmatrix.c (revision c8a9c62294a3909d4cb0ca2bea407dacfcbcf516)
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 inserting entries of Jacobian matrix.
195 
196    Logically Collective on MatFDColoring
197 
198    Input Parameters:
199 +  coloring - the coloring context
200 .  brows - number of rows in the block
201 -  bcols - number of columns in the block
202 
203    Level: intermediate
204 
205 .keywords: Mat, coloring
206 
207 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
208 
209 @*/
210 PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
211 {
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
214   PetscValidLogicalCollectiveInt(matfd,brows,2);
215   PetscValidLogicalCollectiveInt(matfd,bcols,3);
216   if (brows != PETSC_DEFAULT) matfd->brows = brows;
217   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
218   PetscFunctionReturn(0);
219 }
220 
221 #undef __FUNCT__
222 #define __FUNCT__ "MatFDColoringSetUp"
223 /*@
224    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
225 
226    Collective on Mat
227 
228    Input Parameters:
229 +  mat - the matrix containing the nonzero structure of the Jacobian
230 .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
231 -  color - the matrix coloring context
232 
233    Level: beginner
234 
235 .keywords: MatFDColoring, setup
236 
237 .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
238 @*/
239 PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
240 {
241   PetscErrorCode ierr;
242 
243   PetscFunctionBegin;
244   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
245   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
246   if (color->setupcalled) PetscFunctionReturn(0);
247 
248   if (mat->ops->fdcoloringsetup) {
249     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
250   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
251 
252   color->setupcalled = PETSC_TRUE;
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "MatFDColoringGetFunction"
258 /*@C
259    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
260 
261    Not Collective
262 
263    Input Parameters:
264 .  coloring - the coloring context
265 
266    Output Parameters:
267 +  f - the function
268 -  fctx - the optional user-defined function context
269 
270    Level: intermediate
271 
272 .keywords: Mat, Jacobian, finite differences, set, function
273 
274 .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
275 
276 @*/
277 PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
278 {
279   PetscFunctionBegin;
280   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
281   if (f) *f = matfd->f;
282   if (fctx) *fctx = matfd->fctx;
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "MatFDColoringSetFunction"
288 /*@C
289    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
290 
291    Logically Collective on MatFDColoring
292 
293    Input Parameters:
294 +  coloring - the coloring context
295 .  f - the function
296 -  fctx - the optional user-defined function context
297 
298    Calling sequence of (*f) function:
299     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
300     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
301 
302    Level: advanced
303 
304    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
305      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
306      calling MatFDColoringApply()
307 
308    Fortran Notes:
309     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
310   be used without SNES or within the SNES solvers.
311 
312 .keywords: Mat, Jacobian, finite differences, set, function
313 
314 .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
315 
316 @*/
317 PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
318 {
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
321   matfd->f    = f;
322   matfd->fctx = fctx;
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "MatFDColoringSetFromOptions"
328 /*@
329    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
330    the options database.
331 
332    Collective on MatFDColoring
333 
334    The Jacobian, F'(u), is estimated with the differencing approximation
335 .vb
336        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
337        h = error_rel*u[i]                 if  abs(u[i]) > umin
338          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
339        dx_{i} = (0, ... 1, .... 0)
340 .ve
341 
342    Input Parameter:
343 .  coloring - the coloring context
344 
345    Options Database Keys:
346 +  -mat_fd_coloring_err <err> - Sets <err> (square root
347            of relative error in the function)
348 .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
349 .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
350 .  -mat_fd_coloring_view - Activates basic viewing
351 .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
352 -  -mat_fd_coloring_view draw - Activates drawing
353 
354     Level: intermediate
355 
356 .keywords: Mat, finite differences, parameters
357 
358 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
359 
360 @*/
361 PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
362 {
363   PetscErrorCode ierr;
364   PetscBool      flg;
365   char           value[3];
366 
367   PetscFunctionBegin;
368   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
369 
370   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
371   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
372   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
373   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
374   if (flg) {
375     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
376     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
377     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
378   }
379   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
380   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
381   if (flg && matfd->bcols > matfd->ncolors) {
382     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
383     matfd->bcols = matfd->ncolors;
384   }
385 
386   /* process any options handlers added with PetscObjectAddOptionsHandler() */
387   ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
388   PetscOptionsEnd();CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "MatFDColoringViewFromOptions"
394 PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
395 {
396   PetscErrorCode    ierr;
397   PetscBool         flg;
398   PetscViewer       viewer;
399   PetscViewerFormat format;
400 
401   PetscFunctionBegin;
402   if (prefix) {
403     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
404   } else {
405     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
406   }
407   if (flg) {
408     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
409     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
410     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
411     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
412   }
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "MatFDColoringCreate"
418 /*@
419    MatFDColoringCreate - Creates a matrix coloring context for finite difference
420    computation of Jacobians.
421 
422    Collective on Mat
423 
424    Input Parameters:
425 +  mat - the matrix containing the nonzero structure of the Jacobian
426 -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
427 
428     Output Parameter:
429 .   color - the new coloring context
430 
431     Level: intermediate
432 
433 .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
434           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
435           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
436 @*/
437 PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
438 {
439   MatFDColoring  c;
440   MPI_Comm       comm;
441   PetscErrorCode ierr;
442   PetscInt       M,N;
443 
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
446   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
447   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
448   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
449   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
450   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
451   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
452 
453   c->ctype = iscoloring->ctype;
454 
455   if (mat->ops->fdcoloringcreate) {
456     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
457   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
458 
459   ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
460   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
461   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
462   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
463 
464   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
465   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
466   c->currentcolor = -1;
467   c->htype        = "wp";
468   c->fset         = PETSC_FALSE;
469   c->setupcalled  = PETSC_FALSE;
470 
471   *color = c;
472   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
473   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "MatFDColoringDestroy"
479 /*@
480     MatFDColoringDestroy - Destroys a matrix coloring context that was created
481     via MatFDColoringCreate().
482 
483     Collective on MatFDColoring
484 
485     Input Parameter:
486 .   c - coloring context
487 
488     Level: intermediate
489 
490 .seealso: MatFDColoringCreate()
491 @*/
492 PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
493 {
494   PetscErrorCode ierr;
495   PetscInt       i;
496 
497   PetscFunctionBegin;
498   if (!*c) PetscFunctionReturn(0);
499   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
500 
501   for (i=0; i<(*c)->ncolors; i++) {
502     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
503   }
504   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
505   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
506   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
507   ierr = PetscFree((*c)->matentry);CHKERRQ(ierr);
508   ierr = PetscFree((*c)->dy);CHKERRQ(ierr);
509   if ((*c)->vscale) {ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);}
510   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
511   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
512   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
513   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
514   PetscFunctionReturn(0);
515 }
516 
517 #undef __FUNCT__
518 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
519 /*@C
520     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
521       that are currently being perturbed.
522 
523     Not Collective
524 
525     Input Parameters:
526 .   coloring - coloring context created with MatFDColoringCreate()
527 
528     Output Parameters:
529 +   n - the number of local columns being perturbed
530 -   cols - the column indices, in global numbering
531 
532    Level: intermediate
533 
534 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
535 
536 .keywords: coloring, Jacobian, finite differences
537 @*/
538 PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
539 {
540   PetscFunctionBegin;
541   if (coloring->currentcolor >= 0) {
542     *n    = coloring->ncolumns[coloring->currentcolor];
543     *cols = coloring->columns[coloring->currentcolor];
544   } else {
545     *n = 0;
546   }
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "MatFDColoringApply"
552 /*@
553     MatFDColoringApply - Given a matrix for which a MatFDColoring context
554     has been created, computes the Jacobian for a function via finite differences.
555 
556     Collective on MatFDColoring
557 
558     Input Parameters:
559 +   mat - location to store Jacobian
560 .   coloring - coloring context created with MatFDColoringCreate()
561 .   x1 - location at which Jacobian is to be computed
562 -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
563 
564     Options Database Keys:
565 +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
566 .    -mat_fd_coloring_view - Activates basic viewing or coloring
567 .    -mat_fd_coloring_view draw - Activates drawing of coloring
568 -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
569 
570     Level: intermediate
571 
572 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
573 
574 .keywords: coloring, Jacobian, finite differences
575 @*/
576 PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
577 {
578   PetscErrorCode ierr;
579   PetscBool      flg = PETSC_FALSE;
580 
581   PetscFunctionBegin;
582   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
583   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
584   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
585   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
586   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
587   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
588 
589   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
590   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
591   if (flg) {
592     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
593   } else {
594     PetscBool assembled;
595     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
596     if (assembled) {
597       ierr = MatZeroEntries(J);CHKERRQ(ierr);
598     }
599   }
600 
601   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
602   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
603   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
604   PetscFunctionReturn(0);
605 }
606