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