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