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