xref: /petsc/src/mat/matfd/fdmatrix.c (revision 4e7234bf172ce7019c4948834d73bd3002d026ba)
173f4d377SMatthew Knepley /*$Id: fdmatrix.c,v 1.92 2001/08/21 21:03:06 bsmith Exp $*/
2bbf0e169SBarry Smith 
3bbf0e169SBarry Smith /*
4639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
5639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
6bbf0e169SBarry Smith */
7bbf0e169SBarry Smith 
8e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
9bbf0e169SBarry Smith 
108ba1e511SMatthew Knepley /* Logging support */
118ba1e511SMatthew Knepley int MAT_FDCOLORING_COOKIE;
128ba1e511SMatthew Knepley 
134a2ae208SSatish Balay #undef __FUNCT__
143a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
153a7fca6bSBarry Smith int MatFDColoringSetF(MatFDColoring fd,Vec F)
163a7fca6bSBarry Smith {
173a7fca6bSBarry Smith   PetscFunctionBegin;
183a7fca6bSBarry Smith   fd->F = F;
193a7fca6bSBarry Smith   PetscFunctionReturn(0);
203a7fca6bSBarry Smith }
213a7fca6bSBarry Smith 
223a7fca6bSBarry Smith #undef __FUNCT__
234a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
24b0a32e0cSBarry Smith static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
259194eea9SBarry Smith {
269194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
279194eea9SBarry Smith   int           ierr,i,j;
289194eea9SBarry Smith   PetscReal     x,y;
299194eea9SBarry Smith 
309194eea9SBarry Smith   PetscFunctionBegin;
319194eea9SBarry Smith 
329194eea9SBarry Smith   /* loop over colors  */
339194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
349194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
359194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
369194eea9SBarry Smith       x = fd->columnsforrow[i][j];
37b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
389194eea9SBarry Smith     }
399194eea9SBarry Smith   }
409194eea9SBarry Smith   PetscFunctionReturn(0);
419194eea9SBarry Smith }
429194eea9SBarry Smith 
434a2ae208SSatish Balay #undef __FUNCT__
444a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
45b0a32e0cSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
46005c665bSBarry Smith {
479194eea9SBarry Smith   int         ierr;
48005c665bSBarry Smith   PetscTruth  isnull;
49b0a32e0cSBarry Smith   PetscDraw   draw;
5054d96fa3SBarry Smith   PetscReal   xr,yr,xl,yl,h,w;
51005c665bSBarry Smith 
523a40ed3dSBarry Smith   PetscFunctionBegin;
53b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
54b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
559194eea9SBarry Smith 
569194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
57005c665bSBarry Smith 
58005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
59005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
60b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
61b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
629194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
633a40ed3dSBarry Smith   PetscFunctionReturn(0);
64005c665bSBarry Smith }
65005c665bSBarry Smith 
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
68bbf0e169SBarry Smith /*@C
69639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
70bbf0e169SBarry Smith 
714c49b128SBarry Smith    Collective on MatFDColoring
72fee21e36SBarry Smith 
73ef5ee4d1SLois Curfman McInnes    Input  Parameters:
74ef5ee4d1SLois Curfman McInnes +  c - the coloring context
75ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
76ef5ee4d1SLois Curfman McInnes 
7715091d37SBarry Smith    Level: intermediate
7815091d37SBarry Smith 
79b4fc646aSLois Curfman McInnes    Notes:
80b4fc646aSLois Curfman McInnes    The available visualization contexts include
81b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
82b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
83ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
84ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
85ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
86b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
87bbf0e169SBarry Smith 
889194eea9SBarry Smith    Notes:
899194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
90b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
919194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
929194eea9SBarry Smith 
93639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
94005c665bSBarry Smith 
95b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
96bbf0e169SBarry Smith @*/
97b0a32e0cSBarry Smith int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
98bbf0e169SBarry Smith {
99fb9695e5SSatish Balay   int               i,j,ierr;
1006831982aSBarry Smith   PetscTruth        isdraw,isascii;
101f3ef73ceSBarry Smith   PetscViewerFormat format;
102bbf0e169SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
104b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
105b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
106b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
1076831982aSBarry Smith   PetscCheckSameComm(c,viewer);
108bbf0e169SBarry Smith 
109fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
110b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
1110f5bd95cSBarry Smith   if (isdraw) {
112b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1130f5bd95cSBarry Smith   } else if (isascii) {
114b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
115b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
116b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
117b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
118ae09f205SBarry Smith 
119b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
120fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
121b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
122b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
123b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
124b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
125b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
126639f9d9dSBarry Smith         }
127b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
128b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
129b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         }
131bbf0e169SBarry Smith       }
132bbf0e169SBarry Smith     }
133b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1345cd90555SBarry Smith   } else {
13529bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
136bbf0e169SBarry Smith   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
138639f9d9dSBarry Smith }
139639f9d9dSBarry Smith 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
142639f9d9dSBarry Smith /*@
143b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
145639f9d9dSBarry Smith 
146ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
147ef5ee4d1SLois Curfman McInnes 
148ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
149ef5ee4d1SLois Curfman McInnes .vb
15065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
152f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
153ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
154ef5ee4d1SLois Curfman McInnes .ve
155639f9d9dSBarry Smith 
156639f9d9dSBarry Smith    Input Parameters:
157ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
158639f9d9dSBarry Smith .  error_rel - relative error
159f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
160fee21e36SBarry Smith 
16115091d37SBarry Smith    Level: advanced
16215091d37SBarry Smith 
163b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
164b4fc646aSLois Curfman McInnes 
165b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
166639f9d9dSBarry Smith @*/
167329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
168639f9d9dSBarry Smith {
1693a40ed3dSBarry Smith   PetscFunctionBegin;
170639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
171639f9d9dSBarry Smith 
172639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
173639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175639f9d9dSBarry Smith }
176639f9d9dSBarry Smith 
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency"
179005c665bSBarry Smith /*@
180e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
181e0907662SLois Curfman McInnes    matrices.
182005c665bSBarry Smith 
183fee21e36SBarry Smith    Collective on MatFDColoring
184fee21e36SBarry Smith 
185ef5ee4d1SLois Curfman McInnes    Input Parameters:
186ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
187ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
188ef5ee4d1SLois Curfman McInnes 
18915091d37SBarry Smith    Options Database Keys:
19015091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
19115091d37SBarry Smith 
19215091d37SBarry Smith    Level: advanced
19315091d37SBarry Smith 
194e0907662SLois Curfman McInnes    Notes:
195e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
196e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
197e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
198e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
199e0907662SLois Curfman McInnes 
200b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
201ef5ee4d1SLois Curfman McInnes 
20240964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
203005c665bSBarry Smith @*/
204005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
205005c665bSBarry Smith {
2063a40ed3dSBarry Smith   PetscFunctionBegin;
207005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
208005c665bSBarry Smith 
209005c665bSBarry Smith   matfd->freq = freq;
2103a40ed3dSBarry Smith   PetscFunctionReturn(0);
211005c665bSBarry Smith }
212005c665bSBarry Smith 
2134a2ae208SSatish Balay #undef __FUNCT__
2144a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency"
215ff0cfa39SBarry Smith /*@
216ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
217ff0cfa39SBarry Smith    matrices.
218ff0cfa39SBarry Smith 
219ef5ee4d1SLois Curfman McInnes    Not Collective
220ef5ee4d1SLois Curfman McInnes 
221ff0cfa39SBarry Smith    Input Parameters:
222ff0cfa39SBarry Smith .  coloring - the coloring context
223ff0cfa39SBarry Smith 
224ff0cfa39SBarry Smith    Output Parameters:
225ff0cfa39SBarry Smith .  freq - frequency (default is 1)
226ff0cfa39SBarry Smith 
22715091d37SBarry Smith    Options Database Keys:
22815091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
22915091d37SBarry Smith 
23015091d37SBarry Smith    Level: advanced
23115091d37SBarry Smith 
232ff0cfa39SBarry Smith    Notes:
233ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
234ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
235ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
236ff0cfa39SBarry Smith    <freq> nonlinear iterations.
237ff0cfa39SBarry Smith 
238ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
239ef5ee4d1SLois Curfman McInnes 
240ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
241ff0cfa39SBarry Smith @*/
242ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
243ff0cfa39SBarry Smith {
2443a40ed3dSBarry Smith   PetscFunctionBegin;
245ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
246ff0cfa39SBarry Smith 
247ff0cfa39SBarry Smith   *freq = matfd->freq;
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
249ff0cfa39SBarry Smith }
250ff0cfa39SBarry Smith 
2514a2ae208SSatish Balay #undef __FUNCT__
2524a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
253d64ed03dSBarry Smith /*@C
254005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
255005c665bSBarry Smith 
256fee21e36SBarry Smith    Collective on MatFDColoring
257fee21e36SBarry Smith 
258ef5ee4d1SLois Curfman McInnes    Input Parameters:
259ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
260ef5ee4d1SLois Curfman McInnes .  f - the function
261ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
262ef5ee4d1SLois Curfman McInnes 
26315091d37SBarry Smith    Level: intermediate
26415091d37SBarry Smith 
265f881d145SBarry Smith    Notes:
266f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
267f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
268f881d145SBarry Smith   with the TS solvers.
269f881d145SBarry Smith 
270b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
271005c665bSBarry Smith @*/
272840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
273005c665bSBarry Smith {
2743a40ed3dSBarry Smith   PetscFunctionBegin;
275005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
276005c665bSBarry Smith 
277005c665bSBarry Smith   matfd->f    = f;
278005c665bSBarry Smith   matfd->fctx = fctx;
279005c665bSBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
281005c665bSBarry Smith }
282005c665bSBarry Smith 
2834a2ae208SSatish Balay #undef __FUNCT__
2844a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
285639f9d9dSBarry Smith /*@
286b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
287639f9d9dSBarry Smith    the options database.
288639f9d9dSBarry Smith 
289fee21e36SBarry Smith    Collective on MatFDColoring
290fee21e36SBarry Smith 
29165f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
292ef5ee4d1SLois Curfman McInnes .vb
29365f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
294f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
295f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
296ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
297ef5ee4d1SLois Curfman McInnes .ve
298ef5ee4d1SLois Curfman McInnes 
299ef5ee4d1SLois Curfman McInnes    Input Parameter:
300ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
301ef5ee4d1SLois Curfman McInnes 
302b4fc646aSLois Curfman McInnes    Options Database Keys:
303d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
304ef5ee4d1SLois Curfman McInnes            of relative error in the function)
305f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
306ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
307ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
308ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
309ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
310639f9d9dSBarry Smith 
31115091d37SBarry Smith     Level: intermediate
31215091d37SBarry Smith 
313b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
314d1c39f55SBarry Smith 
315d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
316d1c39f55SBarry Smith 
317639f9d9dSBarry Smith @*/
318639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
319639f9d9dSBarry Smith {
320186905e3SBarry Smith   int        ierr;
3213a40ed3dSBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
323639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
324639f9d9dSBarry Smith 
325b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
32687828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
32787828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
328b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
329186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
330b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
331b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
332b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
333b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3343a40ed3dSBarry Smith   PetscFunctionReturn(0);
335005c665bSBarry Smith }
336005c665bSBarry Smith 
3374a2ae208SSatish Balay #undef __FUNCT__
3384a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
339005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
340005c665bSBarry Smith {
341f1af5d2fSBarry Smith   int        ierr;
342f1af5d2fSBarry Smith   PetscTruth flg;
343005c665bSBarry Smith 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
345b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
346005c665bSBarry Smith   if (flg) {
347b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
348005c665bSBarry Smith   }
349b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
350ae09f205SBarry Smith   if (flg) {
351fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
352b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
353b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
354ae09f205SBarry Smith   }
355b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
356005c665bSBarry Smith   if (flg) {
357b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
358b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
359005c665bSBarry Smith   }
3603a40ed3dSBarry Smith   PetscFunctionReturn(0);
361bbf0e169SBarry Smith }
362bbf0e169SBarry Smith 
3634a2ae208SSatish Balay #undef __FUNCT__
3644a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
365bbf0e169SBarry Smith /*@C
366639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
367639f9d9dSBarry Smith    computation of Jacobians.
368bbf0e169SBarry Smith 
369ef5ee4d1SLois Curfman McInnes    Collective on Mat
370ef5ee4d1SLois Curfman McInnes 
371639f9d9dSBarry Smith    Input Parameters:
372ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
373ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
374bbf0e169SBarry Smith 
375bbf0e169SBarry Smith     Output Parameter:
376639f9d9dSBarry Smith .   color - the new coloring context
377bbf0e169SBarry Smith 
37815091d37SBarry Smith     Level: intermediate
37915091d37SBarry Smith 
3806831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
381d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
382d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
383d1c39f55SBarry Smith           MatFDColoringSetParameters()
384bbf0e169SBarry Smith @*/
385639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
386bbf0e169SBarry Smith {
387639f9d9dSBarry Smith   MatFDColoring c;
388639f9d9dSBarry Smith   MPI_Comm      comm;
389639f9d9dSBarry Smith   int           ierr,M,N;
390639f9d9dSBarry Smith 
3913a40ed3dSBarry Smith   PetscFunctionBegin;
392d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
393639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
39429bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
395639f9d9dSBarry Smith 
396f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3979194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
398b0a32e0cSBarry Smith   PetscLogObjectCreate(c);
399639f9d9dSBarry Smith 
400f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
401f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
402639f9d9dSBarry Smith   } else {
40329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
404639f9d9dSBarry Smith   }
405639f9d9dSBarry Smith 
40677d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
40777d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
408005c665bSBarry Smith   c->freq              = 1;
40940964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
41040964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
41149b058dcSBarry Smith   c->currentcolor      = -1;
412639f9d9dSBarry Smith 
413639f9d9dSBarry Smith   *color = c;
414d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4153a40ed3dSBarry Smith   PetscFunctionReturn(0);
416bbf0e169SBarry Smith }
417bbf0e169SBarry Smith 
4184a2ae208SSatish Balay #undef __FUNCT__
4194a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
420bbf0e169SBarry Smith /*@C
421639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
422639f9d9dSBarry Smith     via MatFDColoringCreate().
423bbf0e169SBarry Smith 
424ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
425ef5ee4d1SLois Curfman McInnes 
426b4fc646aSLois Curfman McInnes     Input Parameter:
427639f9d9dSBarry Smith .   c - coloring context
428bbf0e169SBarry Smith 
42915091d37SBarry Smith     Level: intermediate
43015091d37SBarry Smith 
431639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
432bbf0e169SBarry Smith @*/
433639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
434bbf0e169SBarry Smith {
435263760aaSBarry Smith   int i,ierr;
436bbf0e169SBarry Smith 
4373a40ed3dSBarry Smith   PetscFunctionBegin;
4383a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
439d4bb536fSBarry Smith 
440639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
441606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
442606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
443606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
44430b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
445bbf0e169SBarry Smith   }
446606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
447606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
448606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
449606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
450606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
45130b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4524f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
453005c665bSBarry Smith   if (c->w1) {
454005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
455005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
456005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
457005c665bSBarry Smith   }
458b0a32e0cSBarry Smith   PetscLogObjectDestroy(c);
459639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4603a40ed3dSBarry Smith   PetscFunctionReturn(0);
461bbf0e169SBarry Smith }
46243a90d84SBarry Smith 
4634a2ae208SSatish Balay #undef __FUNCT__
46449b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
46549b058dcSBarry Smith /*@C
46649b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
46749b058dcSBarry Smith       that are currently being perturbed.
46849b058dcSBarry Smith 
46949b058dcSBarry Smith     Not Collective
47049b058dcSBarry Smith 
47149b058dcSBarry Smith     Input Parameters:
47249b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
47349b058dcSBarry Smith 
47449b058dcSBarry Smith     Output Parameters:
47549b058dcSBarry Smith +   n - the number of local columns being perturbed
47649b058dcSBarry Smith -   cols - the column indices, in global numbering
47749b058dcSBarry Smith 
47849b058dcSBarry Smith    Level: intermediate
47949b058dcSBarry Smith 
48049b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
48149b058dcSBarry Smith 
48249b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
48349b058dcSBarry Smith @*/
484*4e7234bfSBarry Smith EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int *cols[])
48549b058dcSBarry Smith {
48649b058dcSBarry Smith   PetscFunctionBegin;
48749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
48849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
48949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
49049b058dcSBarry Smith   } else {
49149b058dcSBarry Smith     *n = 0;
49249b058dcSBarry Smith   }
49349b058dcSBarry Smith   PetscFunctionReturn(0);
49449b058dcSBarry Smith }
49549b058dcSBarry Smith 
49649b058dcSBarry Smith #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
49843a90d84SBarry Smith /*@
499e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
500e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
50143a90d84SBarry Smith 
502fee21e36SBarry Smith     Collective on MatFDColoring
503fee21e36SBarry Smith 
504ef5ee4d1SLois Curfman McInnes     Input Parameters:
505ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
506ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
507ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
508ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
509ef5ee4d1SLois Curfman McInnes 
5108bba8e72SBarry Smith     Options Database Keys:
511b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
512b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
513b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
514b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5158bba8e72SBarry Smith 
51615091d37SBarry Smith     Level: intermediate
51715091d37SBarry Smith 
51843a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
51943a90d84SBarry Smith 
52043a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
52143a90d84SBarry Smith @*/
522005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
52343a90d84SBarry Smith {
5245904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
5253a7fca6bSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
526ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
527ea709b57SSatish Balay   PetscScalar   *vscale_array;
528329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
529005c665bSBarry Smith   Vec           w1,w2,w3;
530005c665bSBarry Smith   void          *fctx = coloring->fctx;
531f1af5d2fSBarry Smith   PetscTruth    flg;
532005c665bSBarry Smith 
533a2e34c3dSBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
535e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
536e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
537e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
538e0907662SLois Curfman McInnes 
53940964ad5SBarry Smith   if (coloring->usersetsrecompute) {
54040964ad5SBarry Smith     if (!coloring->recompute) {
54140964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
54276d8deaeSBarry Smith       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
54340964ad5SBarry Smith       PetscFunctionReturn(0);
54440964ad5SBarry Smith     } else {
54540964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
54640964ad5SBarry Smith     }
54740964ad5SBarry Smith   }
54840964ad5SBarry Smith 
549d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
55076d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
55176d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
55276d8deaeSBarry Smith   } else {
55376d8deaeSBarry Smith 
554005c665bSBarry Smith     if (!coloring->w1) {
555005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
556b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w1);
557005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
558b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w2);
559005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
560b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w3);
561005c665bSBarry Smith     }
562005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
56343a90d84SBarry Smith 
5644c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
565b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
566f1af5d2fSBarry Smith     if (flg) {
567b0a32e0cSBarry Smith       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
568e0907662SLois Curfman McInnes     } else {
56943a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
570e0907662SLois Curfman McInnes     }
57143a90d84SBarry Smith 
57243a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
57343a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
574be2fbe1fSBarry Smith 
5753a7fca6bSBarry Smith     /*
5763a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5773a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5783a7fca6bSBarry Smith     */
5793a7fca6bSBarry Smith     if (coloring->F) {
5803a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5813a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5823a7fca6bSBarry Smith       if (m1 != m2) {
5833a7fca6bSBarry Smith 	coloring->F = 0;
5843a7fca6bSBarry Smith       }
5853a7fca6bSBarry Smith     }
5863a7fca6bSBarry Smith 
587be2fbe1fSBarry Smith     if (coloring->F) {
588be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
589be2fbe1fSBarry Smith       coloring->F = 0;
590be2fbe1fSBarry Smith     } else {
59166f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
59243a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
59366f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
594be2fbe1fSBarry Smith     }
59543a90d84SBarry Smith 
59643a90d84SBarry Smith     /*
5979782cf98SBarry Smith        Compute all the scale factors and share with other processors
5989782cf98SBarry Smith     */
5999782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
6004f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
6019782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
6029782cf98SBarry Smith       /*
6039782cf98SBarry Smith 	Loop over each column associated with color adding the
6049782cf98SBarry Smith 	perturbation to the vector w3.
6059782cf98SBarry Smith       */
6069782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
6079782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
6089782cf98SBarry Smith 	dx  = xx[col];
6099782cf98SBarry Smith 	if (dx == 0.0) dx = 1.0;
6109782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6119782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
6129782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
6139782cf98SBarry Smith #else
6149782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6159782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6169782cf98SBarry Smith #endif
6179782cf98SBarry Smith 	dx                *= epsilon;
61830b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
6199782cf98SBarry Smith       }
6209782cf98SBarry Smith     }
621a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
62230b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62330b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6249782cf98SBarry Smith 
625ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
626ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
627b0a32e0cSBarry Smith 
62830b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
62930b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
63030b34957SBarry Smith 
63130b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6329782cf98SBarry Smith     /*
63343a90d84SBarry Smith       Loop over each color
63443a90d84SBarry Smith     */
63543a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
63649b058dcSBarry Smith       coloring->currentcolor = k;
63743a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
63842460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
63943a90d84SBarry Smith       /*
64043a90d84SBarry Smith 	Loop over each column associated with color adding the
64143a90d84SBarry Smith 	perturbation to the vector w3.
64243a90d84SBarry Smith       */
64343a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
64443a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
64542460c72SBarry Smith 	dx  = xx[col];
646ae09f205SBarry Smith 	if (dx == 0.0) dx = 1.0;
647aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
648ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
649ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
65043a90d84SBarry Smith #else
651329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
652329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
65343a90d84SBarry Smith #endif
65443a90d84SBarry Smith 	dx            *= epsilon;
65529bbc08cSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
65642460c72SBarry Smith 	w3_array[col] += dx;
65743a90d84SBarry Smith       }
65842460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6593b28642cSBarry Smith 
66043a90d84SBarry Smith       /*
661e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
66243a90d84SBarry Smith       */
663a2e34c3dSBarry Smith 
66466f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66543a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
66666f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66743a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6689782cf98SBarry Smith 
66943a90d84SBarry Smith       /*
670e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
67143a90d84SBarry Smith       */
6723b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
67343a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
67443a90d84SBarry Smith 	row    = coloring->rows[k][l];
67543a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6765904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
67743a90d84SBarry Smith 	srow   = row + start;
67843a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
67943a90d84SBarry Smith       }
6803b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
68143a90d84SBarry Smith     }
68249b058dcSBarry Smith     coloring->currentcolor = -1;
68330b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
68442460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
68543a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68643a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68776d8deaeSBarry Smith   }
688d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
689a2e34c3dSBarry Smith 
690b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
691a2e34c3dSBarry Smith   if (flg) {
692a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
693a2e34c3dSBarry Smith   }
694b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
695b9722508SLois Curfman McInnes 
6963a40ed3dSBarry Smith   PetscFunctionReturn(0);
69743a90d84SBarry Smith }
698840b8ebdSBarry Smith 
6994a2ae208SSatish Balay #undef __FUNCT__
7004a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
701840b8ebdSBarry Smith /*@
702840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
703840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
704840b8ebdSBarry Smith 
705fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
706fee21e36SBarry Smith 
707ef5ee4d1SLois Curfman McInnes     Input Parameters:
7083b28642cSBarry Smith +   mat - location to store Jacobian
709ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
710ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
711ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
712ef5ee4d1SLois Curfman McInnes 
713840b8ebdSBarry Smith    Options Database Keys:
714ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
715840b8ebdSBarry Smith 
71615091d37SBarry Smith    Level: intermediate
71715091d37SBarry Smith 
718840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
719840b8ebdSBarry Smith 
720840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
721840b8ebdSBarry Smith @*/
722329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
723840b8ebdSBarry Smith {
724329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
7255904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
726ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
727ea709b57SSatish Balay   PetscScalar   *vscale_array;
728329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
729840b8ebdSBarry Smith   Vec           w1,w2,w3;
730840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
731f1af5d2fSBarry Smith   PetscTruth    flg;
732840b8ebdSBarry Smith 
7333a40ed3dSBarry Smith   PetscFunctionBegin;
734840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
735840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
736840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
737840b8ebdSBarry Smith 
738d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
739840b8ebdSBarry Smith   if (!coloring->w1) {
740840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
741b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w1);
742840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
743b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w2);
744840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
745b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w3);
746840b8ebdSBarry Smith   }
747840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
748840b8ebdSBarry Smith 
7495904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
750b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
751f1af5d2fSBarry Smith   if (flg) {
752b0a32e0cSBarry Smith     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
753840b8ebdSBarry Smith   } else {
754840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
755840b8ebdSBarry Smith   }
756840b8ebdSBarry Smith 
757840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
758840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
75966f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
760840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
76166f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
762840b8ebdSBarry Smith 
763840b8ebdSBarry Smith   /*
7645904e1b1SBarry Smith       Compute all the scale factors and share with other processors
765840b8ebdSBarry Smith   */
7665904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7675904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
768840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
769840b8ebdSBarry Smith     /*
770840b8ebdSBarry Smith        Loop over each column associated with color adding the
771840b8ebdSBarry Smith        perturbation to the vector w3.
772840b8ebdSBarry Smith     */
773840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
774840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7755904e1b1SBarry Smith       dx  = xx[col];
776840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
777aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
778840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
779840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
780840b8ebdSBarry Smith #else
781329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
782329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
783840b8ebdSBarry Smith #endif
784840b8ebdSBarry Smith       dx                *= epsilon;
7855904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
786840b8ebdSBarry Smith     }
7875904e1b1SBarry Smith   }
7885904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7895904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7905904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7915904e1b1SBarry Smith 
7925904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
7935904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
7945904e1b1SBarry Smith 
7955904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7965904e1b1SBarry Smith   /*
7975904e1b1SBarry Smith       Loop over each color
7985904e1b1SBarry Smith   */
7995904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8005904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8015904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8025904e1b1SBarry Smith     /*
8035904e1b1SBarry Smith        Loop over each column associated with color adding the
8045904e1b1SBarry Smith        perturbation to the vector w3.
8055904e1b1SBarry Smith     */
8065904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8075904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8085904e1b1SBarry Smith       dx  = xx[col];
8095904e1b1SBarry Smith       if (dx == 0.0) dx = 1.0;
8105904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8115904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8125904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8135904e1b1SBarry Smith #else
8145904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8155904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8165904e1b1SBarry Smith #endif
8175904e1b1SBarry Smith       dx            *= epsilon;
8185904e1b1SBarry Smith       w3_array[col] += dx;
8195904e1b1SBarry Smith     }
8205904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8215904e1b1SBarry Smith 
822840b8ebdSBarry Smith     /*
823840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
824840b8ebdSBarry Smith     */
82566f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
826840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
82766f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
828840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
8295904e1b1SBarry Smith 
830840b8ebdSBarry Smith     /*
831840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
832840b8ebdSBarry Smith     */
8333b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
834840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
835840b8ebdSBarry Smith       row    = coloring->rows[k][l];
836840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8375904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
838840b8ebdSBarry Smith       srow   = row + start;
839840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
840840b8ebdSBarry Smith     }
8413b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
842840b8ebdSBarry Smith   }
8435904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8445904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
845840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
846840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
847d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8483a40ed3dSBarry Smith   PetscFunctionReturn(0);
849840b8ebdSBarry Smith }
8503b28642cSBarry Smith 
8513b28642cSBarry Smith 
8524a2ae208SSatish Balay #undef __FUNCT__
8534a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
85440964ad5SBarry Smith /*@C
85540964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
85640964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
85740964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
85840964ad5SBarry Smith      called again.
85940964ad5SBarry Smith 
86040964ad5SBarry Smith    Collective on MatFDColoring
86140964ad5SBarry Smith 
86240964ad5SBarry Smith    Input  Parameters:
86340964ad5SBarry Smith .  c - the coloring context
86440964ad5SBarry Smith 
86540964ad5SBarry Smith    Level: intermediate
86640964ad5SBarry Smith 
86740964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
86840964ad5SBarry Smith 
86940964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
87040964ad5SBarry Smith 
87140964ad5SBarry Smith .keywords: Mat, finite differences, coloring
87240964ad5SBarry Smith @*/
87340964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c)
87440964ad5SBarry Smith {
87540964ad5SBarry Smith   PetscFunctionBegin;
87640964ad5SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
87740964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
87840964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
87940964ad5SBarry Smith   PetscFunctionReturn(0);
88040964ad5SBarry Smith }
88140964ad5SBarry Smith 
8823b28642cSBarry Smith 
883