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