xref: /petsc/src/mat/matfd/fdmatrix.c (revision 49b058dc59b472b6c22fb813bfc0819ef692001a)
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 
378b4fc646aSLois Curfman McInnes     Options Database Keys:
379ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
380ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
381ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
382639f9d9dSBarry Smith 
38315091d37SBarry Smith     Level: intermediate
38415091d37SBarry Smith 
3856831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
386d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
387d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
388d1c39f55SBarry Smith           MatFDColoringSetParameters()
389bbf0e169SBarry Smith @*/
390639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
391bbf0e169SBarry Smith {
392639f9d9dSBarry Smith   MatFDColoring c;
393639f9d9dSBarry Smith   MPI_Comm      comm;
394639f9d9dSBarry Smith   int           ierr,M,N;
395639f9d9dSBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
397d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
398639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
39929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
400639f9d9dSBarry Smith 
401f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4029194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
403b0a32e0cSBarry Smith   PetscLogObjectCreate(c);
404639f9d9dSBarry Smith 
405f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
406f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
407639f9d9dSBarry Smith   } else {
40829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
409639f9d9dSBarry Smith   }
410639f9d9dSBarry Smith 
41177d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
41277d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
413005c665bSBarry Smith   c->freq              = 1;
41440964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
41540964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
416*49b058dcSBarry Smith   c->currentcolor      = -1;
417005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
418639f9d9dSBarry Smith 
419639f9d9dSBarry Smith   *color = c;
420d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4213a40ed3dSBarry Smith   PetscFunctionReturn(0);
422bbf0e169SBarry Smith }
423bbf0e169SBarry Smith 
4244a2ae208SSatish Balay #undef __FUNCT__
4254a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
426bbf0e169SBarry Smith /*@C
427639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
428639f9d9dSBarry Smith     via MatFDColoringCreate().
429bbf0e169SBarry Smith 
430ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
431ef5ee4d1SLois Curfman McInnes 
432b4fc646aSLois Curfman McInnes     Input Parameter:
433639f9d9dSBarry Smith .   c - coloring context
434bbf0e169SBarry Smith 
43515091d37SBarry Smith     Level: intermediate
43615091d37SBarry Smith 
437639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
438bbf0e169SBarry Smith @*/
439639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
440bbf0e169SBarry Smith {
441263760aaSBarry Smith   int i,ierr;
442bbf0e169SBarry Smith 
4433a40ed3dSBarry Smith   PetscFunctionBegin;
4443a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
445d4bb536fSBarry Smith 
446639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
447606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
448606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
449606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
45030b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
451bbf0e169SBarry Smith   }
452606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
453606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
454606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
455606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
456606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
45730b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4584f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
459005c665bSBarry Smith   if (c->w1) {
460005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
461005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
462005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
463005c665bSBarry Smith   }
464b0a32e0cSBarry Smith   PetscLogObjectDestroy(c);
465639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
467bbf0e169SBarry Smith }
46843a90d84SBarry Smith 
4694a2ae208SSatish Balay #undef __FUNCT__
470*49b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
471*49b058dcSBarry Smith /*@C
472*49b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
473*49b058dcSBarry Smith       that are currently being perturbed.
474*49b058dcSBarry Smith 
475*49b058dcSBarry Smith     Not Collective
476*49b058dcSBarry Smith 
477*49b058dcSBarry Smith     Input Parameters:
478*49b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
479*49b058dcSBarry Smith 
480*49b058dcSBarry Smith     Output Parameters:
481*49b058dcSBarry Smith +   n - the number of local columns being perturbed
482*49b058dcSBarry Smith -   cols - the column indices, in global numbering
483*49b058dcSBarry Smith 
484*49b058dcSBarry Smith    Level: intermediate
485*49b058dcSBarry Smith 
486*49b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
487*49b058dcSBarry Smith 
488*49b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
489*49b058dcSBarry Smith @*/
490*49b058dcSBarry Smith EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int **cols)
491*49b058dcSBarry Smith {
492*49b058dcSBarry Smith   PetscFunctionBegin;
493*49b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
494*49b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
495*49b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
496*49b058dcSBarry Smith   } else {
497*49b058dcSBarry Smith     *n = 0;
498*49b058dcSBarry Smith   }
499*49b058dcSBarry Smith   PetscFunctionReturn(0);
500*49b058dcSBarry Smith }
501*49b058dcSBarry Smith 
502*49b058dcSBarry Smith #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
50443a90d84SBarry Smith /*@
505e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
506e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
50743a90d84SBarry Smith 
508fee21e36SBarry Smith     Collective on MatFDColoring
509fee21e36SBarry Smith 
510ef5ee4d1SLois Curfman McInnes     Input Parameters:
511ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
512ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
513ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
514ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
515ef5ee4d1SLois Curfman McInnes 
5168bba8e72SBarry Smith    Options Database Keys:
517ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
5188bba8e72SBarry Smith 
51915091d37SBarry Smith    Level: intermediate
52015091d37SBarry Smith 
52143a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
52243a90d84SBarry Smith 
52343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
52443a90d84SBarry Smith @*/
525005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
52643a90d84SBarry Smith {
5275904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
5283a7fca6bSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
529ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
530ea709b57SSatish Balay   PetscScalar   *vscale_array;
531329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
532005c665bSBarry Smith   Vec           w1,w2,w3;
533005c665bSBarry Smith   void          *fctx = coloring->fctx;
534f1af5d2fSBarry Smith   PetscTruth    flg;
535005c665bSBarry Smith 
536a2e34c3dSBarry Smith 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
538e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
539e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
540e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
541e0907662SLois Curfman McInnes 
54240964ad5SBarry Smith   if (coloring->usersetsrecompute) {
54340964ad5SBarry Smith     if (!coloring->recompute) {
54440964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
54576d8deaeSBarry Smith       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
54640964ad5SBarry Smith       PetscFunctionReturn(0);
54740964ad5SBarry Smith     } else {
54840964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
54940964ad5SBarry Smith     }
55040964ad5SBarry Smith   }
55140964ad5SBarry Smith 
552d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
55376d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
55476d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
55576d8deaeSBarry Smith   } else {
55676d8deaeSBarry Smith 
557005c665bSBarry Smith     if (!coloring->w1) {
558005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
559b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w1);
560005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
561b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w2);
562005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
563b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w3);
564005c665bSBarry Smith     }
565005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
56643a90d84SBarry Smith 
5674c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
568b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
569f1af5d2fSBarry Smith     if (flg) {
570b0a32e0cSBarry Smith       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
571e0907662SLois Curfman McInnes     } else {
57243a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
573e0907662SLois Curfman McInnes     }
57443a90d84SBarry Smith 
57543a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
57643a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
577be2fbe1fSBarry Smith 
5783a7fca6bSBarry Smith     /*
5793a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5803a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5813a7fca6bSBarry Smith     */
5823a7fca6bSBarry Smith     if (coloring->F) {
5833a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5843a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5853a7fca6bSBarry Smith       if (m1 != m2) {
5863a7fca6bSBarry Smith 	coloring->F = 0;
5873a7fca6bSBarry Smith       }
5883a7fca6bSBarry Smith     }
5893a7fca6bSBarry Smith 
590be2fbe1fSBarry Smith     if (coloring->F) {
591be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
592be2fbe1fSBarry Smith       coloring->F = 0;
593be2fbe1fSBarry Smith     } else {
59443a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
595be2fbe1fSBarry Smith     }
59643a90d84SBarry Smith 
59743a90d84SBarry Smith     /*
5989782cf98SBarry Smith        Compute all the scale factors and share with other processors
5999782cf98SBarry Smith     */
6009782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
6014f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
6029782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
6039782cf98SBarry Smith       /*
6049782cf98SBarry Smith 	Loop over each column associated with color adding the
6059782cf98SBarry Smith 	perturbation to the vector w3.
6069782cf98SBarry Smith       */
6079782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
6089782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
6099782cf98SBarry Smith 	dx  = xx[col];
6109782cf98SBarry Smith 	if (dx == 0.0) dx = 1.0;
6119782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6129782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
6139782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
6149782cf98SBarry Smith #else
6159782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6169782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6179782cf98SBarry Smith #endif
6189782cf98SBarry Smith 	dx                *= epsilon;
61930b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
6209782cf98SBarry Smith       }
6219782cf98SBarry Smith     }
622a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
62330b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62430b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6259782cf98SBarry Smith 
626ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
627ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
628b0a32e0cSBarry Smith 
62930b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
63030b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
63130b34957SBarry Smith 
63230b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6339782cf98SBarry Smith     /*
63443a90d84SBarry Smith       Loop over each color
63543a90d84SBarry Smith     */
63643a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
637*49b058dcSBarry Smith       coloring->currentcolor = k;
63843a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
63942460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
64043a90d84SBarry Smith       /*
64143a90d84SBarry Smith 	Loop over each column associated with color adding the
64243a90d84SBarry Smith 	perturbation to the vector w3.
64343a90d84SBarry Smith       */
64443a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
64543a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
64642460c72SBarry Smith 	dx  = xx[col];
647ae09f205SBarry Smith 	if (dx == 0.0) dx = 1.0;
648aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
649ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
650ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
65143a90d84SBarry Smith #else
652329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
653329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
65443a90d84SBarry Smith #endif
65543a90d84SBarry Smith 	dx            *= epsilon;
65629bbc08cSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
65742460c72SBarry Smith 	w3_array[col] += dx;
65843a90d84SBarry Smith       }
65942460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6603b28642cSBarry Smith 
66143a90d84SBarry Smith       /*
662e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
66343a90d84SBarry Smith       */
664a2e34c3dSBarry Smith 
66543a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
66643a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6679782cf98SBarry Smith 
66843a90d84SBarry Smith       /*
669e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
67043a90d84SBarry Smith       */
6713b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
67243a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
67343a90d84SBarry Smith 	row    = coloring->rows[k][l];
67443a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6755904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
67643a90d84SBarry Smith 	srow   = row + start;
67743a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
67843a90d84SBarry Smith       }
6793b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
68043a90d84SBarry Smith     }
681*49b058dcSBarry Smith     coloring->currentcolor = -1;
68230b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
68342460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
68443a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68543a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68676d8deaeSBarry Smith   }
687d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
688a2e34c3dSBarry Smith 
689b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
690a2e34c3dSBarry Smith   if (flg) {
691a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
692a2e34c3dSBarry Smith   }
6933a40ed3dSBarry Smith   PetscFunctionReturn(0);
69443a90d84SBarry Smith }
695840b8ebdSBarry Smith 
6964a2ae208SSatish Balay #undef __FUNCT__
6974a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
698840b8ebdSBarry Smith /*@
699840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
700840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
701840b8ebdSBarry Smith 
702fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
703fee21e36SBarry Smith 
704ef5ee4d1SLois Curfman McInnes     Input Parameters:
7053b28642cSBarry Smith +   mat - location to store Jacobian
706ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
707ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
708ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
709ef5ee4d1SLois Curfman McInnes 
710840b8ebdSBarry Smith    Options Database Keys:
711ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
712840b8ebdSBarry Smith 
71315091d37SBarry Smith    Level: intermediate
71415091d37SBarry Smith 
715840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
716840b8ebdSBarry Smith 
717840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
718840b8ebdSBarry Smith @*/
719329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
720840b8ebdSBarry Smith {
721329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
7225904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
723ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
724ea709b57SSatish Balay   PetscScalar   *vscale_array;
725329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
726840b8ebdSBarry Smith   Vec           w1,w2,w3;
727840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
728f1af5d2fSBarry Smith   PetscTruth    flg;
729840b8ebdSBarry Smith 
7303a40ed3dSBarry Smith   PetscFunctionBegin;
731840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
732840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
733840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
734840b8ebdSBarry Smith 
735d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
736840b8ebdSBarry Smith   if (!coloring->w1) {
737840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
738b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w1);
739840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
740b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w2);
741840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
742b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w3);
743840b8ebdSBarry Smith   }
744840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
745840b8ebdSBarry Smith 
7465904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
747b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
748f1af5d2fSBarry Smith   if (flg) {
749b0a32e0cSBarry Smith     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
750840b8ebdSBarry Smith   } else {
751840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
752840b8ebdSBarry Smith   }
753840b8ebdSBarry Smith 
754840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
755840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
756840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
757840b8ebdSBarry Smith 
758840b8ebdSBarry Smith   /*
7595904e1b1SBarry Smith       Compute all the scale factors and share with other processors
760840b8ebdSBarry Smith   */
7615904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7625904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
763840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
764840b8ebdSBarry Smith     /*
765840b8ebdSBarry Smith        Loop over each column associated with color adding the
766840b8ebdSBarry Smith        perturbation to the vector w3.
767840b8ebdSBarry Smith     */
768840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
769840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7705904e1b1SBarry Smith       dx  = xx[col];
771840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
772aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
773840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
774840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
775840b8ebdSBarry Smith #else
776329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
777329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
778840b8ebdSBarry Smith #endif
779840b8ebdSBarry Smith       dx                *= epsilon;
7805904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
781840b8ebdSBarry Smith     }
7825904e1b1SBarry Smith   }
7835904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7845904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7855904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7865904e1b1SBarry Smith 
7875904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
7885904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
7895904e1b1SBarry Smith 
7905904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7915904e1b1SBarry Smith   /*
7925904e1b1SBarry Smith       Loop over each color
7935904e1b1SBarry Smith   */
7945904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
7955904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
7965904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
7975904e1b1SBarry Smith     /*
7985904e1b1SBarry Smith        Loop over each column associated with color adding the
7995904e1b1SBarry Smith        perturbation to the vector w3.
8005904e1b1SBarry Smith     */
8015904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8025904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8035904e1b1SBarry Smith       dx  = xx[col];
8045904e1b1SBarry Smith       if (dx == 0.0) dx = 1.0;
8055904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8065904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8075904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8085904e1b1SBarry Smith #else
8095904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8105904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8115904e1b1SBarry Smith #endif
8125904e1b1SBarry Smith       dx            *= epsilon;
8135904e1b1SBarry Smith       w3_array[col] += dx;
8145904e1b1SBarry Smith     }
8155904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8165904e1b1SBarry Smith 
817840b8ebdSBarry Smith     /*
818840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
819840b8ebdSBarry Smith     */
820840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
821840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
8225904e1b1SBarry Smith 
823840b8ebdSBarry Smith     /*
824840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
825840b8ebdSBarry Smith     */
8263b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
827840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
828840b8ebdSBarry Smith       row    = coloring->rows[k][l];
829840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8305904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
831840b8ebdSBarry Smith       srow   = row + start;
832840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
833840b8ebdSBarry Smith     }
8343b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
835840b8ebdSBarry Smith   }
8365904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8375904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
838840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
839840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
840d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8413a40ed3dSBarry Smith   PetscFunctionReturn(0);
842840b8ebdSBarry Smith }
8433b28642cSBarry Smith 
8443b28642cSBarry Smith 
8454a2ae208SSatish Balay #undef __FUNCT__
8464a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
84740964ad5SBarry Smith /*@C
84840964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
84940964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
85040964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
85140964ad5SBarry Smith      called again.
85240964ad5SBarry Smith 
85340964ad5SBarry Smith    Collective on MatFDColoring
85440964ad5SBarry Smith 
85540964ad5SBarry Smith    Input  Parameters:
85640964ad5SBarry Smith .  c - the coloring context
85740964ad5SBarry Smith 
85840964ad5SBarry Smith    Level: intermediate
85940964ad5SBarry Smith 
86040964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
86140964ad5SBarry Smith 
86240964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
86340964ad5SBarry Smith 
86440964ad5SBarry Smith .keywords: Mat, finite differences, coloring
86540964ad5SBarry Smith @*/
86640964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c)
86740964ad5SBarry Smith {
86840964ad5SBarry Smith   PetscFunctionBegin;
86940964ad5SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
87040964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
87140964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
87240964ad5SBarry Smith   PetscFunctionReturn(0);
87340964ad5SBarry Smith }
87440964ad5SBarry Smith 
8753b28642cSBarry Smith 
876