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