xref: /petsc/src/mat/matfd/fdmatrix.c (revision fee21e364a2af8f69c0e7984443fdef19f844ae9)
1840b8ebdSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*fee21e36SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.31 1998/04/09 04:12:18 bsmith Exp bsmith $";
4bbf0e169SBarry Smith #endif
5bbf0e169SBarry Smith 
6bbf0e169SBarry Smith /*
7639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
8639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
9bbf0e169SBarry Smith */
10bbf0e169SBarry Smith 
11bbf0e169SBarry Smith #include "petsc.h"
12bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
13bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
14bbf0e169SBarry Smith #include "pinclude/pviewer.h"
15bbf0e169SBarry Smith 
165615d1e5SSatish Balay #undef __FUNC__
17d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw"
18005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
19005c665bSBarry Smith {
20005c665bSBarry Smith   int         ierr,i,j,pause;
21005c665bSBarry Smith   PetscTruth  isnull;
22005c665bSBarry Smith   Draw        draw;
23d4bb536fSBarry Smith   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
24005c665bSBarry Smith   DrawButton  button;
25005c665bSBarry Smith 
263a40ed3dSBarry Smith   PetscFunctionBegin;
27005c665bSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
283a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
292bdab257SBarry Smith   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
30005c665bSBarry Smith 
31005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
32005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
33005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
34005c665bSBarry Smith 
35005c665bSBarry Smith   /* loop over colors  */
36005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
37005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
38005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
39005c665bSBarry Smith       x = fd->columnsforrow[i][j];
405cd90555SBarry Smith       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); CHKERRQ(ierr);
41005c665bSBarry Smith     }
42005c665bSBarry Smith   }
432bdab257SBarry Smith   ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr);
44005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
453a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
465cd90555SBarry Smith   ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
475cd90555SBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);  CHKERRQ(ierr);
48005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
492bdab257SBarry Smith     ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
50005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
51005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
52005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
53005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
54005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
55005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
56005c665bSBarry Smith     w *= scale; h *= scale;
57005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
58005c665bSBarry Smith     /* loop over colors  */
59005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
60005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
61005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
62005c665bSBarry Smith         x = fd->columnsforrow[i][j];
635cd90555SBarry Smith         ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); CHKERRQ(ierr);
64005c665bSBarry Smith       }
65005c665bSBarry Smith     }
665cd90555SBarry Smith     ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
675cd90555SBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);  CHKERRQ(ierr);
68005c665bSBarry Smith   }
69005c665bSBarry Smith 
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71005c665bSBarry Smith }
72005c665bSBarry Smith 
73005c665bSBarry Smith #undef __FUNC__
74d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView"
75bbf0e169SBarry Smith /*@C
76639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
77bbf0e169SBarry Smith 
78b4fc646aSLois Curfman McInnes    Input  Parameters:
79b4fc646aSLois Curfman McInnes .  c - the coloring context
80b4fc646aSLois Curfman McInnes .  viewer - visualization context
81b4fc646aSLois Curfman McInnes 
82*fee21e36SBarry Smith    Collective on MatFDColoring unless Viewer is VIEWER_STDOUT_SELF
83*fee21e36SBarry Smith 
84b4fc646aSLois Curfman McInnes    Notes:
85b4fc646aSLois Curfman McInnes    The available visualization contexts include
86b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_SELF - standard output (default)
87b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_WORLD - synchronized standard
88b4fc646aSLois Curfman McInnes $       output where only the first processor opens
89b4fc646aSLois Curfman McInnes $       the file.  All other processors send their
90b4fc646aSLois Curfman McInnes $       data to the first processor to print.
91b4fc646aSLois Curfman McInnes $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
92bbf0e169SBarry Smith 
93639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
94005c665bSBarry Smith 
95b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
96bbf0e169SBarry Smith @*/
97b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
98bbf0e169SBarry Smith {
99005c665bSBarry Smith   ViewerType vtype;
100639f9d9dSBarry Smith   int        i,j,format,ierr;
101b4fc646aSLois Curfman McInnes   FILE       *fd;
102bbf0e169SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
104b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
105b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
106b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
107bbf0e169SBarry Smith 
108005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
109005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
110b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
1113a40ed3dSBarry Smith     PetscFunctionReturn(0);
1125cd90555SBarry Smith   } else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
113b4fc646aSLois Curfman McInnes     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
114b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
115b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
116b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
117b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
118ae09f205SBarry Smith 
119ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
120ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
121b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
122b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
123b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
124b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
125b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
126639f9d9dSBarry Smith         }
127b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
128b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
129b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
130b4fc646aSLois Curfman McInnes         }
131bbf0e169SBarry Smith       }
132bbf0e169SBarry Smith     }
1335cd90555SBarry Smith   } else {
1345cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
135bbf0e169SBarry Smith   }
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
137639f9d9dSBarry Smith }
138639f9d9dSBarry Smith 
1395615d1e5SSatish Balay #undef __FUNC__
1405615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
141639f9d9dSBarry Smith /*@
142b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
143b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
144639f9d9dSBarry Smith 
145639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
146639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
147639f9d9dSBarry Smith $          = error_rel*umin                    else
148639f9d9dSBarry Smith $
149639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
150639f9d9dSBarry Smith 
151639f9d9dSBarry Smith    Input Parameters:
152b4fc646aSLois Curfman McInnes .  coloring - the coloring context
153639f9d9dSBarry Smith .  error_rel - relative error
154639f9d9dSBarry Smith .  umin - minimum allowable u-value
155639f9d9dSBarry Smith 
156*fee21e36SBarry Smith    Collective on MatFDColoring
157*fee21e36SBarry Smith 
158b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
159b4fc646aSLois Curfman McInnes 
160b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
161639f9d9dSBarry Smith @*/
162639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
163639f9d9dSBarry Smith {
1643a40ed3dSBarry Smith   PetscFunctionBegin;
165639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
166639f9d9dSBarry Smith 
167639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
168639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170639f9d9dSBarry Smith }
171639f9d9dSBarry Smith 
1725615d1e5SSatish Balay #undef __FUNC__
173005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
174005c665bSBarry Smith /*@
175e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
176e0907662SLois Curfman McInnes    matrices.
177005c665bSBarry Smith 
178005c665bSBarry Smith    Input Parameters:
179b4fc646aSLois Curfman McInnes .  coloring - the coloring context
180e0907662SLois Curfman McInnes .  freq - frequency (default is 1)
181e0907662SLois Curfman McInnes 
182*fee21e36SBarry Smith    Collective on MatFDColoring
183*fee21e36SBarry Smith 
184e0907662SLois Curfman McInnes    Notes:
185e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
186e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
187e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
188e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
189e0907662SLois Curfman McInnes 
190e0907662SLois Curfman McInnes    Options Database Keys:
191e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq>
192005c665bSBarry Smith 
193b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
194005c665bSBarry Smith @*/
195005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
196005c665bSBarry Smith {
1973a40ed3dSBarry Smith   PetscFunctionBegin;
198005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
199005c665bSBarry Smith 
200005c665bSBarry Smith   matfd->freq = freq;
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
202005c665bSBarry Smith }
203005c665bSBarry Smith 
204005c665bSBarry Smith #undef __FUNC__
205ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
206ff0cfa39SBarry Smith /*@
207ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
208ff0cfa39SBarry Smith    matrices.
209ff0cfa39SBarry Smith 
210ff0cfa39SBarry Smith    Input Parameters:
211ff0cfa39SBarry Smith .  coloring - the coloring context
212ff0cfa39SBarry Smith 
213ff0cfa39SBarry Smith    Output Parameters:
214ff0cfa39SBarry Smith .  freq - frequency (default is 1)
215ff0cfa39SBarry Smith 
216*fee21e36SBarry Smith    Not Collective
217*fee21e36SBarry Smith 
218ff0cfa39SBarry Smith    Notes:
219ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
220ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
221ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
222ff0cfa39SBarry Smith    <freq> nonlinear iterations.
223ff0cfa39SBarry Smith 
224ff0cfa39SBarry Smith    Options Database Keys:
225ff0cfa39SBarry Smith $  -mat_fd_coloring_freq <freq>
226ff0cfa39SBarry Smith 
227ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
228ff0cfa39SBarry Smith @*/
229ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
230ff0cfa39SBarry Smith {
2313a40ed3dSBarry Smith   PetscFunctionBegin;
232ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
233ff0cfa39SBarry Smith 
234ff0cfa39SBarry Smith   *freq = matfd->freq;
2353a40ed3dSBarry Smith   PetscFunctionReturn(0);
236ff0cfa39SBarry Smith }
237ff0cfa39SBarry Smith 
238ff0cfa39SBarry Smith #undef __FUNC__
239005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
240d64ed03dSBarry Smith /*@C
241005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
242005c665bSBarry Smith 
243005c665bSBarry Smith    Input Parameters:
244b4fc646aSLois Curfman McInnes .  coloring - the coloring context
245005c665bSBarry Smith .  f - the function
246a10b0f90SLois Curfman McInnes .  fctx - the optional user-defined function context
247005c665bSBarry Smith 
248*fee21e36SBarry Smith    Collective on MatFDColoring
249*fee21e36SBarry Smith 
250b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
251005c665bSBarry Smith @*/
252840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
253005c665bSBarry Smith {
2543a40ed3dSBarry Smith   PetscFunctionBegin;
255005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
256005c665bSBarry Smith 
257005c665bSBarry Smith   matfd->f    = f;
258005c665bSBarry Smith   matfd->fctx = fctx;
259005c665bSBarry Smith 
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261005c665bSBarry Smith }
262005c665bSBarry Smith 
263005c665bSBarry Smith #undef __FUNC__
264d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
265639f9d9dSBarry Smith /*@
266b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
267639f9d9dSBarry Smith    the options database.
268639f9d9dSBarry Smith 
269b4fc646aSLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
270639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
271639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
272ae09f205SBarry Smith $          = error_rel*umin                      else
273639f9d9dSBarry Smith $
274639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
275639f9d9dSBarry Smith 
276639f9d9dSBarry Smith    Input Parameters:
277b4fc646aSLois Curfman McInnes .  coloring - the coloring context
278639f9d9dSBarry Smith 
279*fee21e36SBarry Smith    Collective on MatFDColoring
280*fee21e36SBarry Smith 
281b4fc646aSLois Curfman McInnes    Options Database Keys:
282b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_error <err>, where <err> is the square root
283b4fc646aSLois Curfman McInnes $           of relative error in the function
284b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_umin  <umin>, where umin is described above
285e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
286e0907662SLois Curfman McInnes $           computing a new Jacobian
287ca71c51bSBarry Smith $  -mat_fd_coloring_view
288ca71c51bSBarry Smith $  -mat_fd_coloring_view_info
289ca71c51bSBarry Smith $  -mat_fd_coloring_view_draw
290639f9d9dSBarry Smith 
291b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
292639f9d9dSBarry Smith @*/
293639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
294639f9d9dSBarry Smith {
295005c665bSBarry Smith   int    ierr,flag,freq = 1;
296639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
2973a40ed3dSBarry Smith 
2983a40ed3dSBarry Smith   PetscFunctionBegin;
299639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
300639f9d9dSBarry Smith 
301639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
302639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
303639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
304005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
305005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
306005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
307639f9d9dSBarry Smith   if (flag) {
308639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
309639f9d9dSBarry Smith   }
3103a40ed3dSBarry Smith   PetscFunctionReturn(0);
311639f9d9dSBarry Smith }
312639f9d9dSBarry Smith 
3135615d1e5SSatish Balay #undef __FUNC__
314d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
315639f9d9dSBarry Smith /*@
316639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
317639f9d9dSBarry Smith     using coloring.
318639f9d9dSBarry Smith 
319639f9d9dSBarry Smith     Input Parameter:
320639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
321639f9d9dSBarry Smith 
322*fee21e36SBarry Smith    Collective on MatFDColoring
323*fee21e36SBarry Smith 
324639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
325639f9d9dSBarry Smith @*/
326639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
327639f9d9dSBarry Smith {
3283a40ed3dSBarry Smith   PetscFunctionBegin;
329639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
330639f9d9dSBarry Smith 
33176be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
33276be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
33376be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
33476be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");
33576be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");
33676be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");
3373a40ed3dSBarry Smith   PetscFunctionReturn(0);
338005c665bSBarry Smith }
339005c665bSBarry Smith 
340005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
341005c665bSBarry Smith {
342005c665bSBarry Smith   int ierr,flg;
343005c665bSBarry Smith 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
345005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
346005c665bSBarry Smith   if (flg) {
347f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
348005c665bSBarry Smith   }
349ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
350ae09f205SBarry Smith   if (flg) {
351f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
352f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
353f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
354ae09f205SBarry Smith   }
355005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
356005c665bSBarry Smith   if (flg) {
357005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
358005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
359005c665bSBarry Smith   }
3603a40ed3dSBarry Smith   PetscFunctionReturn(0);
361bbf0e169SBarry Smith }
362bbf0e169SBarry Smith 
3635615d1e5SSatish Balay #undef __FUNC__
3645615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
365bbf0e169SBarry Smith /*@C
366639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
367639f9d9dSBarry Smith    computation of Jacobians.
368bbf0e169SBarry Smith 
369639f9d9dSBarry Smith    Input Parameters:
370639f9d9dSBarry Smith .  mat - the matrix containing the nonzero structure of the Jacobian
371639f9d9dSBarry Smith .  iscoloring - the coloring of the matrix
372bbf0e169SBarry Smith 
373bbf0e169SBarry Smith     Output Parameter:
374639f9d9dSBarry Smith .   color - the new coloring context
375bbf0e169SBarry Smith 
376*fee21e36SBarry Smith    Collective on Mat
377*fee21e36SBarry Smith 
378b4fc646aSLois Curfman McInnes     Options Database Keys:
379b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view
380b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_draw
381b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_info
382639f9d9dSBarry Smith 
383639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
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;
392639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
393e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
394639f9d9dSBarry Smith 
395639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
396f830108cSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
397639f9d9dSBarry Smith   PLogObjectCreate(c);
398639f9d9dSBarry Smith 
399f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
400f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
401639f9d9dSBarry Smith   } else {
402e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
403639f9d9dSBarry Smith   }
404639f9d9dSBarry Smith 
405639f9d9dSBarry Smith   c->error_rel = 1.e-8;
406ae09f205SBarry Smith   c->umin      = 1.e-6;
407005c665bSBarry Smith   c->freq      = 1;
408005c665bSBarry Smith 
409005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
410639f9d9dSBarry Smith 
411639f9d9dSBarry Smith   *color = c;
412639f9d9dSBarry Smith 
4133a40ed3dSBarry Smith   PetscFunctionReturn(0);
414bbf0e169SBarry Smith }
415bbf0e169SBarry Smith 
4165615d1e5SSatish Balay #undef __FUNC__
417d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
418bbf0e169SBarry Smith /*@C
419639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
420639f9d9dSBarry Smith     via MatFDColoringCreate().
421bbf0e169SBarry Smith 
422b4fc646aSLois Curfman McInnes     Input Parameter:
423639f9d9dSBarry Smith .   c - coloring context
424bbf0e169SBarry Smith 
425*fee21e36SBarry Smith     Collective on MatFDColoring
426*fee21e36SBarry Smith 
427639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
428bbf0e169SBarry Smith @*/
429639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
430bbf0e169SBarry Smith {
431263760aaSBarry Smith   int i,ierr;
432bbf0e169SBarry Smith 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
4343a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
435d4bb536fSBarry Smith 
436639f9d9dSBarry Smith 
437639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
438639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
439639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
440639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
441bbf0e169SBarry Smith   }
442639f9d9dSBarry Smith   PetscFree(c->ncolumns);
443639f9d9dSBarry Smith   PetscFree(c->columns);
444639f9d9dSBarry Smith   PetscFree(c->nrows);
445639f9d9dSBarry Smith   PetscFree(c->rows);
446639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
447639f9d9dSBarry Smith   PetscFree(c->scale);
448005c665bSBarry Smith   if (c->w1) {
449005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
450005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
451005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
452005c665bSBarry Smith   }
453639f9d9dSBarry Smith   PLogObjectDestroy(c);
454639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4553a40ed3dSBarry Smith   PetscFunctionReturn(0);
456bbf0e169SBarry Smith }
45743a90d84SBarry Smith 
458005c665bSBarry Smith #include "snes.h"
459005c665bSBarry Smith 
4605615d1e5SSatish Balay #undef __FUNC__
4615615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
46243a90d84SBarry Smith /*@
463e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
464e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
46543a90d84SBarry Smith 
46643a90d84SBarry Smith     Input Parameters:
46743a90d84SBarry Smith .   mat - location to store Jacobian
46843a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
46943a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
470005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
47143a90d84SBarry Smith 
472*fee21e36SBarry Smith    Collective on MatFDColoring
473*fee21e36SBarry Smith 
4748bba8e72SBarry Smith    Options Database Keys:
4758bba8e72SBarry Smith $  -mat_fd_coloring_freq <freq>
4768bba8e72SBarry Smith 
47743a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
47843a90d84SBarry Smith 
47943a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
48043a90d84SBarry Smith @*/
481005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
48243a90d84SBarry Smith {
483e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
48443a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
48543a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
48643a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
487005c665bSBarry Smith   Vec           w1,w2,w3;
488840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
489005c665bSBarry Smith   void          *fctx = coloring->fctx;
490005c665bSBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
492e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
493e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
494e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
495e0907662SLois Curfman McInnes 
496005c665bSBarry Smith 
497005c665bSBarry Smith   if (!coloring->w1) {
498005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
499005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
500005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
501005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
502005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
503005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
504005c665bSBarry Smith   }
505005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
50643a90d84SBarry Smith 
507e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
508e0907662SLois Curfman McInnes   if (fg) {
509e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
510e0907662SLois Curfman McInnes   } else {
51143a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
512e0907662SLois Curfman McInnes   }
51343a90d84SBarry Smith 
51443a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
51543a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
51643a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
51743a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
51843a90d84SBarry Smith 
51943a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
52043a90d84SBarry Smith   /*
52143a90d84SBarry Smith       Loop over each color
52243a90d84SBarry Smith   */
52343a90d84SBarry Smith 
52443a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
52543a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
52643a90d84SBarry Smith     /*
52743a90d84SBarry Smith        Loop over each column associated with color adding the
52843a90d84SBarry Smith        perturbation to the vector w3.
52943a90d84SBarry Smith     */
53043a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
53143a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
53243a90d84SBarry Smith       dx  = xx[col-start];
533ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
5343a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
535ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
536ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
53743a90d84SBarry Smith #else
538ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
539ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
54043a90d84SBarry Smith #endif
54143a90d84SBarry Smith       dx          *= epsilon;
54243a90d84SBarry Smith       wscale[col] = 1.0/dx;
54343a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
54443a90d84SBarry Smith     }
54543a90d84SBarry Smith     VecRestoreArray(x1,&xx);
54643a90d84SBarry Smith     /*
547e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
54843a90d84SBarry Smith     */
54943a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
55043a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
55143a90d84SBarry Smith     /* Communicate scale to all processors */
5523a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
553ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
55443a90d84SBarry Smith #else
555ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
55643a90d84SBarry Smith #endif
55743a90d84SBarry Smith     /*
558e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
55943a90d84SBarry Smith     */
56043a90d84SBarry Smith     VecGetArray(w2,&y);
56143a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
56243a90d84SBarry Smith       row    = coloring->rows[k][l];
56343a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
56443a90d84SBarry Smith       y[row] *= scale[col];
56543a90d84SBarry Smith       srow   = row + start;
56643a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
56743a90d84SBarry Smith     }
56843a90d84SBarry Smith     VecRestoreArray(w2,&y);
56943a90d84SBarry Smith   }
57043a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
57143a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5723a40ed3dSBarry Smith   PetscFunctionReturn(0);
57343a90d84SBarry Smith }
574840b8ebdSBarry Smith 
575840b8ebdSBarry Smith #include "ts.h"
576840b8ebdSBarry Smith 
577840b8ebdSBarry Smith #undef __FUNC__
578840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
579840b8ebdSBarry Smith /*@
580840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
581840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
582840b8ebdSBarry Smith 
583840b8ebdSBarry Smith     Input Parameters:
584840b8ebdSBarry Smith .   mat - location to store Jacobian
585840b8ebdSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
586840b8ebdSBarry Smith .   x1 - location at which Jacobian is to be computed
587840b8ebdSBarry Smith .   sctx - optional context required by function (actually a SNES context)
588840b8ebdSBarry Smith 
589*fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
590*fee21e36SBarry Smith 
591840b8ebdSBarry Smith    Options Database Keys:
592840b8ebdSBarry Smith $  -mat_fd_coloring_freq <freq>
593840b8ebdSBarry Smith 
594840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
595840b8ebdSBarry Smith 
596840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
597840b8ebdSBarry Smith @*/
598840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
599840b8ebdSBarry Smith {
600840b8ebdSBarry Smith   int           k,fg,ierr,N,start,end,l,row,col,srow;
601840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
602840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
603840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
604840b8ebdSBarry Smith   Vec           w1,w2,w3;
605840b8ebdSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
606840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
607840b8ebdSBarry Smith 
6083a40ed3dSBarry Smith   PetscFunctionBegin;
609840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
610840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
611840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
612840b8ebdSBarry Smith 
613840b8ebdSBarry Smith 
614840b8ebdSBarry Smith   if (!coloring->w1) {
615840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
616840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
617840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
618840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
619840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
620840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
621840b8ebdSBarry Smith   }
622840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
623840b8ebdSBarry Smith 
624840b8ebdSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
625840b8ebdSBarry Smith   if (fg) {
626840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
627840b8ebdSBarry Smith   } else {
628840b8ebdSBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
629840b8ebdSBarry Smith   }
630840b8ebdSBarry Smith 
631840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
632840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
633840b8ebdSBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
634840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr);
635840b8ebdSBarry Smith 
636840b8ebdSBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
637840b8ebdSBarry Smith   /*
638840b8ebdSBarry Smith       Loop over each color
639840b8ebdSBarry Smith   */
640840b8ebdSBarry Smith 
641840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
642840b8ebdSBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
643840b8ebdSBarry Smith     /*
644840b8ebdSBarry Smith        Loop over each column associated with color adding the
645840b8ebdSBarry Smith        perturbation to the vector w3.
646840b8ebdSBarry Smith     */
647840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
648840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
649840b8ebdSBarry Smith       dx  = xx[col-start];
650840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
6513a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
652840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
653840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
654840b8ebdSBarry Smith #else
655840b8ebdSBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
656840b8ebdSBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
657840b8ebdSBarry Smith #endif
658840b8ebdSBarry Smith       dx          *= epsilon;
659840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
660840b8ebdSBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
661840b8ebdSBarry Smith     }
662840b8ebdSBarry Smith     VecRestoreArray(x1,&xx);
663840b8ebdSBarry Smith     /*
664840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
665840b8ebdSBarry Smith     */
666840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr);
667840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
668840b8ebdSBarry Smith     /* Communicate scale to all processors */
6693a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
670ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
671840b8ebdSBarry Smith #else
672ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
673840b8ebdSBarry Smith #endif
674840b8ebdSBarry Smith     /*
675840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
676840b8ebdSBarry Smith     */
677840b8ebdSBarry Smith     VecGetArray(w2,&y);
678840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
679840b8ebdSBarry Smith       row    = coloring->rows[k][l];
680840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
681840b8ebdSBarry Smith       y[row] *= scale[col];
682840b8ebdSBarry Smith       srow   = row + start;
683840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
684840b8ebdSBarry Smith     }
685840b8ebdSBarry Smith     VecRestoreArray(w2,&y);
686840b8ebdSBarry Smith   }
687840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
688840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6893a40ed3dSBarry Smith   PetscFunctionReturn(0);
690840b8ebdSBarry Smith }
691