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