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