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