xref: /petsc/src/mat/matfd/fdmatrix.c (revision 5cd905551efc76d042fcd2cc746dd1c50567705b)
1840b8ebdSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*5cd90555SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.30 1998/04/08 16:06:13 curfman 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];
40*5cd90555SBarry 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);}
46*5cd90555SBarry Smith   ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
47*5cd90555SBarry 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];
63*5cd90555SBarry Smith         ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); CHKERRQ(ierr);
64005c665bSBarry Smith       }
65005c665bSBarry Smith     }
66*5cd90555SBarry Smith     ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
67*5cd90555SBarry 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 
82b4fc646aSLois Curfman McInnes    Notes:
83b4fc646aSLois Curfman McInnes    The available visualization contexts include
84b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_SELF - standard output (default)
85b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_WORLD - synchronized standard
86b4fc646aSLois Curfman McInnes $       output where only the first processor opens
87b4fc646aSLois Curfman McInnes $       the file.  All other processors send their
88b4fc646aSLois Curfman McInnes $       data to the first processor to print.
89b4fc646aSLois Curfman McInnes $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
90bbf0e169SBarry Smith 
91639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
92005c665bSBarry Smith 
93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
94bbf0e169SBarry Smith @*/
95b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
96bbf0e169SBarry Smith {
97005c665bSBarry Smith   ViewerType vtype;
98639f9d9dSBarry Smith   int        i,j,format,ierr;
99b4fc646aSLois Curfman McInnes   FILE       *fd;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
102b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
103b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
104b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
105bbf0e169SBarry Smith 
106005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
107005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
108b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
1093a40ed3dSBarry Smith     PetscFunctionReturn(0);
110*5cd90555SBarry Smith   } else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
111b4fc646aSLois Curfman McInnes     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
112b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
113b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
114b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
115b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
116ae09f205SBarry Smith 
117ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
118ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
119b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
120b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
121b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
122b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
123b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
124639f9d9dSBarry Smith         }
125b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
126b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
127b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
128b4fc646aSLois Curfman McInnes         }
129bbf0e169SBarry Smith       }
130bbf0e169SBarry Smith     }
131*5cd90555SBarry Smith   } else {
132*5cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
133bbf0e169SBarry Smith   }
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135639f9d9dSBarry Smith }
136639f9d9dSBarry Smith 
1375615d1e5SSatish Balay #undef __FUNC__
1385615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
139639f9d9dSBarry Smith /*@
140b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
141b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
142639f9d9dSBarry Smith 
143639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
144639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
145639f9d9dSBarry Smith $          = error_rel*umin                    else
146639f9d9dSBarry Smith $
147639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
148639f9d9dSBarry Smith 
149639f9d9dSBarry Smith    Input Parameters:
150b4fc646aSLois Curfman McInnes .  coloring - the coloring context
151639f9d9dSBarry Smith .  error_rel - relative error
152639f9d9dSBarry Smith .  umin - minimum allowable u-value
153639f9d9dSBarry Smith 
154b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
155b4fc646aSLois Curfman McInnes 
156b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
157639f9d9dSBarry Smith @*/
158639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
159639f9d9dSBarry Smith {
1603a40ed3dSBarry Smith   PetscFunctionBegin;
161639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
162639f9d9dSBarry Smith 
163639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
164639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1653a40ed3dSBarry Smith   PetscFunctionReturn(0);
166639f9d9dSBarry Smith }
167639f9d9dSBarry Smith 
1685615d1e5SSatish Balay #undef __FUNC__
169005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
170005c665bSBarry Smith /*@
171e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
172e0907662SLois Curfman McInnes    matrices.
173005c665bSBarry Smith 
174005c665bSBarry Smith    Input Parameters:
175b4fc646aSLois Curfman McInnes .  coloring - the coloring context
176e0907662SLois Curfman McInnes .  freq - frequency (default is 1)
177e0907662SLois Curfman McInnes 
178e0907662SLois Curfman McInnes    Notes:
179e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
180e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
181e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
182e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
183e0907662SLois Curfman McInnes 
184e0907662SLois Curfman McInnes    Options Database Keys:
185e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq>
186005c665bSBarry Smith 
187b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
188005c665bSBarry Smith @*/
189005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
190005c665bSBarry Smith {
1913a40ed3dSBarry Smith   PetscFunctionBegin;
192005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
193005c665bSBarry Smith 
194005c665bSBarry Smith   matfd->freq = freq;
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
196005c665bSBarry Smith }
197005c665bSBarry Smith 
198005c665bSBarry Smith #undef __FUNC__
199ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
200ff0cfa39SBarry Smith /*@
201ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
202ff0cfa39SBarry Smith    matrices.
203ff0cfa39SBarry Smith 
204ff0cfa39SBarry Smith    Input Parameters:
205ff0cfa39SBarry Smith .  coloring - the coloring context
206ff0cfa39SBarry Smith 
207ff0cfa39SBarry Smith    Output Parameters:
208ff0cfa39SBarry Smith .  freq - frequency (default is 1)
209ff0cfa39SBarry Smith 
210ff0cfa39SBarry Smith    Notes:
211ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
212ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
213ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
214ff0cfa39SBarry Smith    <freq> nonlinear iterations.
215ff0cfa39SBarry Smith 
216ff0cfa39SBarry Smith    Options Database Keys:
217ff0cfa39SBarry Smith $  -mat_fd_coloring_freq <freq>
218ff0cfa39SBarry Smith 
219ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
220ff0cfa39SBarry Smith @*/
221ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
222ff0cfa39SBarry Smith {
2233a40ed3dSBarry Smith   PetscFunctionBegin;
224ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
225ff0cfa39SBarry Smith 
226ff0cfa39SBarry Smith   *freq = matfd->freq;
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
228ff0cfa39SBarry Smith }
229ff0cfa39SBarry Smith 
230ff0cfa39SBarry Smith #undef __FUNC__
231005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
232d64ed03dSBarry Smith /*@C
233005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
234005c665bSBarry Smith 
235005c665bSBarry Smith    Input Parameters:
236b4fc646aSLois Curfman McInnes .  coloring - the coloring context
237005c665bSBarry Smith .  f - the function
238a10b0f90SLois Curfman McInnes .  fctx - the optional user-defined function context
239005c665bSBarry Smith 
240b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
241005c665bSBarry Smith @*/
242840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
243005c665bSBarry Smith {
2443a40ed3dSBarry Smith   PetscFunctionBegin;
245005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
246005c665bSBarry Smith 
247005c665bSBarry Smith   matfd->f    = f;
248005c665bSBarry Smith   matfd->fctx = fctx;
249005c665bSBarry Smith 
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
251005c665bSBarry Smith }
252005c665bSBarry Smith 
253005c665bSBarry Smith #undef __FUNC__
254d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
255639f9d9dSBarry Smith /*@
256b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
257639f9d9dSBarry Smith    the options database.
258639f9d9dSBarry Smith 
259b4fc646aSLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
260639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
261639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
262ae09f205SBarry Smith $          = error_rel*umin                      else
263639f9d9dSBarry Smith $
264639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
265639f9d9dSBarry Smith 
266639f9d9dSBarry Smith    Input Parameters:
267b4fc646aSLois Curfman McInnes .  coloring - the coloring context
268639f9d9dSBarry Smith 
269b4fc646aSLois Curfman McInnes    Options Database Keys:
270b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_error <err>, where <err> is the square root
271b4fc646aSLois Curfman McInnes $           of relative error in the function
272b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_umin  <umin>, where umin is described above
273e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
274e0907662SLois Curfman McInnes $           computing a new Jacobian
275ca71c51bSBarry Smith $  -mat_fd_coloring_view
276ca71c51bSBarry Smith $  -mat_fd_coloring_view_info
277ca71c51bSBarry Smith $  -mat_fd_coloring_view_draw
278639f9d9dSBarry Smith 
279b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
280639f9d9dSBarry Smith @*/
281639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
282639f9d9dSBarry Smith {
283005c665bSBarry Smith   int    ierr,flag,freq = 1;
284639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
2853a40ed3dSBarry Smith 
2863a40ed3dSBarry Smith   PetscFunctionBegin;
287639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
288639f9d9dSBarry Smith 
289639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
290639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
291639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
292005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
293005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
294005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
295639f9d9dSBarry Smith   if (flag) {
296639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
297639f9d9dSBarry Smith   }
2983a40ed3dSBarry Smith   PetscFunctionReturn(0);
299639f9d9dSBarry Smith }
300639f9d9dSBarry Smith 
3015615d1e5SSatish Balay #undef __FUNC__
302d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
303639f9d9dSBarry Smith /*@
304639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
305639f9d9dSBarry Smith     using coloring.
306639f9d9dSBarry Smith 
307639f9d9dSBarry Smith     Input Parameter:
308639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
309639f9d9dSBarry Smith 
310639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
311639f9d9dSBarry Smith @*/
312639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
313639f9d9dSBarry Smith {
3143a40ed3dSBarry Smith   PetscFunctionBegin;
315639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
316639f9d9dSBarry Smith 
31776be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
31876be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
31976be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
32076be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");
32176be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");
32276be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
324005c665bSBarry Smith }
325005c665bSBarry Smith 
326005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
327005c665bSBarry Smith {
328005c665bSBarry Smith   int ierr,flg;
329005c665bSBarry Smith 
3303a40ed3dSBarry Smith   PetscFunctionBegin;
331005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
332005c665bSBarry Smith   if (flg) {
333f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
334005c665bSBarry Smith   }
335ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
336ae09f205SBarry Smith   if (flg) {
337f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
338f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
339f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
340ae09f205SBarry Smith   }
341005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
342005c665bSBarry Smith   if (flg) {
343005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
344005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
345005c665bSBarry Smith   }
3463a40ed3dSBarry Smith   PetscFunctionReturn(0);
347bbf0e169SBarry Smith }
348bbf0e169SBarry Smith 
3495615d1e5SSatish Balay #undef __FUNC__
3505615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
351bbf0e169SBarry Smith /*@C
352639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
353639f9d9dSBarry Smith    computation of Jacobians.
354bbf0e169SBarry Smith 
355639f9d9dSBarry Smith    Input Parameters:
356639f9d9dSBarry Smith .  mat - the matrix containing the nonzero structure of the Jacobian
357639f9d9dSBarry Smith .  iscoloring - the coloring of the matrix
358bbf0e169SBarry Smith 
359bbf0e169SBarry Smith     Output Parameter:
360639f9d9dSBarry Smith .   color - the new coloring context
361bbf0e169SBarry Smith 
362b4fc646aSLois Curfman McInnes     Options Database Keys:
363b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view
364b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_draw
365b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_info
366639f9d9dSBarry Smith 
367639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
368bbf0e169SBarry Smith @*/
369639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
370bbf0e169SBarry Smith {
371639f9d9dSBarry Smith   MatFDColoring c;
372639f9d9dSBarry Smith   MPI_Comm      comm;
373639f9d9dSBarry Smith   int           ierr,M,N;
374639f9d9dSBarry Smith 
3753a40ed3dSBarry Smith   PetscFunctionBegin;
376639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
377e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
378639f9d9dSBarry Smith 
379639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
380f830108cSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
381639f9d9dSBarry Smith   PLogObjectCreate(c);
382639f9d9dSBarry Smith 
383f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
384f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
385639f9d9dSBarry Smith   } else {
386e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
387639f9d9dSBarry Smith   }
388639f9d9dSBarry Smith 
389639f9d9dSBarry Smith   c->error_rel = 1.e-8;
390ae09f205SBarry Smith   c->umin      = 1.e-6;
391005c665bSBarry Smith   c->freq      = 1;
392005c665bSBarry Smith 
393005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
394639f9d9dSBarry Smith 
395639f9d9dSBarry Smith   *color = c;
396639f9d9dSBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionReturn(0);
398bbf0e169SBarry Smith }
399bbf0e169SBarry Smith 
4005615d1e5SSatish Balay #undef __FUNC__
401d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
402bbf0e169SBarry Smith /*@C
403639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
404639f9d9dSBarry Smith     via MatFDColoringCreate().
405bbf0e169SBarry Smith 
406b4fc646aSLois Curfman McInnes     Input Parameter:
407639f9d9dSBarry Smith .   c - coloring context
408bbf0e169SBarry Smith 
409639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
410bbf0e169SBarry Smith @*/
411639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
412bbf0e169SBarry Smith {
413263760aaSBarry Smith   int i,ierr;
414bbf0e169SBarry Smith 
4153a40ed3dSBarry Smith   PetscFunctionBegin;
4163a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
417d4bb536fSBarry Smith 
418639f9d9dSBarry Smith 
419639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
420639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
421639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
422639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
423bbf0e169SBarry Smith   }
424639f9d9dSBarry Smith   PetscFree(c->ncolumns);
425639f9d9dSBarry Smith   PetscFree(c->columns);
426639f9d9dSBarry Smith   PetscFree(c->nrows);
427639f9d9dSBarry Smith   PetscFree(c->rows);
428639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
429639f9d9dSBarry Smith   PetscFree(c->scale);
430005c665bSBarry Smith   if (c->w1) {
431005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
432005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
433005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
434005c665bSBarry Smith   }
435639f9d9dSBarry Smith   PLogObjectDestroy(c);
436639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4373a40ed3dSBarry Smith   PetscFunctionReturn(0);
438bbf0e169SBarry Smith }
43943a90d84SBarry Smith 
440005c665bSBarry Smith #include "snes.h"
441005c665bSBarry Smith 
4425615d1e5SSatish Balay #undef __FUNC__
4435615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
44443a90d84SBarry Smith /*@
445e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
446e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
44743a90d84SBarry Smith 
44843a90d84SBarry Smith     Input Parameters:
44943a90d84SBarry Smith .   mat - location to store Jacobian
45043a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
45143a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
452005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
45343a90d84SBarry Smith 
4548bba8e72SBarry Smith    Options Database Keys:
4558bba8e72SBarry Smith $  -mat_fd_coloring_freq <freq>
4568bba8e72SBarry Smith 
45743a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
45843a90d84SBarry Smith 
45943a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
46043a90d84SBarry Smith @*/
461005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
46243a90d84SBarry Smith {
463e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
46443a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
46543a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
46643a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
467005c665bSBarry Smith   Vec           w1,w2,w3;
468840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
469005c665bSBarry Smith   void          *fctx = coloring->fctx;
470005c665bSBarry Smith 
4713a40ed3dSBarry Smith   PetscFunctionBegin;
472e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
473e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
474e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
475e0907662SLois Curfman McInnes 
476005c665bSBarry Smith 
477005c665bSBarry Smith   if (!coloring->w1) {
478005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
479005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
480005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
481005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
482005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
483005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
484005c665bSBarry Smith   }
485005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
48643a90d84SBarry Smith 
487e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
488e0907662SLois Curfman McInnes   if (fg) {
489e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
490e0907662SLois Curfman McInnes   } else {
49143a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
492e0907662SLois Curfman McInnes   }
49343a90d84SBarry Smith 
49443a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
49543a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
49643a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
49743a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
49843a90d84SBarry Smith 
49943a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
50043a90d84SBarry Smith   /*
50143a90d84SBarry Smith       Loop over each color
50243a90d84SBarry Smith   */
50343a90d84SBarry Smith 
50443a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
50543a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
50643a90d84SBarry Smith     /*
50743a90d84SBarry Smith        Loop over each column associated with color adding the
50843a90d84SBarry Smith        perturbation to the vector w3.
50943a90d84SBarry Smith     */
51043a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
51143a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
51243a90d84SBarry Smith       dx  = xx[col-start];
513ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
5143a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
515ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
516ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
51743a90d84SBarry Smith #else
518ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
519ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
52043a90d84SBarry Smith #endif
52143a90d84SBarry Smith       dx          *= epsilon;
52243a90d84SBarry Smith       wscale[col] = 1.0/dx;
52343a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
52443a90d84SBarry Smith     }
52543a90d84SBarry Smith     VecRestoreArray(x1,&xx);
52643a90d84SBarry Smith     /*
527e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
52843a90d84SBarry Smith     */
52943a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
53043a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
53143a90d84SBarry Smith     /* Communicate scale to all processors */
5323a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
533ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
53443a90d84SBarry Smith #else
535ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
53643a90d84SBarry Smith #endif
53743a90d84SBarry Smith     /*
538e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
53943a90d84SBarry Smith     */
54043a90d84SBarry Smith     VecGetArray(w2,&y);
54143a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
54243a90d84SBarry Smith       row    = coloring->rows[k][l];
54343a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
54443a90d84SBarry Smith       y[row] *= scale[col];
54543a90d84SBarry Smith       srow   = row + start;
54643a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
54743a90d84SBarry Smith     }
54843a90d84SBarry Smith     VecRestoreArray(w2,&y);
54943a90d84SBarry Smith   }
55043a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
55143a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5523a40ed3dSBarry Smith   PetscFunctionReturn(0);
55343a90d84SBarry Smith }
554840b8ebdSBarry Smith 
555840b8ebdSBarry Smith #include "ts.h"
556840b8ebdSBarry Smith 
557840b8ebdSBarry Smith #undef __FUNC__
558840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
559840b8ebdSBarry Smith /*@
560840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
561840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
562840b8ebdSBarry Smith 
563840b8ebdSBarry Smith     Input Parameters:
564840b8ebdSBarry Smith .   mat - location to store Jacobian
565840b8ebdSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
566840b8ebdSBarry Smith .   x1 - location at which Jacobian is to be computed
567840b8ebdSBarry Smith .   sctx - optional context required by function (actually a SNES context)
568840b8ebdSBarry Smith 
569840b8ebdSBarry Smith    Options Database Keys:
570840b8ebdSBarry Smith $  -mat_fd_coloring_freq <freq>
571840b8ebdSBarry Smith 
572840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
573840b8ebdSBarry Smith 
574840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
575840b8ebdSBarry Smith @*/
576840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
577840b8ebdSBarry Smith {
578840b8ebdSBarry Smith   int           k,fg,ierr,N,start,end,l,row,col,srow;
579840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
580840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
581840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
582840b8ebdSBarry Smith   Vec           w1,w2,w3;
583840b8ebdSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
584840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
585840b8ebdSBarry Smith 
5863a40ed3dSBarry Smith   PetscFunctionBegin;
587840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
588840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
589840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
590840b8ebdSBarry Smith 
591840b8ebdSBarry Smith 
592840b8ebdSBarry Smith   if (!coloring->w1) {
593840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
594840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
595840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
596840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
597840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
598840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
599840b8ebdSBarry Smith   }
600840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
601840b8ebdSBarry Smith 
602840b8ebdSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
603840b8ebdSBarry Smith   if (fg) {
604840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
605840b8ebdSBarry Smith   } else {
606840b8ebdSBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
607840b8ebdSBarry Smith   }
608840b8ebdSBarry Smith 
609840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
610840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
611840b8ebdSBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
612840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr);
613840b8ebdSBarry Smith 
614840b8ebdSBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
615840b8ebdSBarry Smith   /*
616840b8ebdSBarry Smith       Loop over each color
617840b8ebdSBarry Smith   */
618840b8ebdSBarry Smith 
619840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
620840b8ebdSBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
621840b8ebdSBarry Smith     /*
622840b8ebdSBarry Smith        Loop over each column associated with color adding the
623840b8ebdSBarry Smith        perturbation to the vector w3.
624840b8ebdSBarry Smith     */
625840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
626840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
627840b8ebdSBarry Smith       dx  = xx[col-start];
628840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
6293a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
630840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
631840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
632840b8ebdSBarry Smith #else
633840b8ebdSBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
634840b8ebdSBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
635840b8ebdSBarry Smith #endif
636840b8ebdSBarry Smith       dx          *= epsilon;
637840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
638840b8ebdSBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
639840b8ebdSBarry Smith     }
640840b8ebdSBarry Smith     VecRestoreArray(x1,&xx);
641840b8ebdSBarry Smith     /*
642840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
643840b8ebdSBarry Smith     */
644840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr);
645840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
646840b8ebdSBarry Smith     /* Communicate scale to all processors */
6473a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
648ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
649840b8ebdSBarry Smith #else
650ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
651840b8ebdSBarry Smith #endif
652840b8ebdSBarry Smith     /*
653840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
654840b8ebdSBarry Smith     */
655840b8ebdSBarry Smith     VecGetArray(w2,&y);
656840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
657840b8ebdSBarry Smith       row    = coloring->rows[k][l];
658840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
659840b8ebdSBarry Smith       y[row] *= scale[col];
660840b8ebdSBarry Smith       srow   = row + start;
661840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
662840b8ebdSBarry Smith     }
663840b8ebdSBarry Smith     VecRestoreArray(w2,&y);
664840b8ebdSBarry Smith   }
665840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
666840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6673a40ed3dSBarry Smith   PetscFunctionReturn(0);
668840b8ebdSBarry Smith }
669