xref: /petsc/src/mat/matfd/fdmatrix.c (revision 5cd905551efc76d042fcd2cc746dd1c50567705b)
1 
2 #ifdef PETSC_RCS_HEADER
3 static char vcid[] = "$Id: fdmatrix.c,v 1.30 1998/04/08 16:06:13 curfman Exp bsmith $";
4 #endif
5 
6 /*
7    This is where the abstract matrix operations are defined that are
8   used for finite difference computations of Jacobians using coloring.
9 */
10 
11 #include "petsc.h"
12 #include "src/mat/matimpl.h"        /*I "mat.h" I*/
13 #include "src/vec/vecimpl.h"
14 #include "pinclude/pviewer.h"
15 
16 #undef __FUNC__
17 #define __FUNC__ "MatFDColoringView_Draw"
18 static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
19 {
20   int         ierr,i,j,pause;
21   PetscTruth  isnull;
22   Draw        draw;
23   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
24   DrawButton  button;
25 
26   PetscFunctionBegin;
27   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
28   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
29   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
30 
31   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
32   xr += w;    yr += h;  xl = -w;     yl = -h;
33   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
34 
35   /* loop over colors  */
36   for (i=0; i<fd->ncolors; i++ ) {
37     for ( j=0; j<fd->nrows[i]; j++ ) {
38       y = fd->M - fd->rows[i][j] - fd->rstart;
39       x = fd->columnsforrow[i][j];
40       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); CHKERRQ(ierr);
41     }
42   }
43   ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr);
44   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
45   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
46   ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
47   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);  CHKERRQ(ierr);
48   while (button != BUTTON_RIGHT) {
49     ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
50     if (button == BUTTON_LEFT) scale = .5;
51     else if (button == BUTTON_CENTER) scale = 2.;
52     xl = scale*(xl + w - xc) + xc - w*scale;
53     xr = scale*(xr - w - xc) + xc + w*scale;
54     yl = scale*(yl + h - yc) + yc - h*scale;
55     yr = scale*(yr - h - yc) + yc + h*scale;
56     w *= scale; h *= scale;
57     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
58     /* loop over colors  */
59     for (i=0; i<fd->ncolors; i++ ) {
60       for ( j=0; j<fd->nrows[i]; j++ ) {
61         y = fd->M - fd->rows[i][j] - fd->rstart;
62         x = fd->columnsforrow[i][j];
63         ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); CHKERRQ(ierr);
64       }
65     }
66     ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
67     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);  CHKERRQ(ierr);
68   }
69 
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNC__
74 #define __FUNC__ "MatFDColoringView"
75 /*@C
76    MatFDColoringView - Views a finite difference coloring context.
77 
78    Input  Parameters:
79 .  c - the coloring context
80 .  viewer - visualization context
81 
82    Notes:
83    The available visualization contexts include
84 $     VIEWER_STDOUT_SELF - standard output (default)
85 $     VIEWER_STDOUT_WORLD - synchronized standard
86 $       output where only the first processor opens
87 $       the file.  All other processors send their
88 $       data to the first processor to print.
89 $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
90 
91 .seealso: MatFDColoringCreate()
92 
93 .keywords: Mat, finite differences, coloring, view
94 @*/
95 int MatFDColoringView(MatFDColoring c,Viewer viewer)
96 {
97   ViewerType vtype;
98   int        i,j,format,ierr;
99   FILE       *fd;
100 
101   PetscFunctionBegin;
102   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
103   if (viewer) {PetscValidHeader(viewer);}
104   else {viewer = VIEWER_STDOUT_SELF;}
105 
106   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
107   if (vtype == DRAW_VIEWER) {
108     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
109     PetscFunctionReturn(0);
110   } else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
111     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
112     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
113     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
114     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
115     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
116 
117     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
118     if (format != VIEWER_FORMAT_ASCII_INFO) {
119       for ( i=0; i<c->ncolors; i++ ) {
120         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
121         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
122         for ( j=0; j<c->ncolumns[i]; j++ ) {
123           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
124         }
125         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
126         for ( j=0; j<c->nrows[i]; j++ ) {
127           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
128         }
129       }
130     }
131   } else {
132     SETERRQ(1,1,"Viewer type not supported for this object");
133   }
134   PetscFunctionReturn(0);
135 }
136 
137 #undef __FUNC__
138 #define __FUNC__ "MatFDColoringSetParameters"
139 /*@
140    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
141    a Jacobian matrix using finite differences.
142 
143 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
144 $        h = error_rel*u[i]                    if  u[i] > umin
145 $          = error_rel*umin                    else
146 $
147 $   dx_{i} = (0, ... 1, .... 0)
148 
149    Input Parameters:
150 .  coloring - the coloring context
151 .  error_rel - relative error
152 .  umin - minimum allowable u-value
153 
154 .keywords: Mat, finite differences, coloring, set, parameters
155 
156 .seealso: MatFDColoringCreate()
157 @*/
158 int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
159 {
160   PetscFunctionBegin;
161   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
162 
163   if (error != PETSC_DEFAULT) matfd->error_rel = error;
164   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNC__
169 #define __FUNC__ "MatFDColoringSetFrequency"
170 /*@
171    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
172    matrices.
173 
174    Input Parameters:
175 .  coloring - the coloring context
176 .  freq - frequency (default is 1)
177 
178    Notes:
179    Using a modified Newton strategy, where the Jacobian remains fixed for several
180    iterations, can be cost effective in terms of overall nonlinear solution
181    efficiency.  This parameter indicates that a new Jacobian will be computed every
182    <freq> nonlinear iterations.
183 
184    Options Database Keys:
185 $  -mat_fd_coloring_freq <freq>
186 
187 .keywords: Mat, finite differences, coloring, set, frequency
188 @*/
189 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
190 {
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
193 
194   matfd->freq = freq;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNC__
199 #define __FUNC__ "MatFDColoringGetFrequency"
200 /*@
201    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
202    matrices.
203 
204    Input Parameters:
205 .  coloring - the coloring context
206 
207    Output Parameters:
208 .  freq - frequency (default is 1)
209 
210    Notes:
211    Using a modified Newton strategy, where the Jacobian remains fixed for several
212    iterations, can be cost effective in terms of overall nonlinear solution
213    efficiency.  This parameter indicates that a new Jacobian will be computed every
214    <freq> nonlinear iterations.
215 
216    Options Database Keys:
217 $  -mat_fd_coloring_freq <freq>
218 
219 .keywords: Mat, finite differences, coloring, get, frequency
220 @*/
221 int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
222 {
223   PetscFunctionBegin;
224   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
225 
226   *freq = matfd->freq;
227   PetscFunctionReturn(0);
228 }
229 
230 #undef __FUNC__
231 #define __FUNC__ "MatFDColoringSetFunction"
232 /*@C
233    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
234 
235    Input Parameters:
236 .  coloring - the coloring context
237 .  f - the function
238 .  fctx - the optional user-defined function context
239 
240 .keywords: Mat, Jacobian, finite differences, set, function
241 @*/
242 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
243 {
244   PetscFunctionBegin;
245   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
246 
247   matfd->f    = f;
248   matfd->fctx = fctx;
249 
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNC__
254 #define __FUNC__ "MatFDColoringSetFromOptions"
255 /*@
256    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
257    the options database.
258 
259    The Jacobian is estimated with the differencing approximation
260 $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
261 $        h = error_rel*u[i]                    if  u[i] > umin
262 $          = error_rel*umin                      else
263 $
264 $   dx_{i} = (0, ... 1, .... 0)
265 
266    Input Parameters:
267 .  coloring - the coloring context
268 
269    Options Database Keys:
270 $  -mat_fd_coloring_error <err>, where <err> is the square root
271 $           of relative error in the function
272 $  -mat_fd_coloring_umin  <umin>, where umin is described above
273 $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
274 $           computing a new Jacobian
275 $  -mat_fd_coloring_view
276 $  -mat_fd_coloring_view_info
277 $  -mat_fd_coloring_view_draw
278 
279 .keywords: Mat, finite differences, parameters
280 @*/
281 int MatFDColoringSetFromOptions(MatFDColoring matfd)
282 {
283   int    ierr,flag,freq = 1;
284   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
285 
286   PetscFunctionBegin;
287   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
288 
289   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
290   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
291   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
292   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
293   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
294   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
295   if (flag) {
296     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
297   }
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNC__
302 #define __FUNC__ "MatFDColoringPrintHelp"
303 /*@
304     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
305     using coloring.
306 
307     Input Parameter:
308 .   fdcoloring - the MatFDColoring context
309 
310 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
311 @*/
312 int MatFDColoringPrintHelp(MatFDColoring fd)
313 {
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
316 
317   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
318   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
319   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
320   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");
321   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");
322   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");
323   PetscFunctionReturn(0);
324 }
325 
326 int MatFDColoringView_Private(MatFDColoring fd)
327 {
328   int ierr,flg;
329 
330   PetscFunctionBegin;
331   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
332   if (flg) {
333     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
334   }
335   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
336   if (flg) {
337     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
338     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
339     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
340   }
341   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
342   if (flg) {
343     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
344     ierr = ViewerFlush(VIEWER_DRAWX_(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    Input Parameters:
356 .  mat - the matrix containing the nonzero structure of the Jacobian
357 .  iscoloring - the coloring of the matrix
358 
359     Output Parameter:
360 .   color - the new coloring context
361 
362     Options Database Keys:
363 $    -mat_fd_coloring_view
364 $    -mat_fd_coloring_view_draw
365 $    -mat_fd_coloring_view_info
366 
367 .seealso: MatFDColoringDestroy()
368 @*/
369 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
370 {
371   MatFDColoring c;
372   MPI_Comm      comm;
373   int           ierr,M,N;
374 
375   PetscFunctionBegin;
376   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
377   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
378 
379   PetscObjectGetComm((PetscObject)mat,&comm);
380   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
381   PLogObjectCreate(c);
382 
383   if (mat->ops->fdcoloringcreate) {
384     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
385   } else {
386     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
387   }
388 
389   c->error_rel = 1.e-8;
390   c->umin      = 1.e-6;
391   c->freq      = 1;
392 
393   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
394 
395   *color = c;
396 
397   PetscFunctionReturn(0);
398 }
399 
400 #undef __FUNC__
401 #define __FUNC__ "MatFDColoringDestroy"
402 /*@C
403     MatFDColoringDestroy - Destroys a matrix coloring context that was created
404     via MatFDColoringCreate().
405 
406     Input Parameter:
407 .   c - coloring context
408 
409 .seealso: MatFDColoringCreate()
410 @*/
411 int MatFDColoringDestroy(MatFDColoring c)
412 {
413   int i,ierr;
414 
415   PetscFunctionBegin;
416   if (--c->refct > 0) PetscFunctionReturn(0);
417 
418 
419   for ( i=0; i<c->ncolors; i++ ) {
420     if (c->columns[i])       PetscFree(c->columns[i]);
421     if (c->rows[i])          PetscFree(c->rows[i]);
422     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
423   }
424   PetscFree(c->ncolumns);
425   PetscFree(c->columns);
426   PetscFree(c->nrows);
427   PetscFree(c->rows);
428   PetscFree(c->columnsforrow);
429   PetscFree(c->scale);
430   if (c->w1) {
431     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
432     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
433     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
434   }
435   PLogObjectDestroy(c);
436   PetscHeaderDestroy(c);
437   PetscFunctionReturn(0);
438 }
439 
440 #include "snes.h"
441 
442 #undef __FUNC__
443 #define __FUNC__ "MatFDColoringApply"
444 /*@
445     MatFDColoringApply - Given a matrix for which a MatFDColoring context
446     has been created, computes the Jacobian for a function via finite differences.
447 
448     Input Parameters:
449 .   mat - location to store Jacobian
450 .   coloring - coloring context created with MatFDColoringCreate()
451 .   x1 - location at which Jacobian is to be computed
452 .   sctx - optional context required by function (actually a SNES context)
453 
454    Options Database Keys:
455 $  -mat_fd_coloring_freq <freq>
456 
457 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
458 
459 .keywords: coloring, Jacobian, finite differences
460 @*/
461 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
462 {
463   int           k,fg,ierr,N,start,end,l,row,col,srow;
464   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
465   double        epsilon = coloring->error_rel, umin = coloring->umin;
466   MPI_Comm      comm = coloring->comm;
467   Vec           w1,w2,w3;
468   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
469   void          *fctx = coloring->fctx;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(J,MAT_COOKIE);
473   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
474   PetscValidHeaderSpecific(x1,VEC_COOKIE);
475 
476 
477   if (!coloring->w1) {
478     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
479     PLogObjectParent(coloring,coloring->w1);
480     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
481     PLogObjectParent(coloring,coloring->w2);
482     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
483     PLogObjectParent(coloring,coloring->w3);
484   }
485   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
486 
487   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
488   if (fg) {
489     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
490   } else {
491     ierr = MatZeroEntries(J); CHKERRQ(ierr);
492   }
493 
494   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
495   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
496   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
497   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
498 
499   PetscMemzero(wscale,N*sizeof(Scalar));
500   /*
501       Loop over each color
502   */
503 
504   for (k=0; k<coloring->ncolors; k++) {
505     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
506     /*
507        Loop over each column associated with color adding the
508        perturbation to the vector w3.
509     */
510     for (l=0; l<coloring->ncolumns[k]; l++) {
511       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
512       dx  = xx[col-start];
513       if (dx == 0.0) dx = 1.0;
514 #if !defined(USE_PETSC_COMPLEX)
515       if (dx < umin && dx >= 0.0)      dx = umin;
516       else if (dx < 0.0 && dx > -umin) dx = -umin;
517 #else
518       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
519       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
520 #endif
521       dx          *= epsilon;
522       wscale[col] = 1.0/dx;
523       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
524     }
525     VecRestoreArray(x1,&xx);
526     /*
527        Evaluate function at x1 + dx (here dx is a vector of perturbations)
528     */
529     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
530     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
531     /* Communicate scale to all processors */
532 #if !defined(USE_PETSC_COMPLEX)
533     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
534 #else
535     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
536 #endif
537     /*
538        Loop over rows of vector, putting results into Jacobian matrix
539     */
540     VecGetArray(w2,&y);
541     for (l=0; l<coloring->nrows[k]; l++) {
542       row    = coloring->rows[k][l];
543       col    = coloring->columnsforrow[k][l];
544       y[row] *= scale[col];
545       srow   = row + start;
546       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
547     }
548     VecRestoreArray(w2,&y);
549   }
550   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
551   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
552   PetscFunctionReturn(0);
553 }
554 
555 #include "ts.h"
556 
557 #undef __FUNC__
558 #define __FUNC__ "MatFDColoringApplyTS"
559 /*@
560     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
561     has been created, computes the Jacobian for a function via finite differences.
562 
563     Input Parameters:
564 .   mat - location to store Jacobian
565 .   coloring - coloring context created with MatFDColoringCreate()
566 .   x1 - location at which Jacobian is to be computed
567 .   sctx - optional context required by function (actually a SNES context)
568 
569    Options Database Keys:
570 $  -mat_fd_coloring_freq <freq>
571 
572 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
573 
574 .keywords: coloring, Jacobian, finite differences
575 @*/
576 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
577 {
578   int           k,fg,ierr,N,start,end,l,row,col,srow;
579   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
580   double        epsilon = coloring->error_rel, umin = coloring->umin;
581   MPI_Comm      comm = coloring->comm;
582   Vec           w1,w2,w3;
583   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
584   void          *fctx = coloring->fctx;
585 
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(J,MAT_COOKIE);
588   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
589   PetscValidHeaderSpecific(x1,VEC_COOKIE);
590 
591 
592   if (!coloring->w1) {
593     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
594     PLogObjectParent(coloring,coloring->w1);
595     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
596     PLogObjectParent(coloring,coloring->w2);
597     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
598     PLogObjectParent(coloring,coloring->w3);
599   }
600   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
601 
602   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
603   if (fg) {
604     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
605   } else {
606     ierr = MatZeroEntries(J); CHKERRQ(ierr);
607   }
608 
609   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
610   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
611   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
612   ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr);
613 
614   PetscMemzero(wscale,N*sizeof(Scalar));
615   /*
616       Loop over each color
617   */
618 
619   for (k=0; k<coloring->ncolors; k++) {
620     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
621     /*
622        Loop over each column associated with color adding the
623        perturbation to the vector w3.
624     */
625     for (l=0; l<coloring->ncolumns[k]; l++) {
626       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
627       dx  = xx[col-start];
628       if (dx == 0.0) dx = 1.0;
629 #if !defined(USE_PETSC_COMPLEX)
630       if (dx < umin && dx >= 0.0)      dx = umin;
631       else if (dx < 0.0 && dx > -umin) dx = -umin;
632 #else
633       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
634       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
635 #endif
636       dx          *= epsilon;
637       wscale[col] = 1.0/dx;
638       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
639     }
640     VecRestoreArray(x1,&xx);
641     /*
642        Evaluate function at x1 + dx (here dx is a vector of perturbations)
643     */
644     ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr);
645     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
646     /* Communicate scale to all processors */
647 #if !defined(USE_PETSC_COMPLEX)
648     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
649 #else
650     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
651 #endif
652     /*
653        Loop over rows of vector, putting results into Jacobian matrix
654     */
655     VecGetArray(w2,&y);
656     for (l=0; l<coloring->nrows[k]; l++) {
657       row    = coloring->rows[k][l];
658       col    = coloring->columnsforrow[k][l];
659       y[row] *= scale[col];
660       srow   = row + start;
661       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
662     }
663     VecRestoreArray(w2,&y);
664   }
665   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
666   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
667   PetscFunctionReturn(0);
668 }
669