xref: /petsc/src/mat/matfd/fdmatrix.c (revision 3ab47c1a6b8969b08cfdbd2a4b0fa9b1ca65d2f8)
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   PetscMPIInt    size;
427 
428   PetscFunctionBegin;
429   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
430   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
431   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
432 
433   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
434   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
435 
436   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
437   c->ctype = iscoloring->ctype;
438 
439   if (mat->ops->fdcoloringcreate) {
440     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
441   } else {
442     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
443   }
444 
445   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
446   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
447   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
448   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
449 
450   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
451   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
452   c->freq              = 1;
453   c->usersetsrecompute = PETSC_FALSE;
454   c->recompute         = PETSC_FALSE;
455   c->currentcolor      = -1;
456   c->htype             = "wp";
457 
458   *color = c;
459   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
460   PetscFunctionReturn(0);
461 }
462 
463 #undef __FUNCT__
464 #define __FUNCT__ "MatFDColoringDestroy"
465 /*@
466     MatFDColoringDestroy - Destroys a matrix coloring context that was created
467     via MatFDColoringCreate().
468 
469     Collective on MatFDColoring
470 
471     Input Parameter:
472 .   c - coloring context
473 
474     Level: intermediate
475 
476 .seealso: MatFDColoringCreate()
477 @*/
478 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
479 {
480   PetscErrorCode ierr;
481   PetscInt       i;
482 
483   PetscFunctionBegin;
484   if (--c->refct > 0) PetscFunctionReturn(0);
485 
486   for (i=0; i<c->ncolors; i++) {
487     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
488     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
489     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
490     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
491   }
492   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
493   ierr = PetscFree(c->columns);CHKERRQ(ierr);
494   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
495   ierr = PetscFree(c->rows);CHKERRQ(ierr);
496   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
497   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
498   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
499   if (c->w1) {
500     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
501     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
502   }
503   if (c->w3){
504     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
505   }
506   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
512 /*@C
513     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
514       that are currently being perturbed.
515 
516     Not Collective
517 
518     Input Parameters:
519 .   coloring - coloring context created with MatFDColoringCreate()
520 
521     Output Parameters:
522 +   n - the number of local columns being perturbed
523 -   cols - the column indices, in global numbering
524 
525    Level: intermediate
526 
527 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
528 
529 .keywords: coloring, Jacobian, finite differences
530 @*/
531 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
532 {
533   PetscFunctionBegin;
534   if (coloring->currentcolor >= 0) {
535     *n    = coloring->ncolumns[coloring->currentcolor];
536     *cols = coloring->columns[coloring->currentcolor];
537   } else {
538     *n = 0;
539   }
540   PetscFunctionReturn(0);
541 }
542 
543 #include "petscda.h"      /*I      "petscda.h"    I*/
544 
545 #undef __FUNCT__
546 #define __FUNCT__ "MatFDColoringApply"
547 /*@
548     MatFDColoringApply - Given a matrix for which a MatFDColoring context
549     has been created, computes the Jacobian for a function via finite differences.
550 
551     Collective on MatFDColoring
552 
553     Input Parameters:
554 +   mat - location to store Jacobian
555 .   coloring - coloring context created with MatFDColoringCreate()
556 .   x1 - location at which Jacobian is to be computed
557 -   sctx - optional context required by function (actually a SNES context)
558 
559     Options Database Keys:
560 +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
561 .    -mat_fd_type - "wp" or "ds"  (see MATSNESMF_WP or MATSNESMF_DS)
562 .    -mat_fd_coloring_view - Activates basic viewing or coloring
563 .    -mat_fd_coloring_view_draw - Activates drawing of coloring
564 -    -mat_fd_coloring_view_info - Activates viewing of coloring info
565 
566     Level: intermediate
567 
568 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
569 
570 .keywords: coloring, Jacobian, finite differences
571 @*/
572 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
573 {
574   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
575   PetscErrorCode ierr;
576   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
577   PetscScalar    dx,*y,*xx,*w3_array;
578   PetscScalar    *vscale_array;
579   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
580   Vec            w1=coloring->w1,w2=coloring->w2,w3;
581   void           *fctx = coloring->fctx;
582   PetscTruth     flg;
583   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
584   Vec            x1_tmp;
585 
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
588   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
589   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
590   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
591 
592   if (coloring->usersetsrecompute) {
593     if (!coloring->recompute) {
594       *flag = SAME_PRECONDITIONER;
595       ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr);
596       PetscFunctionReturn(0);
597     } else {
598       coloring->recompute = PETSC_FALSE;
599     }
600   }
601 
602   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
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   x1_tmp = x1;
616   if (!coloring->vscale){
617     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
618   }
619 
620   /*
621     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
622     coloring->F for the coarser grids from the finest
623   */
624   if (coloring->F) {
625     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
626     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
627     if (m1 != m2) {
628       coloring->F = 0;
629       }
630     }
631 
632   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
633     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
634   }
635   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
636 
637   /* Set w1 = F(x1) */
638   if (coloring->F) {
639     w1          = coloring->F; /* use already computed value of function */
640     coloring->F = 0;
641   } else {
642     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
643     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
644     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
645   }
646 
647   if (!coloring->w3) {
648     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
649     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
650   }
651   w3 = coloring->w3;
652 
653     /* Compute all the local scale factors, including ghost points */
654   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
655   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
656   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
657   if (ctype == IS_COLORING_GHOSTED){
658     col_start = 0; col_end = N;
659   } else if (ctype == IS_COLORING_LOCAL){
660     xx = xx - start;
661     vscale_array = vscale_array - start;
662     col_start = start; col_end = N + start;
663   }
664   for (col=col_start; col<col_end; col++){
665     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
666     if (coloring->htype[0] == 'w') {
667       dx = 1.0 + unorm;
668     } else {
669       dx  = xx[col];
670     }
671     if (dx == 0.0) dx = 1.0;
672 #if !defined(PETSC_USE_COMPLEX)
673     if (dx < umin && dx >= 0.0)      dx = umin;
674     else if (dx < 0.0 && dx > -umin) dx = -umin;
675 #else
676     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
677     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
678 #endif
679     dx               *= epsilon;
680     vscale_array[col] = 1.0/dx;
681   }
682   if (ctype == IS_COLORING_LOCAL)  vscale_array = vscale_array + start;
683   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
684   if (ctype == IS_COLORING_LOCAL){
685     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
686     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
687   }
688 
689   if (coloring->vscaleforrow) {
690     vscaleforrow = coloring->vscaleforrow;
691   } else {
692     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
693   }
694 
695   /*
696     Loop over each color
697   */
698   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
699   for (k=0; k<coloring->ncolors; k++) {
700     coloring->currentcolor = k;
701     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
702     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
703     if (ctype == IS_COLORING_LOCAL) w3_array = w3_array - start;
704     /*
705       Loop over each column associated with color
706       adding the perturbation to the vector w3.
707     */
708     for (l=0; l<coloring->ncolumns[k]; l++) {
709       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
710       if (coloring->htype[0] == 'w') {
711         dx = 1.0 + unorm;
712       } else {
713         dx  = xx[col];
714       }
715       if (dx == 0.0) dx = 1.0;
716 #if !defined(PETSC_USE_COMPLEX)
717       if (dx < umin && dx >= 0.0)      dx = umin;
718       else if (dx < 0.0 && dx > -umin) dx = -umin;
719 #else
720       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
721       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
722 #endif
723       dx            *= epsilon;
724       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
725       w3_array[col] += dx;
726     }
727     if (ctype == IS_COLORING_LOCAL) w3_array = w3_array + start;
728     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
729 
730     /*
731       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
732                            w2 = F(x1 + dx) - F(x1)
733     */
734     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
735     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
736     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
737     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
738 
739     /*
740       Loop over rows of vector, putting results into Jacobian matrix
741     */
742     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
743     for (l=0; l<coloring->nrows[k]; l++) {
744       row    = coloring->rows[k][l];             /* local row index */
745       col    = coloring->columnsforrow[k][l];    /* global column index */
746       y[row] *= vscale_array[vscaleforrow[k][l]];
747       srow   = row + start;
748       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
749     }
750     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
751   } /* endof for each color */
752   if (ctype == IS_COLORING_LOCAL) xx = xx + start;
753   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
754   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
755 
756   coloring->currentcolor = -1;
757   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
759   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
760 
761   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
762   if (flg) {
763     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
764   }
765   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "MatFDColoringApplyTS"
771 /*@
772     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
773     has been created, computes the Jacobian for a function via finite differences.
774 
775    Collective on Mat, MatFDColoring, and Vec
776 
777     Input Parameters:
778 +   mat - location to store Jacobian
779 .   coloring - coloring context created with MatFDColoringCreate()
780 .   x1 - location at which Jacobian is to be computed
781 -   sctx - optional context required by function (actually a SNES context)
782 
783    Options Database Keys:
784 .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
785 
786    Level: intermediate
787 
788 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
789 
790 .keywords: coloring, Jacobian, finite differences
791 @*/
792 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
793 {
794   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
795   PetscErrorCode ierr;
796   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
797   PetscScalar    dx,*y,*xx,*w3_array;
798   PetscScalar    *vscale_array;
799   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
800   Vec            w1=coloring->w1,w2=coloring->w2,w3;
801   void           *fctx = coloring->fctx;
802   PetscTruth     flg;
803 
804   PetscFunctionBegin;
805   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
806   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
807   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
808 
809   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
810   if (!coloring->w3) {
811     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
812     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
813   }
814   w3 = coloring->w3;
815 
816   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
817   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
818   if (flg) {
819     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
820   } else {
821     PetscTruth assembled;
822     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
823     if (assembled) {
824       ierr = MatZeroEntries(J);CHKERRQ(ierr);
825     }
826   }
827 
828   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
829   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
830   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
831   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
832   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
833 
834   /*
835       Compute all the scale factors and share with other processors
836   */
837   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
838   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
839   for (k=0; k<coloring->ncolors; k++) {
840     /*
841        Loop over each column associated with color adding the
842        perturbation to the vector w3.
843     */
844     for (l=0; l<coloring->ncolumns[k]; l++) {
845       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
846       dx  = xx[col];
847       if (dx == 0.0) dx = 1.0;
848 #if !defined(PETSC_USE_COMPLEX)
849       if (dx < umin && dx >= 0.0)      dx = umin;
850       else if (dx < 0.0 && dx > -umin) dx = -umin;
851 #else
852       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
853       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
854 #endif
855       dx                *= epsilon;
856       vscale_array[col] = 1.0/dx;
857     }
858   }
859   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
860   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
861   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
862 
863   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
864   else                        vscaleforrow = coloring->columnsforrow;
865 
866   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
867   /*
868       Loop over each color
869   */
870   for (k=0; k<coloring->ncolors; k++) {
871     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
872     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
873     /*
874        Loop over each column associated with color adding the
875        perturbation to the vector w3.
876     */
877     for (l=0; l<coloring->ncolumns[k]; l++) {
878       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
879       dx  = xx[col];
880       if (dx == 0.0) dx = 1.0;
881 #if !defined(PETSC_USE_COMPLEX)
882       if (dx < umin && dx >= 0.0)      dx = umin;
883       else if (dx < 0.0 && dx > -umin) dx = -umin;
884 #else
885       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
886       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
887 #endif
888       dx            *= epsilon;
889       w3_array[col] += dx;
890     }
891     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
892 
893     /*
894        Evaluate function at x1 + dx (here dx is a vector of perturbations)
895     */
896     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
897     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
898     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
899     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
900 
901     /*
902        Loop over rows of vector, putting results into Jacobian matrix
903     */
904     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
905     for (l=0; l<coloring->nrows[k]; l++) {
906       row    = coloring->rows[k][l];
907       col    = coloring->columnsforrow[k][l];
908       y[row] *= vscale_array[vscaleforrow[k][l]];
909       srow   = row + start;
910       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
911     }
912     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
913   }
914   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
915   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
916   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
917   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
918   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
919   PetscFunctionReturn(0);
920 }
921 
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "MatFDColoringSetRecompute()"
925 /*@C
926    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
927      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
928      no additional Jacobian's will be computed (the same one will be used) until this is
929      called again.
930 
931    Collective on MatFDColoring
932 
933    Input  Parameters:
934 .  c - the coloring context
935 
936    Level: intermediate
937 
938    Notes: The MatFDColoringSetFrequency() is ignored once this is called
939 
940 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
941 
942 .keywords: Mat, finite differences, coloring
943 @*/
944 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
945 {
946   PetscFunctionBegin;
947   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
948   c->usersetsrecompute = PETSC_TRUE;
949   c->recompute         = PETSC_TRUE;
950   PetscFunctionReturn(0);
951 }
952 
953 
954