xref: /petsc/src/mat/matfd/fdmatrix.c (revision 634064b45b5c838063ae82f97ffb7e99245dcdb5)
1 /*
2    This is where the abstract matrix operations are defined that are
3   used for finite difference computations of Jacobians using coloring.
4 */
5 
6 #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
7 
8 /* Logging support */
9 int MAT_FDCOLORING_COOKIE = 0;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatFDColoringSetF"
13 int MatFDColoringSetF(MatFDColoring fd,Vec F)
14 {
15   PetscFunctionBegin;
16   fd->F = F;
17   PetscFunctionReturn(0);
18 }
19 
20 #undef __FUNCT__
21 #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
22 static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
23 {
24   MatFDColoring fd = (MatFDColoring)Aa;
25   int           ierr,i,j;
26   PetscReal     x,y;
27 
28   PetscFunctionBegin;
29 
30   /* loop over colors  */
31   for (i=0; i<fd->ncolors; i++) {
32     for (j=0; j<fd->nrows[i]; j++) {
33       y = fd->M - fd->rows[i][j] - fd->rstart;
34       x = fd->columnsforrow[i][j];
35       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
36     }
37   }
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "MatFDColoringView_Draw"
43 static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
44 {
45   int         ierr;
46   PetscTruth  isnull;
47   PetscDraw   draw;
48   PetscReal   xr,yr,xl,yl,h,w;
49 
50   PetscFunctionBegin;
51   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
52   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
53 
54   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
55 
56   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
57   xr += w;     yr += h;    xl = -w;     yl = -h;
58   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
59   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
60   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "MatFDColoringView"
66 /*@C
67    MatFDColoringView - Views a finite difference coloring context.
68 
69    Collective on MatFDColoring
70 
71    Input  Parameters:
72 +  c - the coloring context
73 -  viewer - visualization context
74 
75    Level: intermediate
76 
77    Notes:
78    The available visualization contexts include
79 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
80 .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
81         output where only the first processor opens
82         the file.  All other processors send their
83         data to the first processor to print.
84 -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
85 
86    Notes:
87      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
88    involves more than 33 then some seemingly identical colors are displayed making it look
89    like an illegal coloring. This is just a graphical artifact.
90 
91 .seealso: MatFDColoringCreate()
92 
93 .keywords: Mat, finite differences, coloring, view
94 @*/
95 int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
96 {
97   int               i,j,ierr;
98   PetscTruth        isdraw,iascii;
99   PetscViewerFormat format;
100 
101   PetscFunctionBegin;
102   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
103   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
104   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
105   PetscCheckSameComm(c,1,viewer,2);
106 
107   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
108   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
109   if (isdraw) {
110     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
111   } else if (iascii) {
112     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
113     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
114     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
115     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
116 
117     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
118     if (format != PETSC_VIEWER_ASCII_INFO) {
119       for (i=0; i<c->ncolors; i++) {
120         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
121         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
122         for (j=0; j<c->ncolumns[i]; j++) {
123           ierr = PetscViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
124         }
125         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
126         for (j=0; j<c->nrows[i]; j++) {
127           ierr = PetscViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
128         }
129       }
130     }
131     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
132   } else {
133     SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "MatFDColoringSetParameters"
140 /*@
141    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
142    a Jacobian matrix using finite differences.
143 
144    Collective on MatFDColoring
145 
146    The Jacobian is estimated with the differencing approximation
147 .vb
148        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
149        h = error_rel*u[i]                 if  abs(u[i]) > umin
150          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
151        dx_{i} = (0, ... 1, .... 0)
152 .ve
153 
154    Input Parameters:
155 +  coloring - the coloring context
156 .  error_rel - relative error
157 -  umin - minimum allowable u-value magnitude
158 
159    Level: advanced
160 
161 .keywords: Mat, finite differences, coloring, set, parameters
162 
163 .seealso: MatFDColoringCreate()
164 @*/
165 int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
166 {
167   PetscFunctionBegin;
168   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
169 
170   if (error != PETSC_DEFAULT) matfd->error_rel = error;
171   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "MatFDColoringSetFrequency"
177 /*@
178    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
179    matrices.
180 
181    Collective on MatFDColoring
182 
183    Input Parameters:
184 +  coloring - the coloring context
185 -  freq - frequency (default is 1)
186 
187    Options Database Keys:
188 .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
189 
190    Level: advanced
191 
192    Notes:
193    Using a modified Newton strategy, where the Jacobian remains fixed for several
194    iterations, can be cost effective in terms of overall nonlinear solution
195    efficiency.  This parameter indicates that a new Jacobian will be computed every
196    <freq> nonlinear iterations.
197 
198 .keywords: Mat, finite differences, coloring, set, frequency
199 
200 .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
201 @*/
202 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
203 {
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
206 
207   matfd->freq = freq;
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatFDColoringGetFrequency"
213 /*@
214    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
215    matrices.
216 
217    Not Collective
218 
219    Input Parameters:
220 .  coloring - the coloring context
221 
222    Output Parameters:
223 .  freq - frequency (default is 1)
224 
225    Options Database Keys:
226 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
227 
228    Level: advanced
229 
230    Notes:
231    Using a modified Newton strategy, where the Jacobian remains fixed for several
232    iterations, can be cost effective in terms of overall nonlinear solution
233    efficiency.  This parameter indicates that a new Jacobian will be computed every
234    <freq> nonlinear iterations.
235 
236 .keywords: Mat, finite differences, coloring, get, frequency
237 
238 .seealso: MatFDColoringSetFrequency()
239 @*/
240 int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
241 {
242   PetscFunctionBegin;
243   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
244   *freq = matfd->freq;
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "MatFDColoringSetFunction"
250 /*@C
251    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
252 
253    Collective on MatFDColoring
254 
255    Input Parameters:
256 +  coloring - the coloring context
257 .  f - the function
258 -  fctx - the optional user-defined function context
259 
260    Level: intermediate
261 
262    Notes:
263     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
264   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
265   with the TS solvers.
266 
267 .keywords: Mat, Jacobian, finite differences, set, function
268 @*/
269 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
270 {
271   PetscFunctionBegin;
272   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
273   matfd->f    = f;
274   matfd->fctx = fctx;
275   PetscFunctionReturn(0);
276 }
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "MatFDColoringSetFromOptions"
280 /*@
281    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
282    the options database.
283 
284    Collective on MatFDColoring
285 
286    The Jacobian, F'(u), is estimated with the differencing approximation
287 .vb
288        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
289        h = error_rel*u[i]                 if  abs(u[i]) > umin
290          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
291        dx_{i} = (0, ... 1, .... 0)
292 .ve
293 
294    Input Parameter:
295 .  coloring - the coloring context
296 
297    Options Database Keys:
298 +  -mat_fd_coloring_err <err> - Sets <err> (square root
299            of relative error in the function)
300 .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
301 .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
302 .  -mat_fd_coloring_view - Activates basic viewing
303 .  -mat_fd_coloring_view_info - Activates viewing info
304 -  -mat_fd_coloring_view_draw - Activates drawing
305 
306     Level: intermediate
307 
308 .keywords: Mat, finite differences, parameters
309 
310 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
311 
312 @*/
313 int MatFDColoringSetFromOptions(MatFDColoring matfd)
314 {
315   int        ierr;
316 
317   PetscFunctionBegin;
318   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
319 
320   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
321     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
322     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
323     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
324     /* not used here; but so they are presented in the GUI */
325     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
326     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
327     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
328   PetscOptionsEnd();CHKERRQ(ierr);
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "MatFDColoringView_Private"
334 int MatFDColoringView_Private(MatFDColoring fd)
335 {
336   int        ierr;
337   PetscTruth flg;
338 
339   PetscFunctionBegin;
340   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
341   if (flg) {
342     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
343   }
344   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
345   if (flg) {
346     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
347     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
348     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
349   }
350   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
351   if (flg) {
352     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
353     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
354   }
355   PetscFunctionReturn(0);
356 }
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "MatFDColoringCreate"
360 /*@C
361    MatFDColoringCreate - Creates a matrix coloring context for finite difference
362    computation of Jacobians.
363 
364    Collective on Mat
365 
366    Input Parameters:
367 +  mat - the matrix containing the nonzero structure of the Jacobian
368 -  iscoloring - the coloring of the matrix
369 
370     Output Parameter:
371 .   color - the new coloring context
372 
373     Level: intermediate
374 
375 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
376           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
377           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
378           MatFDColoringSetParameters()
379 @*/
380 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
381 {
382   MatFDColoring c;
383   MPI_Comm      comm;
384   int           ierr,M,N;
385 
386   PetscFunctionBegin;
387   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
388   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
389   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
390 
391   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
392   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
393   PetscLogObjectCreate(c);
394 
395   if (mat->ops->fdcoloringcreate) {
396     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
397   } else {
398     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
399   }
400 
401   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
402   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
403   c->freq              = 1;
404   c->usersetsrecompute = PETSC_FALSE;
405   c->recompute         = PETSC_FALSE;
406   c->currentcolor      = -1;
407 
408   *color = c;
409   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "MatFDColoringDestroy"
415 /*@C
416     MatFDColoringDestroy - Destroys a matrix coloring context that was created
417     via MatFDColoringCreate().
418 
419     Collective on MatFDColoring
420 
421     Input Parameter:
422 .   c - coloring context
423 
424     Level: intermediate
425 
426 .seealso: MatFDColoringCreate()
427 @*/
428 int MatFDColoringDestroy(MatFDColoring c)
429 {
430   int i,ierr;
431 
432   PetscFunctionBegin;
433   if (--c->refct > 0) PetscFunctionReturn(0);
434 
435   for (i=0; i<c->ncolors; i++) {
436     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
437     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
438     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
439     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
440   }
441   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
442   ierr = PetscFree(c->columns);CHKERRQ(ierr);
443   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
444   ierr = PetscFree(c->rows);CHKERRQ(ierr);
445   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
446   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
447   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
448   if (c->w1) {
449     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
450     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
451     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
452   }
453   PetscLogObjectDestroy(c);
454   PetscHeaderDestroy(c);
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
460 /*@C
461     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
462       that are currently being perturbed.
463 
464     Not Collective
465 
466     Input Parameters:
467 .   coloring - coloring context created with MatFDColoringCreate()
468 
469     Output Parameters:
470 +   n - the number of local columns being perturbed
471 -   cols - the column indices, in global numbering
472 
473    Level: intermediate
474 
475 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
476 
477 .keywords: coloring, Jacobian, finite differences
478 @*/
479 EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int *cols[])
480 {
481   PetscFunctionBegin;
482   if (coloring->currentcolor >= 0) {
483     *n    = coloring->ncolumns[coloring->currentcolor];
484     *cols = coloring->columns[coloring->currentcolor];
485   } else {
486     *n = 0;
487   }
488   PetscFunctionReturn(0);
489 }
490 
491 #undef __FUNCT__
492 #define __FUNCT__ "MatFDColoringApply"
493 /*@
494     MatFDColoringApply - Given a matrix for which a MatFDColoring context
495     has been created, computes the Jacobian for a function via finite differences.
496 
497     Collective on MatFDColoring
498 
499     Input Parameters:
500 +   mat - location to store Jacobian
501 .   coloring - coloring context created with MatFDColoringCreate()
502 .   x1 - location at which Jacobian is to be computed
503 -   sctx - optional context required by function (actually a SNES context)
504 
505     Options Database Keys:
506 +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
507 .    -mat_fd_coloring_view - Activates basic viewing or coloring
508 .    -mat_fd_coloring_view_draw - Activates drawing of coloring
509 -    -mat_fd_coloring_view_info - Activates viewing of coloring info
510 
511     Level: intermediate
512 
513 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
514 
515 .keywords: coloring, Jacobian, finite differences
516 @*/
517 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
518 {
519   int           (*f)(void*,Vec,Vec,void*) = (int (*)(void*,Vec,Vec,void *))coloring->f;
520   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
521   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
522   PetscScalar   *vscale_array;
523   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
524   Vec           w1,w2,w3;
525   void          *fctx = coloring->fctx;
526   PetscTruth    flg;
527 
528 
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
531   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
532   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
533 
534   if (coloring->usersetsrecompute) {
535     if (!coloring->recompute) {
536       *flag = SAME_PRECONDITIONER;
537       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
538       PetscFunctionReturn(0);
539     } else {
540       coloring->recompute = PETSC_FALSE;
541     }
542   }
543 
544   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
545   if (J->ops->fdcoloringapply) {
546     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
547   } else {
548 
549     if (!coloring->w1) {
550       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
551       PetscLogObjectParent(coloring,coloring->w1);
552       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
553       PetscLogObjectParent(coloring,coloring->w2);
554       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
555       PetscLogObjectParent(coloring,coloring->w3);
556     }
557     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
558 
559     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
560     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
561     if (flg) {
562       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
563     } else {
564       ierr = MatZeroEntries(J);CHKERRQ(ierr);
565     }
566 
567     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
568     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
569 
570     /*
571       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
572       coloring->F for the coarser grids from the finest
573     */
574     if (coloring->F) {
575       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
576       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
577       if (m1 != m2) {
578 	coloring->F = 0;
579       }
580     }
581 
582     if (coloring->F) {
583       w1          = coloring->F; /* use already computed value of function */
584       coloring->F = 0;
585     } else {
586       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
587       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
588       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
589     }
590 
591     /*
592        Compute all the scale factors and share with other processors
593     */
594     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
595     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
596     for (k=0; k<coloring->ncolors; k++) {
597       /*
598 	Loop over each column associated with color adding the
599 	perturbation to the vector w3.
600       */
601       for (l=0; l<coloring->ncolumns[k]; l++) {
602 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
603 	dx  = xx[col];
604 	if (dx == 0.0) dx = 1.0;
605 #if !defined(PETSC_USE_COMPLEX)
606 	if (dx < umin && dx >= 0.0)      dx = umin;
607 	else if (dx < 0.0 && dx > -umin) dx = -umin;
608 #else
609 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
610 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
611 #endif
612 	dx                *= epsilon;
613 	vscale_array[col] = 1.0/dx;
614       }
615     }
616     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
617     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
618     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619 
620     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
621 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
622 
623     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
624     else                        vscaleforrow = coloring->columnsforrow;
625 
626     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
627     /*
628       Loop over each color
629     */
630     for (k=0; k<coloring->ncolors; k++) {
631       coloring->currentcolor = k;
632       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
633       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
634       /*
635 	Loop over each column associated with color adding the
636 	perturbation to the vector w3.
637       */
638       for (l=0; l<coloring->ncolumns[k]; l++) {
639 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
640 	dx  = xx[col];
641 	if (dx == 0.0) dx = 1.0;
642 #if !defined(PETSC_USE_COMPLEX)
643 	if (dx < umin && dx >= 0.0)      dx = umin;
644 	else if (dx < 0.0 && dx > -umin) dx = -umin;
645 #else
646 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
647 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
648 #endif
649 	dx            *= epsilon;
650 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
651 	w3_array[col] += dx;
652       }
653       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
654 
655       /*
656 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
657       */
658 
659       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
660       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
661       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
662       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
663 
664       /*
665 	Loop over rows of vector, putting results into Jacobian matrix
666       */
667       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
668       for (l=0; l<coloring->nrows[k]; l++) {
669 	row    = coloring->rows[k][l];
670 	col    = coloring->columnsforrow[k][l];
671 	y[row] *= vscale_array[vscaleforrow[k][l]];
672 	srow   = row + start;
673 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
674       }
675       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
676     }
677     coloring->currentcolor = -1;
678     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
679     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
680     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
681     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
682   }
683   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
684 
685   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
686   if (flg) {
687     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
688   }
689   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
690 
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "MatFDColoringApplyTS"
696 /*@
697     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
698     has been created, computes the Jacobian for a function via finite differences.
699 
700    Collective on Mat, MatFDColoring, and Vec
701 
702     Input Parameters:
703 +   mat - location to store Jacobian
704 .   coloring - coloring context created with MatFDColoringCreate()
705 .   x1 - location at which Jacobian is to be computed
706 -   sctx - optional context required by function (actually a SNES context)
707 
708    Options Database Keys:
709 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
710 
711    Level: intermediate
712 
713 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
714 
715 .keywords: coloring, Jacobian, finite differences
716 @*/
717 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
718 {
719   int           (*f)(void*,PetscReal,Vec,Vec,void*)=(int (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
720   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
721   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
722   PetscScalar   *vscale_array;
723   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
724   Vec           w1,w2,w3;
725   void          *fctx = coloring->fctx;
726   PetscTruth    flg;
727 
728   PetscFunctionBegin;
729   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
730   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
731   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
732 
733   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
734   if (!coloring->w1) {
735     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
736     PetscLogObjectParent(coloring,coloring->w1);
737     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
738     PetscLogObjectParent(coloring,coloring->w2);
739     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
740     PetscLogObjectParent(coloring,coloring->w3);
741   }
742   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
743 
744   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
745   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
746   if (flg) {
747     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
748   } else {
749     ierr = MatZeroEntries(J);CHKERRQ(ierr);
750   }
751 
752   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
753   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
754   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
755   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
756   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
757 
758   /*
759       Compute all the scale factors and share with other processors
760   */
761   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
762   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
763   for (k=0; k<coloring->ncolors; k++) {
764     /*
765        Loop over each column associated with color adding the
766        perturbation to the vector w3.
767     */
768     for (l=0; l<coloring->ncolumns[k]; l++) {
769       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
770       dx  = xx[col];
771       if (dx == 0.0) dx = 1.0;
772 #if !defined(PETSC_USE_COMPLEX)
773       if (dx < umin && dx >= 0.0)      dx = umin;
774       else if (dx < 0.0 && dx > -umin) dx = -umin;
775 #else
776       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
777       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
778 #endif
779       dx                *= epsilon;
780       vscale_array[col] = 1.0/dx;
781     }
782   }
783   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
784   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
785   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
786 
787   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
788   else                        vscaleforrow = coloring->columnsforrow;
789 
790   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
791   /*
792       Loop over each color
793   */
794   for (k=0; k<coloring->ncolors; k++) {
795     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
796     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
797     /*
798        Loop over each column associated with color adding the
799        perturbation to the vector w3.
800     */
801     for (l=0; l<coloring->ncolumns[k]; l++) {
802       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
803       dx  = xx[col];
804       if (dx == 0.0) dx = 1.0;
805 #if !defined(PETSC_USE_COMPLEX)
806       if (dx < umin && dx >= 0.0)      dx = umin;
807       else if (dx < 0.0 && dx > -umin) dx = -umin;
808 #else
809       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
810       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
811 #endif
812       dx            *= epsilon;
813       w3_array[col] += dx;
814     }
815     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
816 
817     /*
818        Evaluate function at x1 + dx (here dx is a vector of perturbations)
819     */
820     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
821     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
822     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
823     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
824 
825     /*
826        Loop over rows of vector, putting results into Jacobian matrix
827     */
828     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
829     for (l=0; l<coloring->nrows[k]; l++) {
830       row    = coloring->rows[k][l];
831       col    = coloring->columnsforrow[k][l];
832       y[row] *= vscale_array[vscaleforrow[k][l]];
833       srow   = row + start;
834       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
835     }
836     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
837   }
838   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
839   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
840   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
841   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
842   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
843   PetscFunctionReturn(0);
844 }
845 
846 
847 #undef __FUNCT__
848 #define __FUNCT__ "MatFDColoringSetRecompute()"
849 /*@C
850    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
851      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
852      no additional Jacobian's will be computed (the same one will be used) until this is
853      called again.
854 
855    Collective on MatFDColoring
856 
857    Input  Parameters:
858 .  c - the coloring context
859 
860    Level: intermediate
861 
862    Notes: The MatFDColoringSetFrequency() is ignored once this is called
863 
864 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
865 
866 .keywords: Mat, finite differences, coloring
867 @*/
868 int MatFDColoringSetRecompute(MatFDColoring c)
869 {
870   PetscFunctionBegin;
871   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
872   c->usersetsrecompute = PETSC_TRUE;
873   c->recompute         = PETSC_TRUE;
874   PetscFunctionReturn(0);
875 }
876 
877 
878