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