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