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