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