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