xref: /petsc/src/mat/matfd/fdmatrix.c (revision 2cfcea2917e9827eeaf4fa9c2e8dd75b054ced7b)
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__ "MatFDColoringGetFunction"
254 /*@C
255    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
256 
257    Collective on MatFDColoring
258 
259    Input Parameters:
260 .  coloring - the coloring context
261 
262    Output Parameters:
263 +  f - the function
264 -  fctx - the optional user-defined function context
265 
266    Level: intermediate
267 
268 .keywords: Mat, Jacobian, finite differences, set, function
269 @*/
270 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
271 {
272   PetscFunctionBegin;
273   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
274   if (f) *f = matfd->f;
275   if (fctx) *fctx = matfd->fctx;
276   PetscFunctionReturn(0);
277 }
278 
279 #undef __FUNCT__
280 #define __FUNCT__ "MatFDColoringSetFunction"
281 /*@C
282    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
283 
284    Collective on MatFDColoring
285 
286    Input Parameters:
287 +  coloring - the coloring context
288 .  f - the function
289 -  fctx - the optional user-defined function context
290 
291    Level: intermediate
292 
293    Notes:
294     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
295   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
296   with the TS solvers.
297 
298 .keywords: Mat, Jacobian, finite differences, set, function
299 @*/
300 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
301 {
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
304   matfd->f    = f;
305   matfd->fctx = fctx;
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "MatFDColoringSetFromOptions"
311 /*@
312    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
313    the options database.
314 
315    Collective on MatFDColoring
316 
317    The Jacobian, F'(u), is estimated with the differencing approximation
318 .vb
319        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
320        h = error_rel*u[i]                 if  abs(u[i]) > umin
321          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
322        dx_{i} = (0, ... 1, .... 0)
323 .ve
324 
325    Input Parameter:
326 .  coloring - the coloring context
327 
328    Options Database Keys:
329 +  -mat_fd_coloring_err <err> - Sets <err> (square root
330            of relative error in the function)
331 .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
332 .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
333 .  -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
334 .  -mat_fd_coloring_view - Activates basic viewing
335 .  -mat_fd_coloring_view_info - Activates viewing info
336 -  -mat_fd_coloring_view_draw - Activates drawing
337 
338     Level: intermediate
339 
340 .keywords: Mat, finite differences, parameters
341 
342 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
343 
344 @*/
345 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
346 {
347   PetscErrorCode ierr;
348   PetscTruth     flg;
349   char           value[3];
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
353 
354   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
355     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
356     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
357     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
358     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
359     if (flg) {
360       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
361       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
362       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
363     }
364     /* not used here; but so they are presented in the GUI */
365     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
366     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
367     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
368   PetscOptionsEnd();CHKERRQ(ierr);
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "MatFDColoringView_Private"
374 PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
375 {
376   PetscErrorCode ierr;
377   PetscTruth     flg;
378 
379   PetscFunctionBegin;
380   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
381   if (flg) {
382     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
383   }
384   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
385   if (flg) {
386     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
387     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
388     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
389   }
390   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
391   if (flg) {
392     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
393     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
394   }
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "MatFDColoringCreate"
400 /*@
401    MatFDColoringCreate - Creates a matrix coloring context for finite difference
402    computation of Jacobians.
403 
404    Collective on Mat
405 
406    Input Parameters:
407 +  mat - the matrix containing the nonzero structure of the Jacobian
408 -  iscoloring - the coloring of the matrix
409 
410     Output Parameter:
411 .   color - the new coloring context
412 
413     Level: intermediate
414 
415 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
416           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
417           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
418           MatFDColoringSetParameters()
419 @*/
420 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
421 {
422   MatFDColoring  c;
423   MPI_Comm       comm;
424   PetscErrorCode ierr;
425   PetscInt       M,N;
426 
427   PetscFunctionBegin;
428   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
429   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
430   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
431 
432   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
433   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
434 
435   if (mat->ops->fdcoloringcreate) {
436     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
437   } else {
438     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
439   }
440 
441   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
442   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
443   c->freq              = 1;
444   c->usersetsrecompute = PETSC_FALSE;
445   c->recompute         = PETSC_FALSE;
446   c->currentcolor      = -1;
447   c->htype             = "wp";
448 
449   *color = c;
450   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "MatFDColoringDestroy"
456 /*@
457     MatFDColoringDestroy - Destroys a matrix coloring context that was created
458     via MatFDColoringCreate().
459 
460     Collective on MatFDColoring
461 
462     Input Parameter:
463 .   c - coloring context
464 
465     Level: intermediate
466 
467 .seealso: MatFDColoringCreate()
468 @*/
469 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
470 {
471   PetscErrorCode ierr;
472   PetscInt       i;
473 
474   PetscFunctionBegin;
475   if (--c->refct > 0) PetscFunctionReturn(0);
476 
477   for (i=0; i<c->ncolors; i++) {
478     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
479     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
480     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
481     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
482   }
483   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
484   ierr = PetscFree(c->columns);CHKERRQ(ierr);
485   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
486   ierr = PetscFree(c->rows);CHKERRQ(ierr);
487   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
488   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
489   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
490   if (c->w1) {
491     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
492     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
493     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
494   }
495   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
496   PetscFunctionReturn(0);
497 }
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
501 /*@C
502     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
503       that are currently being perturbed.
504 
505     Not Collective
506 
507     Input Parameters:
508 .   coloring - coloring context created with MatFDColoringCreate()
509 
510     Output Parameters:
511 +   n - the number of local columns being perturbed
512 -   cols - the column indices, in global numbering
513 
514    Level: intermediate
515 
516 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
517 
518 .keywords: coloring, Jacobian, finite differences
519 @*/
520 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
521 {
522   PetscFunctionBegin;
523   if (coloring->currentcolor >= 0) {
524     *n    = coloring->ncolumns[coloring->currentcolor];
525     *cols = coloring->columns[coloring->currentcolor];
526   } else {
527     *n = 0;
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "MatFDColoringApply"
534 /*@
535     MatFDColoringApply - Given a matrix for which a MatFDColoring context
536     has been created, computes the Jacobian for a function via finite differences.
537 
538     Collective on MatFDColoring
539 
540     Input Parameters:
541 +   mat - location to store Jacobian
542 .   coloring - coloring context created with MatFDColoringCreate()
543 .   x1 - location at which Jacobian is to be computed
544 -   sctx - optional context required by function (actually a SNES context)
545 
546     Options Database Keys:
547 +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
548 .    -mat_fd_type - "wp" or "ds"  (see MATSNESMF_WP or MATSNESMF_DS)
549 .    -mat_fd_coloring_view - Activates basic viewing or coloring
550 .    -mat_fd_coloring_view_draw - Activates drawing of coloring
551 -    -mat_fd_coloring_view_info - Activates viewing of coloring info
552 
553     Level: intermediate
554 
555 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
556 
557 .keywords: coloring, Jacobian, finite differences
558 @*/
559 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
560 {
561   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
562   PetscErrorCode ierr;
563   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
564   PetscScalar    dx,*y,*xx,*w3_array;
565   PetscScalar    *vscale_array;
566   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
567   Vec            w1,w2,w3;
568   void           *fctx = coloring->fctx;
569   PetscTruth     flg;
570 
571 
572   PetscFunctionBegin;
573   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
574   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
575   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
576   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
577 
578   if (coloring->usersetsrecompute) {
579     if (!coloring->recompute) {
580       *flag = SAME_PRECONDITIONER;
581       ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr);
582       PetscFunctionReturn(0);
583     } else {
584       coloring->recompute = PETSC_FALSE;
585     }
586   }
587 
588   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
589   if (J->ops->fdcoloringapply) {
590     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
591   } else {
592 
593     if (!coloring->w1) {
594       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
595       ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
596       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
597       ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
598       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
599       ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
600     }
601     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
602 
603     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
604     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
605     if (flg) {
606       ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
607     } else {
608       PetscTruth assembled;
609       ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
610       if (assembled) {
611 	ierr = MatZeroEntries(J);CHKERRQ(ierr);
612       }
613     }
614 
615     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
616     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
617 
618     /*
619       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
620       coloring->F for the coarser grids from the finest
621     */
622     if (coloring->F) {
623       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
624       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
625       if (m1 != m2) {
626 	coloring->F = 0;
627       }
628     }
629 
630     if (coloring->F) {
631       w1          = coloring->F; /* use already computed value of function */
632       coloring->F = 0;
633     } else {
634       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
635       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
636       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
637     }
638 
639     if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
640       ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
641     }
642 
643     /*
644        Compute all the scale factors and share with other processors
645     */
646     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
647     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
648     for (k=0; k<coloring->ncolors; k++) {
649       /*
650 	Loop over each column associated with color adding the
651 	perturbation to the vector w3.
652       */
653       for (l=0; l<coloring->ncolumns[k]; l++) {
654 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
655         if (coloring->htype[0] == 'w') {
656           dx = 1.0 + unorm;
657         } else {
658   	  dx  = xx[col];
659         }
660 	if (dx == 0.0) dx = 1.0;
661 #if !defined(PETSC_USE_COMPLEX)
662 	if (dx < umin && dx >= 0.0)      dx = umin;
663 	else if (dx < 0.0 && dx > -umin) dx = -umin;
664 #else
665 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
666 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
667 #endif
668 	dx                *= epsilon;
669 	vscale_array[col] = 1.0/dx;
670       }
671     }
672     vscale_array = vscale_array + start;
673     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
674     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
675     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
676 
677     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
678 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
679 
680     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
681     else                        vscaleforrow = coloring->columnsforrow;
682 
683     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
684     /*
685       Loop over each color
686     */
687     for (k=0; k<coloring->ncolors; k++) {
688       coloring->currentcolor = k;
689       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
690       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
691       /*
692 	Loop over each column associated with color adding the
693 	perturbation to the vector w3.
694       */
695       for (l=0; l<coloring->ncolumns[k]; l++) {
696 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
697         if (coloring->htype[0] == 'w') {
698           dx = 1.0 + unorm;
699         } else {
700   	  dx  = xx[col];
701         }
702 	if (dx == 0.0) dx = 1.0;
703 #if !defined(PETSC_USE_COMPLEX)
704 	if (dx < umin && dx >= 0.0)      dx = umin;
705 	else if (dx < 0.0 && dx > -umin) dx = -umin;
706 #else
707 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
708 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
709 #endif
710 	dx            *= epsilon;
711 	if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
712 	w3_array[col] += dx;
713       }
714       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
715 
716       /*
717 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
718       */
719 
720       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
721       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
722       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
723       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
724 
725       /*
726 	Loop over rows of vector, putting results into Jacobian matrix
727       */
728       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
729       for (l=0; l<coloring->nrows[k]; l++) {
730 	row    = coloring->rows[k][l];
731 	col    = coloring->columnsforrow[k][l];
732 	y[row] *= vscale_array[vscaleforrow[k][l]];
733 	srow   = row + start;
734 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
735       }
736       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
737     }
738     coloring->currentcolor = -1;
739     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
740     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
741     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
742     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
743   }
744   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
745 
746   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
747   if (flg) {
748     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
749   }
750   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
751 
752   PetscFunctionReturn(0);
753 }
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "MatFDColoringApplyTS"
757 /*@
758     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
759     has been created, computes the Jacobian for a function via finite differences.
760 
761    Collective on Mat, MatFDColoring, and Vec
762 
763     Input Parameters:
764 +   mat - location to store Jacobian
765 .   coloring - coloring context created with MatFDColoringCreate()
766 .   x1 - location at which Jacobian is to be computed
767 -   sctx - optional context required by function (actually a SNES context)
768 
769    Options Database Keys:
770 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
771 
772    Level: intermediate
773 
774 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
775 
776 .keywords: coloring, Jacobian, finite differences
777 @*/
778 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
779 {
780   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
781   PetscErrorCode ierr;
782   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
783   PetscScalar    dx,*y,*xx,*w3_array;
784   PetscScalar    *vscale_array;
785   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
786   Vec            w1,w2,w3;
787   void           *fctx = coloring->fctx;
788   PetscTruth     flg;
789 
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
792   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
793   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
794 
795   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
796   if (!coloring->w1) {
797     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
798     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
799     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
800     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
801     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
802     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
803   }
804   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
805 
806   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
807   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
808   if (flg) {
809     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
810   } else {
811     PetscTruth assembled;
812     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
813     if (assembled) {
814       ierr = MatZeroEntries(J);CHKERRQ(ierr);
815     }
816   }
817 
818   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
819   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
820   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
821   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
822   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
823 
824   /*
825       Compute all the scale factors and share with other processors
826   */
827   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
828   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
829   for (k=0; k<coloring->ncolors; k++) {
830     /*
831        Loop over each column associated with color adding the
832        perturbation to the vector w3.
833     */
834     for (l=0; l<coloring->ncolumns[k]; l++) {
835       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
836       dx  = xx[col];
837       if (dx == 0.0) dx = 1.0;
838 #if !defined(PETSC_USE_COMPLEX)
839       if (dx < umin && dx >= 0.0)      dx = umin;
840       else if (dx < 0.0 && dx > -umin) dx = -umin;
841 #else
842       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
843       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
844 #endif
845       dx                *= epsilon;
846       vscale_array[col] = 1.0/dx;
847     }
848   }
849   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
850   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
851   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
852 
853   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
854   else                        vscaleforrow = coloring->columnsforrow;
855 
856   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
857   /*
858       Loop over each color
859   */
860   for (k=0; k<coloring->ncolors; k++) {
861     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
862     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
863     /*
864        Loop over each column associated with color adding the
865        perturbation to the vector w3.
866     */
867     for (l=0; l<coloring->ncolumns[k]; l++) {
868       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
869       dx  = xx[col];
870       if (dx == 0.0) dx = 1.0;
871 #if !defined(PETSC_USE_COMPLEX)
872       if (dx < umin && dx >= 0.0)      dx = umin;
873       else if (dx < 0.0 && dx > -umin) dx = -umin;
874 #else
875       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
876       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
877 #endif
878       dx            *= epsilon;
879       w3_array[col] += dx;
880     }
881     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
882 
883     /*
884        Evaluate function at x1 + dx (here dx is a vector of perturbations)
885     */
886     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
887     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
888     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
889     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
890 
891     /*
892        Loop over rows of vector, putting results into Jacobian matrix
893     */
894     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
895     for (l=0; l<coloring->nrows[k]; l++) {
896       row    = coloring->rows[k][l];
897       col    = coloring->columnsforrow[k][l];
898       y[row] *= vscale_array[vscaleforrow[k][l]];
899       srow   = row + start;
900       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
901     }
902     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
903   }
904   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
905   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
906   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
907   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
908   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
909   PetscFunctionReturn(0);
910 }
911 
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "MatFDColoringSetRecompute()"
915 /*@C
916    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
917      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
918      no additional Jacobian's will be computed (the same one will be used) until this is
919      called again.
920 
921    Collective on MatFDColoring
922 
923    Input  Parameters:
924 .  c - the coloring context
925 
926    Level: intermediate
927 
928    Notes: The MatFDColoringSetFrequency() is ignored once this is called
929 
930 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
931 
932 .keywords: Mat, finite differences, coloring
933 @*/
934 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
935 {
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
938   c->usersetsrecompute = PETSC_TRUE;
939   c->recompute         = PETSC_TRUE;
940   PetscFunctionReturn(0);
941 }
942 
943 
944