xref: /petsc/src/mat/matfd/fdmatrix.c (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
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   PetscBool      flg;
340   char           value[3];
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
344 
345   PetscObjectOptionsBegin((PetscObject)matfd);
346   PetscCall(PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,NULL));
347   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,NULL));
348   PetscCall(PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,sizeof(value),&flg));
349   if (flg) {
350     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
351     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
352     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
353   }
354   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL));
355   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg));
356   if (flg && matfd->bcols > matfd->ncolors) {
357     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
358     matfd->bcols = matfd->ncolors;
359   }
360 
361   /* process any options handlers added with PetscObjectAddOptionsHandler() */
362   PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd));
363   PetscOptionsEnd();
364   PetscFunctionReturn(0);
365 }
366 
367 /*@C
368    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
369 
370    Collective on MatFDColoring
371 
372    Input Parameters:
373 +  coloring - the coloring context
374 -  type - either MATMFFD_WP or MATMFFD_DS
375 
376    Options Database Keys:
377 .  -mat_fd_type - "wp" or "ds"
378 
379    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
380          but the process of computing the entries is the same as as with the MatMFFD operation so we should reuse the names instead of
381          introducing another one.
382 
383    Level: intermediate
384 
385 .seealso: `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
386 
387 @*/
388 PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
389 {
390   PetscFunctionBegin;
391   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
392   /*
393      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
394      and this function is being provided as patch for a release so it shouldn't change the implementaton
395   */
396   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
397   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
398   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
399   PetscFunctionReturn(0);
400 }
401 
402 PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
403 {
404   PetscBool         flg;
405   PetscViewer       viewer;
406   PetscViewerFormat format;
407 
408   PetscFunctionBegin;
409   if (prefix) {
410     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,prefix,optionname,&viewer,&format,&flg));
411   } else {
412     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg));
413   }
414   if (flg) {
415     PetscCall(PetscViewerPushFormat(viewer,format));
416     PetscCall(MatFDColoringView(fd,viewer));
417     PetscCall(PetscViewerPopFormat(viewer));
418     PetscCall(PetscViewerDestroy(&viewer));
419   }
420   PetscFunctionReturn(0);
421 }
422 
423 /*@
424    MatFDColoringCreate - Creates a matrix coloring context for finite difference
425    computation of Jacobians.
426 
427    Collective on Mat
428 
429    Input Parameters:
430 +  mat - the matrix containing the nonzero structure of the Jacobian
431 -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
432 
433     Output Parameter:
434 .   color - the new coloring context
435 
436     Level: intermediate
437 
438 .seealso: `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
439           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
440           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
441 @*/
442 PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
443 {
444   MatFDColoring  c;
445   MPI_Comm       comm;
446   PetscInt       M,N;
447 
448   PetscFunctionBegin;
449   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
450   PetscCheck(mat->assembled,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
451   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0));
452   PetscCall(MatGetSize(mat,&M,&N));
453   PetscCheck(M == N,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
454   PetscCall(PetscObjectGetComm((PetscObject)mat,&comm));
455   PetscCall(PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView));
456 
457   c->ctype = iscoloring->ctype;
458   PetscCall(PetscObjectGetId((PetscObject)mat,&c->matid));
459 
460   if (mat->ops->fdcoloringcreate) {
461     PetscCall((*mat->ops->fdcoloringcreate)(mat,iscoloring,c));
462   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
463 
464   PetscCall(MatCreateVecs(mat,NULL,&c->w1));
465   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
466   PetscCall(VecBindToCPU(c->w1,PETSC_TRUE));
467   PetscCall(PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1));
468   PetscCall(VecDuplicate(c->w1,&c->w2));
469   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
470   PetscCall(VecBindToCPU(c->w2,PETSC_TRUE));
471   PetscCall(PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2));
472 
473   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
474   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
475   c->currentcolor = -1;
476   c->htype        = "wp";
477   c->fset         = PETSC_FALSE;
478   c->setupcalled  = PETSC_FALSE;
479 
480   *color = c;
481   PetscCall(PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c));
482   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0));
483   PetscFunctionReturn(0);
484 }
485 
486 /*@
487     MatFDColoringDestroy - Destroys a matrix coloring context that was created
488     via MatFDColoringCreate().
489 
490     Collective on MatFDColoring
491 
492     Input Parameter:
493 .   c - coloring context
494 
495     Level: intermediate
496 
497 .seealso: `MatFDColoringCreate()`
498 @*/
499 PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
500 {
501   PetscInt       i;
502   MatFDColoring  color = *c;
503 
504   PetscFunctionBegin;
505   if (!*c) PetscFunctionReturn(0);
506   if (--((PetscObject)color)->refct > 0) {*c = NULL; PetscFunctionReturn(0);}
507 
508   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
509   for (i=0; i<color->ncolors; i++) {
510     PetscCall(ISDestroy(&color->isa[i]));
511   }
512   PetscCall(PetscFree(color->isa));
513   PetscCall(PetscFree2(color->ncolumns,color->columns));
514   PetscCall(PetscFree(color->nrows));
515   if (color->htype[0] == 'w') {
516     PetscCall(PetscFree(color->matentry2));
517   } else {
518     PetscCall(PetscFree(color->matentry));
519   }
520   PetscCall(PetscFree(color->dy));
521   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
522   PetscCall(VecDestroy(&color->w1));
523   PetscCall(VecDestroy(&color->w2));
524   PetscCall(VecDestroy(&color->w3));
525   PetscCall(PetscHeaderDestroy(c));
526   PetscFunctionReturn(0);
527 }
528 
529 /*@C
530     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
531       that are currently being perturbed.
532 
533     Not Collective
534 
535     Input Parameter:
536 .   coloring - coloring context created with MatFDColoringCreate()
537 
538     Output Parameters:
539 +   n - the number of local columns being perturbed
540 -   cols - the column indices, in global numbering
541 
542    Level: advanced
543 
544    Note: IF the matrix type is BAIJ, then the block column indices are returned
545 
546    Fortran Note:
547    This routine has a different interface for Fortran
548 $     #include <petsc/finclude/petscmat.h>
549 $          use petscmat
550 $          PetscInt, pointer :: array(:)
551 $          PetscErrorCode  ierr
552 $          MatFDColoring   i
553 $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
554 $      use the entries of array ...
555 $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
556 
557 .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
558 
559 @*/
560 PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[])
561 {
562   PetscFunctionBegin;
563   if (coloring->currentcolor >= 0) {
564     *n    = coloring->ncolumns[coloring->currentcolor];
565     *cols = coloring->columns[coloring->currentcolor];
566   } else {
567     *n = 0;
568   }
569   PetscFunctionReturn(0);
570 }
571 
572 /*@
573     MatFDColoringApply - Given a matrix for which a MatFDColoring context
574     has been created, computes the Jacobian for a function via finite differences.
575 
576     Collective on MatFDColoring
577 
578     Input Parameters:
579 +   mat - location to store Jacobian
580 .   coloring - coloring context created with MatFDColoringCreate()
581 .   x1 - location at which Jacobian is to be computed
582 -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
583 
584     Options Database Keys:
585 +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
586 .    -mat_fd_coloring_view - Activates basic viewing or coloring
587 .    -mat_fd_coloring_view draw - Activates drawing of coloring
588 -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
589 
590     Level: intermediate
591 
592 .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
593 
594 @*/
595 PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
596 {
597   PetscBool      eq;
598 
599   PetscFunctionBegin;
600   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
601   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
602   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
603   PetscCall(PetscObjectCompareId((PetscObject)J,coloring->matid,&eq));
604   PetscCheck(eq,PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
605   PetscCheck(coloring->f,PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
606   PetscCheck(J->ops->fdcoloringapply,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
607   PetscCheck(coloring->setupcalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
608 
609   PetscCall(MatSetUnfactored(J));
610   PetscCall(PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0));
611   PetscCall((*J->ops->fdcoloringapply)(J,coloring,x1,sctx));
612   PetscCall(PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0));
613   if (!coloring->viewed) {
614     PetscCall(MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view"));
615     coloring->viewed = PETSC_TRUE;
616   }
617   PetscFunctionReturn(0);
618 }
619