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