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