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