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