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