xref: /petsc/src/mat/matfd/fdmatrix.c (revision 0e9bae810fdaeb60e2713eaa8ddb89f42e079fd1)
1 
2 /*
3    This is where the abstract matrix operations are defined that are
4   used for finite difference computations of Jacobians using coloring.
5 */
6 
7 #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatFDColoringSetF"
11 PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
12 {
13   PetscErrorCode ierr;
14   PetscFunctionBegin;
15   if (F) {
16     ierr = VecCopy(F,fd->w1);CHKERRQ(ierr);
17     fd->fset = PETSC_TRUE;
18   } else {
19     fd->fset = PETSC_FALSE;
20   }
21   PetscFunctionReturn(0);
22 }
23 
24 #undef __FUNCT__
25 #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
26 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
27 {
28   MatFDColoring  fd = (MatFDColoring)Aa;
29   PetscErrorCode ierr;
30   PetscInt       i,j;
31   PetscReal      x,y;
32 
33   PetscFunctionBegin;
34 
35   /* loop over colors  */
36   for (i=0; i<fd->ncolors; i++) {
37     for (j=0; j<fd->nrows[i]; j++) {
38       y = fd->M - fd->rows[i][j] - fd->rstart;
39       x = fd->columnsforrow[i][j];
40       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
41     }
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "MatFDColoringView_Draw"
48 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
49 {
50   PetscErrorCode ierr;
51   PetscBool      isnull;
52   PetscDraw      draw;
53   PetscReal      xr,yr,xl,yl,h,w;
54 
55   PetscFunctionBegin;
56   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
57   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
58 
59   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
60 
61   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
62   xr += w;     yr += h;    xl = -w;     yl = -h;
63   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
64   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
65   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "MatFDColoringView"
71 /*@C
72    MatFDColoringView - Views a finite difference coloring context.
73 
74    Collective on MatFDColoring
75 
76    Input  Parameters:
77 +  c - the coloring context
78 -  viewer - visualization context
79 
80    Level: intermediate
81 
82    Notes:
83    The available visualization contexts include
84 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
85 .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
86         output where only the first processor opens
87         the file.  All other processors send their
88         data to the first processor to print.
89 -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
90 
91    Notes:
92      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
93    involves more than 33 then some seemingly identical colors are displayed making it look
94    like an illegal coloring. This is just a graphical artifact.
95 
96 .seealso: MatFDColoringCreate()
97 
98 .keywords: Mat, finite differences, coloring, view
99 @*/
100 PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
101 {
102   PetscErrorCode    ierr;
103   PetscInt          i,j;
104   PetscBool         isdraw,iascii;
105   PetscViewerFormat format;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
109   if (!viewer) {
110     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
111   }
112   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
113   PetscCheckSameComm(c,1,viewer,2);
114 
115   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
116   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
117   if (isdraw) {
118     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
119   } else if (iascii) {
120     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer,"MatFDColoring Object");CHKERRQ(ierr);
121     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
122     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
123     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
124 
125     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
126     if (format != PETSC_VIEWER_ASCII_INFO) {
127       for (i=0; i<c->ncolors; i++) {
128         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
129         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
130         for (j=0; j<c->ncolumns[i]; j++) {
131           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
132         }
133         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
134         for (j=0; j<c->nrows[i]; j++) {
135           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
136         }
137       }
138     }
139     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
140   } else {
141     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
142   }
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "MatFDColoringSetParameters"
148 /*@
149    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
150    a Jacobian matrix using finite differences.
151 
152    Logically Collective on MatFDColoring
153 
154    The Jacobian is estimated with the differencing approximation
155 .vb
156        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
157        h = error_rel*u[i]                 if  abs(u[i]) > umin
158          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
159        dx_{i} = (0, ... 1, .... 0)
160 .ve
161 
162    Input Parameters:
163 +  coloring - the coloring context
164 .  error_rel - relative error
165 -  umin - minimum allowable u-value magnitude
166 
167    Level: advanced
168 
169 .keywords: Mat, finite differences, coloring, set, parameters
170 
171 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
172 
173 @*/
174 PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
175 {
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
178   PetscValidLogicalCollectiveReal(matfd,error,2);
179   PetscValidLogicalCollectiveReal(matfd,umin,3);
180   if (error != PETSC_DEFAULT) matfd->error_rel = error;
181   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
182   PetscFunctionReturn(0);
183 }
184 
185 
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "MatFDColoringGetFunction"
189 /*@C
190    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
191 
192    Not Collective
193 
194    Input Parameters:
195 .  coloring - the coloring context
196 
197    Output Parameters:
198 +  f - the function
199 -  fctx - the optional user-defined function context
200 
201    Level: intermediate
202 
203 .keywords: Mat, Jacobian, finite differences, set, function
204 
205 .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
206 
207 @*/
208 PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
209 {
210   PetscFunctionBegin;
211   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
212   if (f) *f = matfd->f;
213   if (fctx) *fctx = matfd->fctx;
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "MatFDColoringSetFunction"
219 /*@C
220    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
221 
222    Logically Collective on MatFDColoring
223 
224    Input Parameters:
225 +  coloring - the coloring context
226 .  f - the function
227 -  fctx - the optional user-defined function context
228 
229    Calling sequence of (*f) function:
230     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
231     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
232 
233    Level: advanced
234 
235    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
236      SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by
237      calling MatFDColoringApply()
238 
239    Fortran Notes:
240     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
241   be used without SNES or within the SNES solvers.
242 
243 .keywords: Mat, Jacobian, finite differences, set, function
244 
245 .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
246 
247 @*/
248 PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
249 {
250   PetscFunctionBegin;
251   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
252   matfd->f    = f;
253   matfd->fctx = fctx;
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "MatFDColoringSetFromOptions"
259 /*@
260    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
261    the options database.
262 
263    Collective on MatFDColoring
264 
265    The Jacobian, F'(u), is estimated with the differencing approximation
266 .vb
267        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
268        h = error_rel*u[i]                 if  abs(u[i]) > umin
269          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
270        dx_{i} = (0, ... 1, .... 0)
271 .ve
272 
273    Input Parameter:
274 .  coloring - the coloring context
275 
276    Options Database Keys:
277 +  -mat_fd_coloring_err <err> - Sets <err> (square root
278            of relative error in the function)
279 .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
280 .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
281 .  -mat_fd_coloring_view - Activates basic viewing
282 .  -mat_fd_coloring_view_info - Activates viewing info
283 -  -mat_fd_coloring_view_draw - Activates drawing
284 
285     Level: intermediate
286 
287 .keywords: Mat, finite differences, parameters
288 
289 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
290 
291 @*/
292 PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
293 {
294   PetscErrorCode ierr;
295   PetscBool      flg;
296   char           value[3];
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
300 
301   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
302     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
303     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
304     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
305     if (flg) {
306       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
307       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
308       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
309     }
310     /* not used here; but so they are presented in the GUI */
311     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
312     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
313     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
314 
315     /* process any options handlers added with PetscObjectAddOptionsHandler() */
316     ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
317   PetscOptionsEnd();CHKERRQ(ierr);
318   PetscFunctionReturn(0);
319 }
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "MatFDColoringView_Private"
323 PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
324 {
325   PetscErrorCode ierr;
326   PetscBool      flg = PETSC_FALSE;
327   PetscViewer    viewer;
328 
329   PetscFunctionBegin;
330   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
331   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
332   if (flg) {
333     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
334   }
335   flg  = PETSC_FALSE;
336   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
337   if (flg) {
338     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
339     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
340     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
341   }
342   flg  = PETSC_FALSE;
343   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
344   if (flg) {
345     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
346     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
347   }
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "MatFDColoringCreate"
353 /*@
354    MatFDColoringCreate - Creates a matrix coloring context for finite difference
355    computation of Jacobians.
356 
357    Collective on Mat
358 
359    Input Parameters:
360 +  mat - the matrix containing the nonzero structure of the Jacobian
361 -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
362 
363     Output Parameter:
364 .   color - the new coloring context
365 
366     Level: intermediate
367 
368 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
369           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
370           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
371 @*/
372 PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
373 {
374   MatFDColoring  c;
375   MPI_Comm       comm;
376   PetscErrorCode ierr;
377   PetscInt       M,N;
378 
379   PetscFunctionBegin;
380   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
381   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
382   if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices");
383   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
384   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
385 
386   c->ctype = iscoloring->ctype;
387 
388   if (mat->ops->fdcoloringcreate) {
389     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
390   } else SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
391 
392   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
393   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
394   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
395   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
396 
397   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
398   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
399   c->currentcolor      = -1;
400   c->htype             = "wp";
401   c->fset              = PETSC_FALSE;
402 
403   *color = c;
404   ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
405   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "MatFDColoringDestroy"
411 /*@
412     MatFDColoringDestroy - Destroys a matrix coloring context that was created
413     via MatFDColoringCreate().
414 
415     Collective on MatFDColoring
416 
417     Input Parameter:
418 .   c - coloring context
419 
420     Level: intermediate
421 
422 .seealso: MatFDColoringCreate()
423 @*/
424 PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
425 {
426   PetscErrorCode ierr;
427   PetscInt       i;
428 
429   PetscFunctionBegin;
430   if (!*c) PetscFunctionReturn(0);
431   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
432 
433   for (i=0; i<(*c)->ncolors; i++) {
434     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
435     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
436     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
437     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
438   }
439   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
440   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
441   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
442   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
443   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
444   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
445   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
446   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
447   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
448   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
449   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
450   PetscFunctionReturn(0);
451 }
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
455 /*@C
456     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
457       that are currently being perturbed.
458 
459     Not Collective
460 
461     Input Parameters:
462 .   coloring - coloring context created with MatFDColoringCreate()
463 
464     Output Parameters:
465 +   n - the number of local columns being perturbed
466 -   cols - the column indices, in global numbering
467 
468    Level: intermediate
469 
470 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
471 
472 .keywords: coloring, Jacobian, finite differences
473 @*/
474 PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
475 {
476   PetscFunctionBegin;
477   if (coloring->currentcolor >= 0) {
478     *n    = coloring->ncolumns[coloring->currentcolor];
479     *cols = coloring->columns[coloring->currentcolor];
480   } else {
481     *n = 0;
482   }
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "MatFDColoringApply"
488 /*@
489     MatFDColoringApply - Given a matrix for which a MatFDColoring context
490     has been created, computes the Jacobian for a function via finite differences.
491 
492     Collective on MatFDColoring
493 
494     Input Parameters:
495 +   mat - location to store Jacobian
496 .   coloring - coloring context created with MatFDColoringCreate()
497 .   x1 - location at which Jacobian is to be computed
498 -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
499 
500     Options Database Keys:
501 +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
502 .    -mat_fd_coloring_view - Activates basic viewing or coloring
503 .    -mat_fd_coloring_view_draw - Activates drawing of coloring
504 -    -mat_fd_coloring_view_info - Activates viewing of coloring info
505 
506     Level: intermediate
507 
508 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
509 
510 .keywords: coloring, Jacobian, finite differences
511 @*/
512 PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
513 {
514   PetscErrorCode ierr;
515 
516   PetscFunctionBegin;
517   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
518   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
519   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
520   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
521   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
522   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 
526 #undef __FUNCT__
527 #define __FUNCT__ "MatFDColoringApply_AIJ"
528 PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
529 {
530   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
531   PetscErrorCode ierr;
532   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow;
533   PetscScalar    dx,*y,*xx,*w3_array;
534   PetscScalar    *vscale_array;
535   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
536   Vec            w1=coloring->w1,w2=coloring->w2,w3;
537   void           *fctx = coloring->fctx;
538   PetscBool      flg = PETSC_FALSE;
539   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
540   Vec            x1_tmp;
541 
542   PetscFunctionBegin;
543   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
544   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
545   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
546   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
547 
548   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
549   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
550   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
551   if (flg) {
552     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
553   } else {
554     PetscBool  assembled;
555     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
556     if (assembled) {
557       ierr = MatZeroEntries(J);CHKERRQ(ierr);
558     }
559   }
560 
561   x1_tmp = x1;
562   if (!coloring->vscale){
563     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
564   }
565 
566   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
567     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
568   }
569   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
570 
571   /* Set w1 = F(x1) */
572   if (!coloring->fset) {
573     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
574     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
575     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
576   } else {
577     coloring->fset = PETSC_FALSE;
578 }
579 
580   if (!coloring->w3) {
581     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
582     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
583   }
584   w3 = coloring->w3;
585 
586     /* Compute all the local scale factors, including ghost points */
587   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
588   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
589   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
590   if (ctype == IS_COLORING_GHOSTED){
591     col_start = 0; col_end = N;
592   } else if (ctype == IS_COLORING_GLOBAL){
593     xx = xx - start;
594     vscale_array = vscale_array - start;
595     col_start = start; col_end = N + start;
596   }
597   for (col=col_start; col<col_end; col++){
598     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
599     if (coloring->htype[0] == 'w') {
600       dx = 1.0 + unorm;
601     } else {
602       dx  = xx[col];
603     }
604     if (dx == (PetscScalar)0.0) dx = 1.0;
605 #if !defined(PETSC_USE_COMPLEX)
606     if (dx < umin && dx >= 0.0)      dx = umin;
607     else if (dx < 0.0 && dx > -umin) dx = -umin;
608 #else
609     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
610     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
611 #endif
612     dx               *= epsilon;
613     vscale_array[col] = (PetscScalar)1.0/dx;
614   }
615   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
616   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
617   if (ctype == IS_COLORING_GLOBAL){
618     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620   }
621 
622   if (coloring->vscaleforrow) {
623     vscaleforrow = coloring->vscaleforrow;
624   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
625 
626   /*
627     Loop over each color
628   */
629   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
630   for (k=0; k<coloring->ncolors; k++) {
631     coloring->currentcolor = k;
632     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
633     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
634     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
635     /*
636       Loop over each column associated with color
637       adding the perturbation to the vector w3.
638     */
639     for (l=0; l<coloring->ncolumns[k]; l++) {
640       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
641       if (coloring->htype[0] == 'w') {
642         dx = 1.0 + unorm;
643       } else {
644         dx  = xx[col];
645       }
646       if (dx == (PetscScalar)0.0) dx = 1.0;
647 #if !defined(PETSC_USE_COMPLEX)
648       if (dx < umin && dx >= 0.0)      dx = umin;
649       else if (dx < 0.0 && dx > -umin) dx = -umin;
650 #else
651       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
652       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
653 #endif
654       dx            *= epsilon;
655       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
656       w3_array[col] += dx;
657     }
658     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
659     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
660 
661     /*
662       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
663                            w2 = F(x1 + dx) - F(x1)
664     */
665     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
666     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
667     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
668     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
669 
670     /*
671       Loop over rows of vector, putting results into Jacobian matrix
672     */
673     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
674     for (l=0; l<coloring->nrows[k]; l++) {
675       row    = coloring->rows[k][l];             /* local row index */
676       col    = coloring->columnsforrow[k][l];    /* global column index */
677       y[row] *= vscale_array[vscaleforrow[k][l]];
678       srow   = row + start;
679       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
680     }
681     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
682   } /* endof for each color */
683   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
684   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
685   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
686 
687   coloring->currentcolor = -1;
688   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
690   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
691 
692   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
693   PetscFunctionReturn(0);
694 }
695