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