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