xref: /petsc/src/mat/matfd/fdmatrix.c (revision 52e6d16ba989af9362d0fcebb12e73dd7c9666ed)
1 /*
2    This is where the abstract matrix operations are defined that are
3   used for finite difference computations of Jacobians using coloring.
4 */
5 
6 #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
7 
8 /* Logging support */
9 PetscCookie MAT_FDCOLORING_COOKIE = 0;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatFDColoringSetF"
13 PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F)
14 {
15   PetscFunctionBegin;
16   fd->F = F;
17   PetscFunctionReturn(0);
18 }
19 
20 #undef __FUNCT__
21 #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
22 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
23 {
24   MatFDColoring  fd = (MatFDColoring)Aa;
25   PetscErrorCode ierr;
26   PetscInt       i,j;
27   PetscReal      x,y;
28 
29   PetscFunctionBegin;
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 = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
37     }
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "MatFDColoringView_Draw"
44 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
45 {
46   PetscErrorCode ierr;
47   PetscTruth     isnull;
48   PetscDraw      draw;
49   PetscReal      xr,yr,xl,yl,h,w;
50 
51   PetscFunctionBegin;
52   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
53   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
54 
55   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
56 
57   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
58   xr += w;     yr += h;    xl = -w;     yl = -h;
59   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
60   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
61   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "MatFDColoringView"
67 /*@C
68    MatFDColoringView - Views a finite difference coloring context.
69 
70    Collective on MatFDColoring
71 
72    Input  Parameters:
73 +  c - the coloring context
74 -  viewer - visualization context
75 
76    Level: intermediate
77 
78    Notes:
79    The available visualization contexts include
80 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
81 .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
82         output where only the first processor opens
83         the file.  All other processors send their
84         data to the first processor to print.
85 -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
86 
87    Notes:
88      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
89    involves more than 33 then some seemingly identical colors are displayed making it look
90    like an illegal coloring. This is just a graphical artifact.
91 
92 .seealso: MatFDColoringCreate()
93 
94 .keywords: Mat, finite differences, coloring, view
95 @*/
96 PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer)
97 {
98   PetscErrorCode    ierr;
99   PetscInt          i,j;
100   PetscTruth        isdraw,iascii;
101   PetscViewerFormat format;
102 
103   PetscFunctionBegin;
104   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
105   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
106   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
107   PetscCheckSameComm(c,1,viewer,2);
108 
109   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
110   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
111   if (isdraw) {
112     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
113   } else if (iascii) {
114     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
115     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
116     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
117     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
118 
119     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
120     if (format != PETSC_VIEWER_ASCII_INFO) {
121       for (i=0; i<c->ncolors; i++) {
122         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
123         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
124         for (j=0; j<c->ncolumns[i]; j++) {
125           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
126         }
127         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
128         for (j=0; j<c->nrows[i]; j++) {
129           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
130         }
131       }
132     }
133     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
134   } else {
135     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
136   }
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "MatFDColoringSetParameters"
142 /*@
143    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144    a Jacobian matrix using finite differences.
145 
146    Collective on MatFDColoring
147 
148    The Jacobian is estimated with the differencing approximation
149 .vb
150        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151        h = error_rel*u[i]                 if  abs(u[i]) > umin
152          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
153        dx_{i} = (0, ... 1, .... 0)
154 .ve
155 
156    Input Parameters:
157 +  coloring - the coloring context
158 .  error_rel - relative error
159 -  umin - minimum allowable u-value magnitude
160 
161    Level: advanced
162 
163 .keywords: Mat, finite differences, coloring, set, parameters
164 
165 .seealso: MatFDColoringCreate()
166 @*/
167 PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
168 {
169   PetscFunctionBegin;
170   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
171 
172   if (error != PETSC_DEFAULT) matfd->error_rel = error;
173   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "MatFDColoringSetFrequency"
179 /*@
180    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
181    matrices.
182 
183    Collective on MatFDColoring
184 
185    Input Parameters:
186 +  coloring - the coloring context
187 -  freq - frequency (default is 1)
188 
189    Options Database Keys:
190 .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
191 
192    Level: advanced
193 
194    Notes:
195    Using a modified Newton strategy, where the Jacobian remains fixed for several
196    iterations, can be cost effective in terms of overall nonlinear solution
197    efficiency.  This parameter indicates that a new Jacobian will be computed every
198    <freq> nonlinear iterations.
199 
200 .keywords: Mat, finite differences, coloring, set, frequency
201 
202 .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
203 @*/
204 PetscErrorCode MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq)
205 {
206   PetscFunctionBegin;
207   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
208 
209   matfd->freq = freq;
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "MatFDColoringGetFrequency"
215 /*@
216    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
217    matrices.
218 
219    Not Collective
220 
221    Input Parameters:
222 .  coloring - the coloring context
223 
224    Output Parameters:
225 .  freq - frequency (default is 1)
226 
227    Options Database Keys:
228 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
229 
230    Level: advanced
231 
232    Notes:
233    Using a modified Newton strategy, where the Jacobian remains fixed for several
234    iterations, can be cost effective in terms of overall nonlinear solution
235    efficiency.  This parameter indicates that a new Jacobian will be computed every
236    <freq> nonlinear iterations.
237 
238 .keywords: Mat, finite differences, coloring, get, frequency
239 
240 .seealso: MatFDColoringSetFrequency()
241 @*/
242 PetscErrorCode MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq)
243 {
244   PetscFunctionBegin;
245   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
246   *freq = matfd->freq;
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "MatFDColoringSetFunction"
252 /*@C
253    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
254 
255    Collective on MatFDColoring
256 
257    Input Parameters:
258 +  coloring - the coloring context
259 .  f - the function
260 -  fctx - the optional user-defined function context
261 
262    Level: intermediate
263 
264    Notes:
265     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
266   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
267   with the TS solvers.
268 
269 .keywords: Mat, Jacobian, finite differences, set, function
270 @*/
271 PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
272 {
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
275   matfd->f    = f;
276   matfd->fctx = fctx;
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "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_err <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 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
313 
314 @*/
315 PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
321 
322   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
323     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
324     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
325     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
326     /* not used here; but so they are presented in the GUI */
327     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
328     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
329     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
330   PetscOptionsEnd();CHKERRQ(ierr);
331   PetscFunctionReturn(0);
332 }
333 
334 #undef __FUNCT__
335 #define __FUNCT__ "MatFDColoringView_Private"
336 PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
337 {
338   PetscErrorCode ierr;
339   PetscTruth     flg;
340 
341   PetscFunctionBegin;
342   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
343   if (flg) {
344     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
345   }
346   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
347   if (flg) {
348     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
349     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
350     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
351   }
352   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
353   if (flg) {
354     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
355     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
356   }
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "MatFDColoringCreate"
362 /*@C
363    MatFDColoringCreate - Creates a matrix coloring context for finite difference
364    computation of Jacobians.
365 
366    Collective on Mat
367 
368    Input Parameters:
369 +  mat - the matrix containing the nonzero structure of the Jacobian
370 -  iscoloring - the coloring of the matrix
371 
372     Output Parameter:
373 .   color - the new coloring context
374 
375     Level: intermediate
376 
377 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
378           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
379           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
380           MatFDColoringSetParameters()
381 @*/
382 PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
383 {
384   MatFDColoring  c;
385   MPI_Comm       comm;
386   PetscErrorCode ierr;
387   PetscInt       M,N;
388 
389   PetscFunctionBegin;
390   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
391   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
392   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
393 
394   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
395   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
396 
397   if (mat->ops->fdcoloringcreate) {
398     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
399   } else {
400     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
401   }
402 
403   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
404   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
405   c->freq              = 1;
406   c->usersetsrecompute = PETSC_FALSE;
407   c->recompute         = PETSC_FALSE;
408   c->currentcolor      = -1;
409 
410   *color = c;
411   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
412   PetscFunctionReturn(0);
413 }
414 
415 #undef __FUNCT__
416 #define __FUNCT__ "MatFDColoringDestroy"
417 /*@C
418     MatFDColoringDestroy - Destroys a matrix coloring context that was created
419     via MatFDColoringCreate().
420 
421     Collective on MatFDColoring
422 
423     Input Parameter:
424 .   c - coloring context
425 
426     Level: intermediate
427 
428 .seealso: MatFDColoringCreate()
429 @*/
430 PetscErrorCode MatFDColoringDestroy(MatFDColoring c)
431 {
432   PetscErrorCode ierr;
433   PetscInt       i;
434 
435   PetscFunctionBegin;
436   if (--c->refct > 0) PetscFunctionReturn(0);
437 
438   for (i=0; i<c->ncolors; i++) {
439     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
440     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
441     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
442     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
443   }
444   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
445   ierr = PetscFree(c->columns);CHKERRQ(ierr);
446   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
447   ierr = PetscFree(c->rows);CHKERRQ(ierr);
448   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
449   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
450   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
451   if (c->w1) {
452     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
453     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
454     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
455   }
456   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
462 /*@C
463     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
464       that are currently being perturbed.
465 
466     Not Collective
467 
468     Input Parameters:
469 .   coloring - coloring context created with MatFDColoringCreate()
470 
471     Output Parameters:
472 +   n - the number of local columns being perturbed
473 -   cols - the column indices, in global numbering
474 
475    Level: intermediate
476 
477 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
478 
479 .keywords: coloring, Jacobian, finite differences
480 @*/
481 EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
482 {
483   PetscFunctionBegin;
484   if (coloring->currentcolor >= 0) {
485     *n    = coloring->ncolumns[coloring->currentcolor];
486     *cols = coloring->columns[coloring->currentcolor];
487   } else {
488     *n = 0;
489   }
490   PetscFunctionReturn(0);
491 }
492 
493 #undef __FUNCT__
494 #define __FUNCT__ "MatFDColoringApply"
495 /*@
496     MatFDColoringApply - Given a matrix for which a MatFDColoring context
497     has been created, computes the Jacobian for a function via finite differences.
498 
499     Collective on MatFDColoring
500 
501     Input Parameters:
502 +   mat - location to store Jacobian
503 .   coloring - coloring context created with MatFDColoringCreate()
504 .   x1 - location at which Jacobian is to be computed
505 -   sctx - optional context required by function (actually a SNES context)
506 
507     Options Database Keys:
508 +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
509 .    -mat_fd_coloring_view - Activates basic viewing or coloring
510 .    -mat_fd_coloring_view_draw - Activates drawing of coloring
511 -    -mat_fd_coloring_view_info - Activates viewing of coloring info
512 
513     Level: intermediate
514 
515 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
516 
517 .keywords: coloring, Jacobian, finite differences
518 @*/
519 PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
520 {
521   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
522   PetscErrorCode ierr;
523   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
524   PetscScalar    dx,mone = -1.0,*y,*xx,*w3_array;
525   PetscScalar    *vscale_array;
526   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
527   Vec            w1,w2,w3;
528   void           *fctx = coloring->fctx;
529   PetscTruth     flg;
530 
531 
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
534   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
535   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
536 
537   if (coloring->usersetsrecompute) {
538     if (!coloring->recompute) {
539       *flag = SAME_PRECONDITIONER;
540       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
541       PetscFunctionReturn(0);
542     } else {
543       coloring->recompute = PETSC_FALSE;
544     }
545   }
546 
547   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
548   if (J->ops->fdcoloringapply) {
549     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
550   } else {
551 
552     if (!coloring->w1) {
553       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
554       ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
555       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
556       ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
557       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
558       ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
559     }
560     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
561 
562     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
563     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
564     if (flg) {
565       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
566     } else {
567       PetscTruth assembled;
568       ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
569       if (assembled) {
570 	ierr = MatZeroEntries(J);CHKERRQ(ierr);
571       }
572     }
573 
574     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
575     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
576 
577     /*
578       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
579       coloring->F for the coarser grids from the finest
580     */
581     if (coloring->F) {
582       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
583       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
584       if (m1 != m2) {
585 	coloring->F = 0;
586       }
587     }
588 
589     if (coloring->F) {
590       w1          = coloring->F; /* use already computed value of function */
591       coloring->F = 0;
592     } else {
593       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
594       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
595       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
596     }
597 
598     /*
599        Compute all the scale factors and share with other processors
600     */
601     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
602     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
603     for (k=0; k<coloring->ncolors; k++) {
604       /*
605 	Loop over each column associated with color adding the
606 	perturbation to the vector w3.
607       */
608       for (l=0; l<coloring->ncolumns[k]; l++) {
609 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
610 	dx  = xx[col];
611 	if (dx == 0.0) dx = 1.0;
612 #if !defined(PETSC_USE_COMPLEX)
613 	if (dx < umin && dx >= 0.0)      dx = umin;
614 	else if (dx < 0.0 && dx > -umin) dx = -umin;
615 #else
616 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
617 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
618 #endif
619 	dx                *= epsilon;
620 	vscale_array[col] = 1.0/dx;
621       }
622     }
623     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
624     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
625     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626 
627     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
628 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
629 
630     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
631     else                        vscaleforrow = coloring->columnsforrow;
632 
633     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
634     /*
635       Loop over each color
636     */
637     for (k=0; k<coloring->ncolors; k++) {
638       coloring->currentcolor = k;
639       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
640       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
641       /*
642 	Loop over each column associated with color adding the
643 	perturbation to the vector w3.
644       */
645       for (l=0; l<coloring->ncolumns[k]; l++) {
646 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
647 	dx  = xx[col];
648 	if (dx == 0.0) dx = 1.0;
649 #if !defined(PETSC_USE_COMPLEX)
650 	if (dx < umin && dx >= 0.0)      dx = umin;
651 	else if (dx < 0.0 && dx > -umin) dx = -umin;
652 #else
653 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
654 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
655 #endif
656 	dx            *= epsilon;
657 	if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
658 	w3_array[col] += dx;
659       }
660       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
661 
662       /*
663 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
664       */
665 
666       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
667       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
668       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
669       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
670 
671       /*
672 	Loop over rows of vector, putting results into Jacobian matrix
673       */
674       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
675       for (l=0; l<coloring->nrows[k]; l++) {
676 	row    = coloring->rows[k][l];
677 	col    = coloring->columnsforrow[k][l];
678 	y[row] *= vscale_array[vscaleforrow[k][l]];
679 	srow   = row + start;
680 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
681       }
682       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
683     }
684     coloring->currentcolor = -1;
685     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
686     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
687     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
688     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689   }
690   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
691 
692   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
693   if (flg) {
694     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
695   }
696   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
697 
698   PetscFunctionReturn(0);
699 }
700 
701 #undef __FUNCT__
702 #define __FUNCT__ "MatFDColoringApplyTS"
703 /*@
704     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
705     has been created, computes the Jacobian for a function via finite differences.
706 
707    Collective on Mat, MatFDColoring, and Vec
708 
709     Input Parameters:
710 +   mat - location to store Jacobian
711 .   coloring - coloring context created with MatFDColoringCreate()
712 .   x1 - location at which Jacobian is to be computed
713 -   sctx - optional context required by function (actually a SNES context)
714 
715    Options Database Keys:
716 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
717 
718    Level: intermediate
719 
720 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
721 
722 .keywords: coloring, Jacobian, finite differences
723 @*/
724 PetscErrorCode MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
725 {
726   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
727   PetscErrorCode ierr;
728   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
729   PetscScalar    dx,mone = -1.0,*y,*xx,*w3_array;
730   PetscScalar    *vscale_array;
731   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
732   Vec            w1,w2,w3;
733   void           *fctx = coloring->fctx;
734   PetscTruth     flg;
735 
736   PetscFunctionBegin;
737   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
738   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
739   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
740 
741   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
742   if (!coloring->w1) {
743     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
744     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
745     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
746     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
747     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
748     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
749   }
750   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
751 
752   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
753   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
754   if (flg) {
755     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
756   } else {
757     PetscTruth assembled;
758     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
759     if (assembled) {
760       ierr = MatZeroEntries(J);CHKERRQ(ierr);
761     }
762   }
763 
764   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
765   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
766   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
767   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
768   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
769 
770   /*
771       Compute all the scale factors and share with other processors
772   */
773   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
774   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
775   for (k=0; k<coloring->ncolors; k++) {
776     /*
777        Loop over each column associated with color adding the
778        perturbation to the vector w3.
779     */
780     for (l=0; l<coloring->ncolumns[k]; l++) {
781       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
782       dx  = xx[col];
783       if (dx == 0.0) dx = 1.0;
784 #if !defined(PETSC_USE_COMPLEX)
785       if (dx < umin && dx >= 0.0)      dx = umin;
786       else if (dx < 0.0 && dx > -umin) dx = -umin;
787 #else
788       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
789       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
790 #endif
791       dx                *= epsilon;
792       vscale_array[col] = 1.0/dx;
793     }
794   }
795   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
796   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
797   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
798 
799   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
800   else                        vscaleforrow = coloring->columnsforrow;
801 
802   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
803   /*
804       Loop over each color
805   */
806   for (k=0; k<coloring->ncolors; k++) {
807     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
808     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
809     /*
810        Loop over each column associated with color adding the
811        perturbation to the vector w3.
812     */
813     for (l=0; l<coloring->ncolumns[k]; l++) {
814       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
815       dx  = xx[col];
816       if (dx == 0.0) dx = 1.0;
817 #if !defined(PETSC_USE_COMPLEX)
818       if (dx < umin && dx >= 0.0)      dx = umin;
819       else if (dx < 0.0 && dx > -umin) dx = -umin;
820 #else
821       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
822       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
823 #endif
824       dx            *= epsilon;
825       w3_array[col] += dx;
826     }
827     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
828 
829     /*
830        Evaluate function at x1 + dx (here dx is a vector of perturbations)
831     */
832     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
833     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
834     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
835     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
836 
837     /*
838        Loop over rows of vector, putting results into Jacobian matrix
839     */
840     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
841     for (l=0; l<coloring->nrows[k]; l++) {
842       row    = coloring->rows[k][l];
843       col    = coloring->columnsforrow[k][l];
844       y[row] *= vscale_array[vscaleforrow[k][l]];
845       srow   = row + start;
846       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
847     }
848     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
849   }
850   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
851   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
852   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
853   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
854   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "MatFDColoringSetRecompute()"
861 /*@C
862    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
863      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
864      no additional Jacobian's will be computed (the same one will be used) until this is
865      called again.
866 
867    Collective on MatFDColoring
868 
869    Input  Parameters:
870 .  c - the coloring context
871 
872    Level: intermediate
873 
874    Notes: The MatFDColoringSetFrequency() is ignored once this is called
875 
876 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
877 
878 .keywords: Mat, finite differences, coloring
879 @*/
880 PetscErrorCode MatFDColoringSetRecompute(MatFDColoring c)
881 {
882   PetscFunctionBegin;
883   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
884   c->usersetsrecompute = PETSC_TRUE;
885   c->recompute         = PETSC_TRUE;
886   PetscFunctionReturn(0);
887 }
888 
889 
890