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