xref: /petsc/src/mat/matfd/fdmatrix.c (revision 840b8ebde3de12c1c951820c153e7e1ece966aad)
1*840b8ebdSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*840b8ebdSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.23 1997/10/12 23:24:25 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 
26005c665bSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
27005c665bSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
282bdab257SBarry Smith   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
29005c665bSBarry Smith 
30005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
31005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
32005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
33005c665bSBarry Smith 
34005c665bSBarry Smith   /* loop over colors  */
35005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
36005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
37005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
38005c665bSBarry Smith       x = fd->columnsforrow[i][j];
39005c665bSBarry Smith       DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
40005c665bSBarry Smith     }
41005c665bSBarry Smith   }
422bdab257SBarry Smith   ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr);
43005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
44005c665bSBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
45005c665bSBarry Smith   ierr = DrawCheckResizedWindow(draw);
462bdab257SBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
47005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
482bdab257SBarry Smith     ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
49005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
50005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
51005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
52005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
53005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
54005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
55005c665bSBarry Smith     w *= scale; h *= scale;
56005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
57005c665bSBarry Smith     /* loop over colors  */
58005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
59005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
60005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
61005c665bSBarry Smith         x = fd->columnsforrow[i][j];
62005c665bSBarry Smith         DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
63005c665bSBarry Smith       }
64005c665bSBarry Smith     }
65005c665bSBarry Smith     ierr = DrawCheckResizedWindow(draw);
662bdab257SBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
67005c665bSBarry Smith   }
68005c665bSBarry Smith 
69005c665bSBarry Smith   return 0;
70005c665bSBarry Smith }
71005c665bSBarry Smith 
72005c665bSBarry Smith #undef __FUNC__
73d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView"
74bbf0e169SBarry Smith /*@C
75639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
76bbf0e169SBarry Smith 
77b4fc646aSLois Curfman McInnes    Input  Parameters:
78b4fc646aSLois Curfman McInnes .  c - the coloring context
79b4fc646aSLois Curfman McInnes .  viewer - visualization context
80b4fc646aSLois Curfman McInnes 
81b4fc646aSLois Curfman McInnes    Notes:
82b4fc646aSLois Curfman McInnes    The available visualization contexts include
83b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_SELF - standard output (default)
84b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_WORLD - synchronized standard
85b4fc646aSLois Curfman McInnes $       output where only the first processor opens
86b4fc646aSLois Curfman McInnes $       the file.  All other processors send their
87b4fc646aSLois Curfman McInnes $       data to the first processor to print.
88b4fc646aSLois Curfman McInnes $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
89bbf0e169SBarry Smith 
90639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
91005c665bSBarry Smith 
92b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
93bbf0e169SBarry Smith @*/
94b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
95bbf0e169SBarry Smith {
96005c665bSBarry Smith   ViewerType vtype;
97639f9d9dSBarry Smith   int        i,j,format,ierr;
98b4fc646aSLois Curfman McInnes   FILE       *fd;
99bbf0e169SBarry Smith 
100b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
101b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
102b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
103bbf0e169SBarry Smith 
104005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
105005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
106b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
107005c665bSBarry Smith     return 0;
108005c665bSBarry Smith   }
109b4fc646aSLois Curfman McInnes   else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
110b4fc646aSLois Curfman McInnes     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
111b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
112b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
113b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
114b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
115ae09f205SBarry Smith 
116ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
117ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
118b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
119b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
120b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
121b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
122b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
123639f9d9dSBarry Smith         }
124b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
125b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
126b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
127b4fc646aSLois Curfman McInnes         }
128bbf0e169SBarry Smith       }
129bbf0e169SBarry Smith     }
130bbf0e169SBarry Smith   }
131639f9d9dSBarry Smith   return 0;
132639f9d9dSBarry Smith }
133639f9d9dSBarry Smith 
1345615d1e5SSatish Balay #undef __FUNC__
1355615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
136639f9d9dSBarry Smith /*@
137b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
138b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
139639f9d9dSBarry Smith 
140639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
141639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
142639f9d9dSBarry Smith $          = error_rel*umin                    else
143639f9d9dSBarry Smith $
144639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
145639f9d9dSBarry Smith 
146639f9d9dSBarry Smith    Input Parameters:
147b4fc646aSLois Curfman McInnes .  coloring - the coloring context
148639f9d9dSBarry Smith .  error_rel - relative error
149639f9d9dSBarry Smith .  umin - minimum allowable u-value
150639f9d9dSBarry Smith 
151b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
152b4fc646aSLois Curfman McInnes 
153b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
154639f9d9dSBarry Smith @*/
155639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
156639f9d9dSBarry Smith {
157639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
158639f9d9dSBarry Smith 
159639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
160639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
161639f9d9dSBarry Smith   return 0;
162639f9d9dSBarry Smith }
163639f9d9dSBarry Smith 
1645615d1e5SSatish Balay #undef __FUNC__
165005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
166005c665bSBarry Smith /*@
167e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
168e0907662SLois Curfman McInnes    matrices.
169005c665bSBarry Smith 
170005c665bSBarry Smith    Input Parameters:
171b4fc646aSLois Curfman McInnes .  coloring - the coloring context
172e0907662SLois Curfman McInnes .  freq - frequency (default is 1)
173e0907662SLois Curfman McInnes 
174e0907662SLois Curfman McInnes    Notes:
175e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
176e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
177e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
178e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
179e0907662SLois Curfman McInnes 
180e0907662SLois Curfman McInnes    Options Database Keys:
181e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq>
182005c665bSBarry Smith 
183b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
184005c665bSBarry Smith @*/
185005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
186005c665bSBarry Smith {
187005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
188005c665bSBarry Smith 
189005c665bSBarry Smith   matfd->freq = freq;
190005c665bSBarry Smith   return 0;
191005c665bSBarry Smith }
192005c665bSBarry Smith 
193005c665bSBarry Smith #undef __FUNC__
194ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
195ff0cfa39SBarry Smith /*@
196ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
197ff0cfa39SBarry Smith    matrices.
198ff0cfa39SBarry Smith 
199ff0cfa39SBarry Smith    Input Parameters:
200ff0cfa39SBarry Smith .  coloring - the coloring context
201ff0cfa39SBarry Smith 
202ff0cfa39SBarry Smith    Output Parameters:
203ff0cfa39SBarry Smith .  freq - frequency (default is 1)
204ff0cfa39SBarry Smith 
205ff0cfa39SBarry Smith    Notes:
206ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
207ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
208ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
209ff0cfa39SBarry Smith    <freq> nonlinear iterations.
210ff0cfa39SBarry Smith 
211ff0cfa39SBarry Smith    Options Database Keys:
212ff0cfa39SBarry Smith $  -mat_fd_coloring_freq <freq>
213ff0cfa39SBarry Smith 
214ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
215ff0cfa39SBarry Smith @*/
216ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
217ff0cfa39SBarry Smith {
218ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
219ff0cfa39SBarry Smith 
220ff0cfa39SBarry Smith   *freq = matfd->freq;
221ff0cfa39SBarry Smith   return 0;
222ff0cfa39SBarry Smith }
223ff0cfa39SBarry Smith 
224ff0cfa39SBarry Smith #undef __FUNC__
225005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
226005c665bSBarry Smith /*@
227005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
228005c665bSBarry Smith 
229005c665bSBarry Smith    Input Parameters:
230b4fc646aSLois Curfman McInnes .  coloring - the coloring context
231005c665bSBarry Smith .  f - the function
232005c665bSBarry Smith .  fctx - the function context
233005c665bSBarry Smith 
234b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
235005c665bSBarry Smith @*/
236*840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
237005c665bSBarry Smith {
238005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
239005c665bSBarry Smith 
240005c665bSBarry Smith   matfd->f    = f;
241005c665bSBarry Smith   matfd->fctx = fctx;
242005c665bSBarry Smith 
243005c665bSBarry Smith   return 0;
244005c665bSBarry Smith }
245005c665bSBarry Smith 
246005c665bSBarry Smith #undef __FUNC__
247d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
248639f9d9dSBarry Smith /*@
249b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
250639f9d9dSBarry Smith    the options database.
251639f9d9dSBarry Smith 
252b4fc646aSLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
253639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
254639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
255ae09f205SBarry Smith $          = error_rel*umin                      else
256639f9d9dSBarry Smith $
257639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
258639f9d9dSBarry Smith 
259639f9d9dSBarry Smith    Input Parameters:
260b4fc646aSLois Curfman McInnes .  coloring - the coloring context
261639f9d9dSBarry Smith 
262b4fc646aSLois Curfman McInnes    Options Database Keys:
263b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_error <err>, where <err> is the square root
264b4fc646aSLois Curfman McInnes $           of relative error in the function
265b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_umin  <umin>, where umin is described above
266e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
267e0907662SLois Curfman McInnes $           computing a new Jacobian
268ca71c51bSBarry Smith $  -mat_fd_coloring_view
269ca71c51bSBarry Smith $  -mat_fd_coloring_view_info
270ca71c51bSBarry Smith $  -mat_fd_coloring_view_draw
271639f9d9dSBarry Smith 
272b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
273639f9d9dSBarry Smith @*/
274639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
275639f9d9dSBarry Smith {
276005c665bSBarry Smith   int    ierr,flag,freq = 1;
277639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
278639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
279639f9d9dSBarry Smith 
280639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
281639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
282639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
283005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
284005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
285005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
286639f9d9dSBarry Smith   if (flag) {
287639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
288639f9d9dSBarry Smith   }
289639f9d9dSBarry Smith   return 0;
290639f9d9dSBarry Smith }
291639f9d9dSBarry Smith 
2925615d1e5SSatish Balay #undef __FUNC__
293d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
294639f9d9dSBarry Smith /*@
295639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
296639f9d9dSBarry Smith     using coloring.
297639f9d9dSBarry Smith 
298639f9d9dSBarry Smith     Input Parameter:
299639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
300639f9d9dSBarry Smith 
301639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
302639f9d9dSBarry Smith @*/
303639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
304639f9d9dSBarry Smith {
305639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
306639f9d9dSBarry Smith 
307e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
308e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
309e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
310005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
311005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
312ae09f205SBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n");
313005c665bSBarry Smith   return 0;
314005c665bSBarry Smith }
315005c665bSBarry Smith 
316005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
317005c665bSBarry Smith {
318005c665bSBarry Smith   int ierr,flg;
319005c665bSBarry Smith 
320005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
321005c665bSBarry Smith   if (flg) {
322f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
323005c665bSBarry Smith   }
324ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
325ae09f205SBarry Smith   if (flg) {
326f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
327f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
328f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
329ae09f205SBarry Smith   }
330005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
331005c665bSBarry Smith   if (flg) {
332005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
333005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
334005c665bSBarry Smith   }
335bbf0e169SBarry Smith   return 0;
336bbf0e169SBarry Smith }
337bbf0e169SBarry Smith 
3385615d1e5SSatish Balay #undef __FUNC__
3395615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
340bbf0e169SBarry Smith /*@C
341639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
342639f9d9dSBarry Smith    computation of Jacobians.
343bbf0e169SBarry Smith 
344639f9d9dSBarry Smith    Input Parameters:
345639f9d9dSBarry Smith .  mat - the matrix containing the nonzero structure of the Jacobian
346639f9d9dSBarry Smith .  iscoloring - the coloring of the matrix
347bbf0e169SBarry Smith 
348bbf0e169SBarry Smith     Output Parameter:
349639f9d9dSBarry Smith .   color - the new coloring context
350bbf0e169SBarry Smith 
351b4fc646aSLois Curfman McInnes     Options Database Keys:
352b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view
353b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_draw
354b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_info
355639f9d9dSBarry Smith 
356639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
357bbf0e169SBarry Smith @*/
358639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
359bbf0e169SBarry Smith {
360639f9d9dSBarry Smith   MatFDColoring c;
361639f9d9dSBarry Smith   MPI_Comm      comm;
362639f9d9dSBarry Smith   int           ierr,M,N;
363639f9d9dSBarry Smith 
364639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
365e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
366639f9d9dSBarry Smith 
367639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
368d4bb536fSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
369639f9d9dSBarry Smith   PLogObjectCreate(c);
370639f9d9dSBarry Smith 
371639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
372639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
373639f9d9dSBarry Smith   } else {
374e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
375639f9d9dSBarry Smith   }
376639f9d9dSBarry Smith 
377639f9d9dSBarry Smith   c->error_rel = 1.e-8;
378ae09f205SBarry Smith   c->umin      = 1.e-6;
379005c665bSBarry Smith   c->freq      = 1;
380005c665bSBarry Smith 
381005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
382639f9d9dSBarry Smith 
383639f9d9dSBarry Smith   *color = c;
384639f9d9dSBarry Smith 
385bbf0e169SBarry Smith   return 0;
386bbf0e169SBarry Smith }
387bbf0e169SBarry Smith 
3885615d1e5SSatish Balay #undef __FUNC__
389d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
390bbf0e169SBarry Smith /*@C
391639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
392639f9d9dSBarry Smith     via MatFDColoringCreate().
393bbf0e169SBarry Smith 
394b4fc646aSLois Curfman McInnes     Input Parameter:
395639f9d9dSBarry Smith .   c - coloring context
396bbf0e169SBarry Smith 
397639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
398bbf0e169SBarry Smith @*/
399639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
400bbf0e169SBarry Smith {
401263760aaSBarry Smith   int i,ierr;
402bbf0e169SBarry Smith 
403d4bb536fSBarry Smith   if (--c->refct > 0) return 0;
404d4bb536fSBarry Smith 
405639f9d9dSBarry Smith 
406639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
407639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
408639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
409639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
410bbf0e169SBarry Smith   }
411639f9d9dSBarry Smith   PetscFree(c->ncolumns);
412639f9d9dSBarry Smith   PetscFree(c->columns);
413639f9d9dSBarry Smith   PetscFree(c->nrows);
414639f9d9dSBarry Smith   PetscFree(c->rows);
415639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
416639f9d9dSBarry Smith   PetscFree(c->scale);
417005c665bSBarry Smith   if (c->w1) {
418005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
419005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
420005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
421005c665bSBarry Smith   }
422639f9d9dSBarry Smith   PLogObjectDestroy(c);
423639f9d9dSBarry Smith   PetscHeaderDestroy(c);
424bbf0e169SBarry Smith   return 0;
425bbf0e169SBarry Smith }
42643a90d84SBarry Smith 
427005c665bSBarry Smith #include "snes.h"
428005c665bSBarry Smith 
4295615d1e5SSatish Balay #undef __FUNC__
4305615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
43143a90d84SBarry Smith /*@
432e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
433e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
43443a90d84SBarry Smith 
43543a90d84SBarry Smith     Input Parameters:
43643a90d84SBarry Smith .   mat - location to store Jacobian
43743a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
43843a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
439005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
44043a90d84SBarry Smith 
4418bba8e72SBarry Smith    Options Database Keys:
4428bba8e72SBarry Smith $  -mat_fd_coloring_freq <freq>
4438bba8e72SBarry Smith 
44443a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
44543a90d84SBarry Smith 
44643a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
44743a90d84SBarry Smith @*/
448005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
44943a90d84SBarry Smith {
450e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
45143a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
45243a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
45343a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
454005c665bSBarry Smith   Vec           w1,w2,w3;
455*840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
456005c665bSBarry Smith   void          *fctx = coloring->fctx;
457005c665bSBarry Smith 
458e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
459e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
460e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
461e0907662SLois Curfman McInnes 
462005c665bSBarry Smith 
463005c665bSBarry Smith   if (!coloring->w1) {
464005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
465005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
466005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
467005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
468005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
469005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
470005c665bSBarry Smith   }
471005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
47243a90d84SBarry Smith 
473e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
474e0907662SLois Curfman McInnes   if (fg) {
475e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
476e0907662SLois Curfman McInnes   } else {
47743a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
478e0907662SLois Curfman McInnes   }
47943a90d84SBarry Smith 
48043a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
48143a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
48243a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
48343a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
48443a90d84SBarry Smith 
48543a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
48643a90d84SBarry Smith   /*
48743a90d84SBarry Smith       Loop over each color
48843a90d84SBarry Smith   */
48943a90d84SBarry Smith 
49043a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
49143a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
49243a90d84SBarry Smith     /*
49343a90d84SBarry Smith        Loop over each column associated with color adding the
49443a90d84SBarry Smith        perturbation to the vector w3.
49543a90d84SBarry Smith     */
49643a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
49743a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
49843a90d84SBarry Smith       dx  = xx[col-start];
499ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
50043a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
501ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
502ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
50343a90d84SBarry Smith #else
504ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
505ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
50643a90d84SBarry Smith #endif
50743a90d84SBarry Smith       dx          *= epsilon;
50843a90d84SBarry Smith       wscale[col] = 1.0/dx;
50943a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
51043a90d84SBarry Smith     }
51143a90d84SBarry Smith     VecRestoreArray(x1,&xx);
51243a90d84SBarry Smith     /*
513e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
51443a90d84SBarry Smith     */
51543a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
51643a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
51743a90d84SBarry Smith     /* Communicate scale to all processors */
51843a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
51943a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
52043a90d84SBarry Smith #else
52143a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
52243a90d84SBarry Smith #endif
52343a90d84SBarry Smith     /*
524e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
52543a90d84SBarry Smith     */
52643a90d84SBarry Smith     VecGetArray(w2,&y);
52743a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
52843a90d84SBarry Smith       row    = coloring->rows[k][l];
52943a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
53043a90d84SBarry Smith       y[row] *= scale[col];
53143a90d84SBarry Smith       srow   = row + start;
53243a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
53343a90d84SBarry Smith     }
53443a90d84SBarry Smith     VecRestoreArray(w2,&y);
53543a90d84SBarry Smith   }
53643a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
53743a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
53843a90d84SBarry Smith   return 0;
53943a90d84SBarry Smith }
540*840b8ebdSBarry Smith 
541*840b8ebdSBarry Smith #include "ts.h"
542*840b8ebdSBarry Smith 
543*840b8ebdSBarry Smith #undef __FUNC__
544*840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
545*840b8ebdSBarry Smith /*@
546*840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
547*840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
548*840b8ebdSBarry Smith 
549*840b8ebdSBarry Smith     Input Parameters:
550*840b8ebdSBarry Smith .   mat - location to store Jacobian
551*840b8ebdSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
552*840b8ebdSBarry Smith .   x1 - location at which Jacobian is to be computed
553*840b8ebdSBarry Smith .   sctx - optional context required by function (actually a SNES context)
554*840b8ebdSBarry Smith 
555*840b8ebdSBarry Smith    Options Database Keys:
556*840b8ebdSBarry Smith $  -mat_fd_coloring_freq <freq>
557*840b8ebdSBarry Smith 
558*840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
559*840b8ebdSBarry Smith 
560*840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
561*840b8ebdSBarry Smith @*/
562*840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
563*840b8ebdSBarry Smith {
564*840b8ebdSBarry Smith   int           k,fg,ierr,N,start,end,l,row,col,srow;
565*840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
566*840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
567*840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
568*840b8ebdSBarry Smith   Vec           w1,w2,w3;
569*840b8ebdSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
570*840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
571*840b8ebdSBarry Smith 
572*840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
573*840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
574*840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
575*840b8ebdSBarry Smith 
576*840b8ebdSBarry Smith 
577*840b8ebdSBarry Smith   if (!coloring->w1) {
578*840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
579*840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
580*840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
581*840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
582*840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
583*840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
584*840b8ebdSBarry Smith   }
585*840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
586*840b8ebdSBarry Smith 
587*840b8ebdSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
588*840b8ebdSBarry Smith   if (fg) {
589*840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
590*840b8ebdSBarry Smith   } else {
591*840b8ebdSBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
592*840b8ebdSBarry Smith   }
593*840b8ebdSBarry Smith 
594*840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
595*840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
596*840b8ebdSBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
597*840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr);
598*840b8ebdSBarry Smith 
599*840b8ebdSBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
600*840b8ebdSBarry Smith   /*
601*840b8ebdSBarry Smith       Loop over each color
602*840b8ebdSBarry Smith   */
603*840b8ebdSBarry Smith 
604*840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
605*840b8ebdSBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
606*840b8ebdSBarry Smith     /*
607*840b8ebdSBarry Smith        Loop over each column associated with color adding the
608*840b8ebdSBarry Smith        perturbation to the vector w3.
609*840b8ebdSBarry Smith     */
610*840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
611*840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
612*840b8ebdSBarry Smith       dx  = xx[col-start];
613*840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
614*840b8ebdSBarry Smith #if !defined(PETSC_COMPLEX)
615*840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
616*840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
617*840b8ebdSBarry Smith #else
618*840b8ebdSBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
619*840b8ebdSBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
620*840b8ebdSBarry Smith #endif
621*840b8ebdSBarry Smith       dx          *= epsilon;
622*840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
623*840b8ebdSBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
624*840b8ebdSBarry Smith     }
625*840b8ebdSBarry Smith     VecRestoreArray(x1,&xx);
626*840b8ebdSBarry Smith     /*
627*840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
628*840b8ebdSBarry Smith     */
629*840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr);
630*840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
631*840b8ebdSBarry Smith     /* Communicate scale to all processors */
632*840b8ebdSBarry Smith #if !defined(PETSC_COMPLEX)
633*840b8ebdSBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
634*840b8ebdSBarry Smith #else
635*840b8ebdSBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
636*840b8ebdSBarry Smith #endif
637*840b8ebdSBarry Smith     /*
638*840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
639*840b8ebdSBarry Smith     */
640*840b8ebdSBarry Smith     VecGetArray(w2,&y);
641*840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
642*840b8ebdSBarry Smith       row    = coloring->rows[k][l];
643*840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
644*840b8ebdSBarry Smith       y[row] *= scale[col];
645*840b8ebdSBarry Smith       srow   = row + start;
646*840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
647*840b8ebdSBarry Smith     }
648*840b8ebdSBarry Smith     VecRestoreArray(w2,&y);
649*840b8ebdSBarry Smith   }
650*840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
651*840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
652*840b8ebdSBarry Smith   return 0;
653*840b8ebdSBarry Smith }
654