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