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