xref: /petsc/src/mat/matfd/fdmatrix.c (revision 09f3b4e5628a00a1eaf17d80982cfbcc515cc9c1)
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     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
479     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
480     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
481     if (c->vscaleforrow && c->vscaleforrow[i]) {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   if (c->vscaleforrow) {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 
577   if (coloring->usersetsrecompute) {
578     if (!coloring->recompute) {
579       *flag = SAME_PRECONDITIONER;
580       ierr = PetscVerboseInfo((J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"));CHKERRQ(ierr);
581       PetscFunctionReturn(0);
582     } else {
583       coloring->recompute = PETSC_FALSE;
584     }
585   }
586 
587   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
588   if (J->ops->fdcoloringapply) {
589     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
590   } else {
591 
592     if (!coloring->w1) {
593       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
594       ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
595       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
596       ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
597       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
598       ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
599     }
600     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
601 
602     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
603     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
604     if (flg) {
605       ierr = PetscVerboseInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
606     } else {
607       PetscTruth assembled;
608       ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
609       if (assembled) {
610 	ierr = MatZeroEntries(J);CHKERRQ(ierr);
611       }
612     }
613 
614     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
615     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
616 
617     /*
618       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
619       coloring->F for the coarser grids from the finest
620     */
621     if (coloring->F) {
622       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
623       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
624       if (m1 != m2) {
625 	coloring->F = 0;
626       }
627     }
628 
629     if (coloring->F) {
630       w1          = coloring->F; /* use already computed value of function */
631       coloring->F = 0;
632     } else {
633       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
634       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
635       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
636     }
637 
638     if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
639       ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
640     }
641 
642     /*
643        Compute all the scale factors and share with other processors
644     */
645     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
646     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
647     for (k=0; k<coloring->ncolors; k++) {
648       /*
649 	Loop over each column associated with color adding the
650 	perturbation to the vector w3.
651       */
652       for (l=0; l<coloring->ncolumns[k]; l++) {
653 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
654         if (coloring->htype[0] == 'w') {
655           dx = 1.0 + unorm;
656         } else {
657   	  dx  = xx[col];
658         }
659 	if (dx == 0.0) dx = 1.0;
660 #if !defined(PETSC_USE_COMPLEX)
661 	if (dx < umin && dx >= 0.0)      dx = umin;
662 	else if (dx < 0.0 && dx > -umin) dx = -umin;
663 #else
664 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
665 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
666 #endif
667 	dx                *= epsilon;
668 	vscale_array[col] = 1.0/dx;
669       }
670     }
671     vscale_array = vscale_array + start;
672     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
673     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
675 
676     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
677 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
678 
679     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
680     else                        vscaleforrow = coloring->columnsforrow;
681 
682     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
683     /*
684       Loop over each color
685     */
686     for (k=0; k<coloring->ncolors; k++) {
687       coloring->currentcolor = k;
688       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
689       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
690       /*
691 	Loop over each column associated with color adding the
692 	perturbation to the vector w3.
693       */
694       for (l=0; l<coloring->ncolumns[k]; l++) {
695 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
696         if (coloring->htype[0] == 'w') {
697           dx = 1.0 + unorm;
698         } else {
699   	  dx  = xx[col];
700         }
701 	if (dx == 0.0) dx = 1.0;
702 #if !defined(PETSC_USE_COMPLEX)
703 	if (dx < umin && dx >= 0.0)      dx = umin;
704 	else if (dx < 0.0 && dx > -umin) dx = -umin;
705 #else
706 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
707 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
708 #endif
709 	dx            *= epsilon;
710 	if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
711 	w3_array[col] += dx;
712       }
713       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
714 
715       /*
716 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
717       */
718 
719       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
720       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
721       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
722       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
723 
724       /*
725 	Loop over rows of vector, putting results into Jacobian matrix
726       */
727       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
728       for (l=0; l<coloring->nrows[k]; l++) {
729 	row    = coloring->rows[k][l];
730 	col    = coloring->columnsforrow[k][l];
731 	y[row] *= vscale_array[vscaleforrow[k][l]];
732 	srow   = row + start;
733 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
734       }
735       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
736     }
737     coloring->currentcolor = -1;
738     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
739     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
740     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
741     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
742   }
743   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
744 
745   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
746   if (flg) {
747     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
748   }
749   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
750 
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "MatFDColoringApplyTS"
756 /*@
757     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
758     has been created, computes the Jacobian for a function via finite differences.
759 
760    Collective on Mat, MatFDColoring, and Vec
761 
762     Input Parameters:
763 +   mat - location to store Jacobian
764 .   coloring - coloring context created with MatFDColoringCreate()
765 .   x1 - location at which Jacobian is to be computed
766 -   sctx - optional context required by function (actually a SNES context)
767 
768    Options Database Keys:
769 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
770 
771    Level: intermediate
772 
773 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
774 
775 .keywords: coloring, Jacobian, finite differences
776 @*/
777 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
778 {
779   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
780   PetscErrorCode ierr;
781   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
782   PetscScalar    dx,*y,*xx,*w3_array;
783   PetscScalar    *vscale_array;
784   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
785   Vec            w1,w2,w3;
786   void           *fctx = coloring->fctx;
787   PetscTruth     flg;
788 
789   PetscFunctionBegin;
790   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
791   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
792   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
793 
794   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
795   if (!coloring->w1) {
796     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
797     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
798     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
799     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
800     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
801     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
802   }
803   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
804 
805   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
806   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
807   if (flg) {
808     ierr = PetscVerboseInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
809   } else {
810     PetscTruth assembled;
811     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
812     if (assembled) {
813       ierr = MatZeroEntries(J);CHKERRQ(ierr);
814     }
815   }
816 
817   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
818   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
819   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
820   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
821   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
822 
823   /*
824       Compute all the scale factors and share with other processors
825   */
826   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
827   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
828   for (k=0; k<coloring->ncolors; k++) {
829     /*
830        Loop over each column associated with color adding the
831        perturbation to the vector w3.
832     */
833     for (l=0; l<coloring->ncolumns[k]; l++) {
834       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
835       dx  = xx[col];
836       if (dx == 0.0) dx = 1.0;
837 #if !defined(PETSC_USE_COMPLEX)
838       if (dx < umin && dx >= 0.0)      dx = umin;
839       else if (dx < 0.0 && dx > -umin) dx = -umin;
840 #else
841       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
842       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
843 #endif
844       dx                *= epsilon;
845       vscale_array[col] = 1.0/dx;
846     }
847   }
848   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
849   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
850   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
851 
852   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
853   else                        vscaleforrow = coloring->columnsforrow;
854 
855   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
856   /*
857       Loop over each color
858   */
859   for (k=0; k<coloring->ncolors; k++) {
860     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
861     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
862     /*
863        Loop over each column associated with color adding the
864        perturbation to the vector w3.
865     */
866     for (l=0; l<coloring->ncolumns[k]; l++) {
867       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
868       dx  = xx[col];
869       if (dx == 0.0) dx = 1.0;
870 #if !defined(PETSC_USE_COMPLEX)
871       if (dx < umin && dx >= 0.0)      dx = umin;
872       else if (dx < 0.0 && dx > -umin) dx = -umin;
873 #else
874       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
875       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
876 #endif
877       dx            *= epsilon;
878       w3_array[col] += dx;
879     }
880     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
881 
882     /*
883        Evaluate function at x1 + dx (here dx is a vector of perturbations)
884     */
885     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
886     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
887     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
888     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
889 
890     /*
891        Loop over rows of vector, putting results into Jacobian matrix
892     */
893     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
894     for (l=0; l<coloring->nrows[k]; l++) {
895       row    = coloring->rows[k][l];
896       col    = coloring->columnsforrow[k][l];
897       y[row] *= vscale_array[vscaleforrow[k][l]];
898       srow   = row + start;
899       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
900     }
901     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
902   }
903   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
904   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
905   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
906   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
907   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
908   PetscFunctionReturn(0);
909 }
910 
911 
912 #undef __FUNCT__
913 #define __FUNCT__ "MatFDColoringSetRecompute()"
914 /*@C
915    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
916      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
917      no additional Jacobian's will be computed (the same one will be used) until this is
918      called again.
919 
920    Collective on MatFDColoring
921 
922    Input  Parameters:
923 .  c - the coloring context
924 
925    Level: intermediate
926 
927    Notes: The MatFDColoringSetFrequency() is ignored once this is called
928 
929 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
930 
931 .keywords: Mat, finite differences, coloring
932 @*/
933 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
934 {
935   PetscFunctionBegin;
936   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
937   c->usersetsrecompute = PETSC_TRUE;
938   c->recompute         = PETSC_TRUE;
939   PetscFunctionReturn(0);
940 }
941 
942 
943