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