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