xref: /petsc/src/mat/matfd/fdmatrix.c (revision ff0cfa393765184a9047ba166a67ffec829df148)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*ff0cfa39SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.21 1997/10/10 20:56:30 bsmith Exp bsmith $";
3bbf0e169SBarry Smith #endif
4bbf0e169SBarry Smith 
5bbf0e169SBarry Smith /*
6639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
7639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
8bbf0e169SBarry Smith */
9bbf0e169SBarry Smith 
10bbf0e169SBarry Smith #include "petsc.h"
11bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
12bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
13bbf0e169SBarry Smith #include "pinclude/pviewer.h"
14bbf0e169SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw"
17005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
18005c665bSBarry Smith {
19005c665bSBarry Smith   int         ierr,i,j,pause;
20005c665bSBarry Smith   PetscTruth  isnull;
21005c665bSBarry Smith   Draw        draw;
22d4bb536fSBarry Smith   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
23005c665bSBarry Smith   DrawButton  button;
24005c665bSBarry Smith 
25005c665bSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
26005c665bSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
272bdab257SBarry Smith   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
28005c665bSBarry Smith 
29005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
30005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
31005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
32005c665bSBarry Smith 
33005c665bSBarry Smith   /* loop over colors  */
34005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
35005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
36005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
37005c665bSBarry Smith       x = fd->columnsforrow[i][j];
38005c665bSBarry Smith       DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
39005c665bSBarry Smith     }
40005c665bSBarry Smith   }
412bdab257SBarry Smith   ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr);
42005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr);
43005c665bSBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
44005c665bSBarry Smith   ierr = DrawCheckResizedWindow(draw);
452bdab257SBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
46005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
472bdab257SBarry Smith     ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
48005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
49005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
50005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
51005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
52005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
53005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
54005c665bSBarry Smith     w *= scale; h *= scale;
55005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
56005c665bSBarry Smith     /* loop over colors  */
57005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
58005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
59005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
60005c665bSBarry Smith         x = fd->columnsforrow[i][j];
61005c665bSBarry Smith         DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);
62005c665bSBarry Smith       }
63005c665bSBarry Smith     }
64005c665bSBarry Smith     ierr = DrawCheckResizedWindow(draw);
652bdab257SBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
66005c665bSBarry Smith   }
67005c665bSBarry Smith 
68005c665bSBarry Smith   return 0;
69005c665bSBarry Smith }
70005c665bSBarry Smith 
71005c665bSBarry Smith #undef __FUNC__
72d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView"
73bbf0e169SBarry Smith /*@C
74639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
75bbf0e169SBarry Smith 
76b4fc646aSLois Curfman McInnes    Input  Parameters:
77b4fc646aSLois Curfman McInnes .  c - the coloring context
78b4fc646aSLois Curfman McInnes .  viewer - visualization context
79b4fc646aSLois Curfman McInnes 
80b4fc646aSLois Curfman McInnes    Notes:
81b4fc646aSLois Curfman McInnes    The available visualization contexts include
82b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_SELF - standard output (default)
83b4fc646aSLois Curfman McInnes $     VIEWER_STDOUT_WORLD - synchronized standard
84b4fc646aSLois Curfman McInnes $       output where only the first processor opens
85b4fc646aSLois Curfman McInnes $       the file.  All other processors send their
86b4fc646aSLois Curfman McInnes $       data to the first processor to print.
87b4fc646aSLois Curfman McInnes $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
88bbf0e169SBarry Smith 
89639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
90005c665bSBarry Smith 
91b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
92bbf0e169SBarry Smith @*/
93b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
94bbf0e169SBarry Smith {
95005c665bSBarry Smith   ViewerType vtype;
96639f9d9dSBarry Smith   int        i,j,format,ierr;
97b4fc646aSLois Curfman McInnes   FILE       *fd;
98bbf0e169SBarry Smith 
99b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
100b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
101b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
102bbf0e169SBarry Smith 
103005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
104005c665bSBarry Smith   if (vtype == DRAW_VIEWER) {
105b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
106005c665bSBarry Smith     return 0;
107005c665bSBarry Smith   }
108b4fc646aSLois Curfman McInnes   else if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
109b4fc646aSLois Curfman McInnes     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
110b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
111b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
112b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
113b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
114ae09f205SBarry Smith 
115ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
116ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
117b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
118b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
119b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
120b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
121b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
122639f9d9dSBarry Smith         }
123b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
124b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
125b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
126b4fc646aSLois Curfman McInnes         }
127bbf0e169SBarry Smith       }
128bbf0e169SBarry Smith     }
129bbf0e169SBarry Smith   }
130639f9d9dSBarry Smith   return 0;
131639f9d9dSBarry Smith }
132639f9d9dSBarry Smith 
1335615d1e5SSatish Balay #undef __FUNC__
1345615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
135639f9d9dSBarry Smith /*@
136b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
137b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
138639f9d9dSBarry Smith 
139639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
140639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
141639f9d9dSBarry Smith $          = error_rel*umin                    else
142639f9d9dSBarry Smith $
143639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
144639f9d9dSBarry Smith 
145639f9d9dSBarry Smith    Input Parameters:
146b4fc646aSLois Curfman McInnes .  coloring - the coloring context
147639f9d9dSBarry Smith .  error_rel - relative error
148639f9d9dSBarry Smith .  umin - minimum allowable u-value
149639f9d9dSBarry Smith 
150b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
151b4fc646aSLois Curfman McInnes 
152b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
153639f9d9dSBarry Smith @*/
154639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
155639f9d9dSBarry Smith {
156639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
157639f9d9dSBarry Smith 
158639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
159639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
160639f9d9dSBarry Smith   return 0;
161639f9d9dSBarry Smith }
162639f9d9dSBarry Smith 
1635615d1e5SSatish Balay #undef __FUNC__
164005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
165005c665bSBarry Smith /*@
166e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
167e0907662SLois Curfman McInnes    matrices.
168005c665bSBarry Smith 
169005c665bSBarry Smith    Input Parameters:
170b4fc646aSLois Curfman McInnes .  coloring - the coloring context
171e0907662SLois Curfman McInnes .  freq - frequency (default is 1)
172e0907662SLois Curfman McInnes 
173e0907662SLois Curfman McInnes    Notes:
174e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
175e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
176e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
177e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
178e0907662SLois Curfman McInnes 
179e0907662SLois Curfman McInnes    Options Database Keys:
180e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq>
181005c665bSBarry Smith 
182b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
183005c665bSBarry Smith @*/
184005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
185005c665bSBarry Smith {
186005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
187005c665bSBarry Smith 
188005c665bSBarry Smith   matfd->freq = freq;
189005c665bSBarry Smith   return 0;
190005c665bSBarry Smith }
191005c665bSBarry Smith 
192005c665bSBarry Smith #undef __FUNC__
193*ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
194*ff0cfa39SBarry Smith /*@
195*ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
196*ff0cfa39SBarry Smith    matrices.
197*ff0cfa39SBarry Smith 
198*ff0cfa39SBarry Smith    Input Parameters:
199*ff0cfa39SBarry Smith .  coloring - the coloring context
200*ff0cfa39SBarry Smith 
201*ff0cfa39SBarry Smith    Output Parameters:
202*ff0cfa39SBarry Smith .  freq - frequency (default is 1)
203*ff0cfa39SBarry Smith 
204*ff0cfa39SBarry Smith    Notes:
205*ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
206*ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
207*ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
208*ff0cfa39SBarry Smith    <freq> nonlinear iterations.
209*ff0cfa39SBarry Smith 
210*ff0cfa39SBarry Smith    Options Database Keys:
211*ff0cfa39SBarry Smith $  -mat_fd_coloring_freq <freq>
212*ff0cfa39SBarry Smith 
213*ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
214*ff0cfa39SBarry Smith @*/
215*ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
216*ff0cfa39SBarry Smith {
217*ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
218*ff0cfa39SBarry Smith 
219*ff0cfa39SBarry Smith   *freq = matfd->freq;
220*ff0cfa39SBarry Smith   return 0;
221*ff0cfa39SBarry Smith }
222*ff0cfa39SBarry Smith 
223*ff0cfa39SBarry Smith #undef __FUNC__
224005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
225005c665bSBarry Smith /*@
226005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
227005c665bSBarry Smith 
228005c665bSBarry Smith    Input Parameters:
229b4fc646aSLois Curfman McInnes .  coloring - the coloring context
230005c665bSBarry Smith .  f - the function
231005c665bSBarry Smith .  fctx - the function context
232005c665bSBarry Smith 
233b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
234005c665bSBarry Smith @*/
235005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx)
236005c665bSBarry Smith {
237005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
238005c665bSBarry Smith 
239005c665bSBarry Smith   matfd->f    = f;
240005c665bSBarry Smith   matfd->fctx = fctx;
241005c665bSBarry Smith 
242005c665bSBarry Smith   return 0;
243005c665bSBarry Smith }
244005c665bSBarry Smith 
245005c665bSBarry Smith #undef __FUNC__
246d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
247639f9d9dSBarry Smith /*@
248b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
249639f9d9dSBarry Smith    the options database.
250639f9d9dSBarry Smith 
251b4fc646aSLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
252639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
253639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
254ae09f205SBarry Smith $          = error_rel*umin                      else
255639f9d9dSBarry Smith $
256639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
257639f9d9dSBarry Smith 
258639f9d9dSBarry Smith    Input Parameters:
259b4fc646aSLois Curfman McInnes .  coloring - the coloring context
260639f9d9dSBarry Smith 
261b4fc646aSLois Curfman McInnes    Options Database Keys:
262b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_error <err>, where <err> is the square root
263b4fc646aSLois Curfman McInnes $           of relative error in the function
264b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_umin  <umin>, where umin is described above
265e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
266e0907662SLois Curfman McInnes $           computing a new Jacobian
267ca71c51bSBarry Smith $  -mat_fd_coloring_view
268ca71c51bSBarry Smith $  -mat_fd_coloring_view_info
269ca71c51bSBarry Smith $  -mat_fd_coloring_view_draw
270639f9d9dSBarry Smith 
271b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
272639f9d9dSBarry Smith @*/
273639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
274639f9d9dSBarry Smith {
275005c665bSBarry Smith   int    ierr,flag,freq = 1;
276639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
277639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
278639f9d9dSBarry Smith 
279639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
280639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
281639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
282005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
283005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
284005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
285639f9d9dSBarry Smith   if (flag) {
286639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
287639f9d9dSBarry Smith   }
288639f9d9dSBarry Smith   return 0;
289639f9d9dSBarry Smith }
290639f9d9dSBarry Smith 
2915615d1e5SSatish Balay #undef __FUNC__
292d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
293639f9d9dSBarry Smith /*@
294639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
295639f9d9dSBarry Smith     using coloring.
296639f9d9dSBarry Smith 
297639f9d9dSBarry Smith     Input Parameter:
298639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
299639f9d9dSBarry Smith 
300639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
301639f9d9dSBarry Smith @*/
302639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
303639f9d9dSBarry Smith {
304639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
305639f9d9dSBarry Smith 
306e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
307e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
308e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
309005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
310005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
311ae09f205SBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n");
312005c665bSBarry Smith   return 0;
313005c665bSBarry Smith }
314005c665bSBarry Smith 
315005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
316005c665bSBarry Smith {
317005c665bSBarry Smith   int ierr,flg;
318005c665bSBarry Smith 
319005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
320005c665bSBarry Smith   if (flg) {
321f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
322005c665bSBarry Smith   }
323ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
324ae09f205SBarry Smith   if (flg) {
325f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
326f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
327f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
328ae09f205SBarry Smith   }
329005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
330005c665bSBarry Smith   if (flg) {
331005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
332005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
333005c665bSBarry Smith   }
334bbf0e169SBarry Smith   return 0;
335bbf0e169SBarry Smith }
336bbf0e169SBarry Smith 
3375615d1e5SSatish Balay #undef __FUNC__
3385615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
339bbf0e169SBarry Smith /*@C
340639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
341639f9d9dSBarry Smith    computation of Jacobians.
342bbf0e169SBarry Smith 
343639f9d9dSBarry Smith    Input Parameters:
344639f9d9dSBarry Smith .  mat - the matrix containing the nonzero structure of the Jacobian
345639f9d9dSBarry Smith .  iscoloring - the coloring of the matrix
346bbf0e169SBarry Smith 
347bbf0e169SBarry Smith     Output Parameter:
348639f9d9dSBarry Smith .   color - the new coloring context
349bbf0e169SBarry Smith 
350b4fc646aSLois Curfman McInnes     Options Database Keys:
351b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view
352b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_draw
353b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_info
354639f9d9dSBarry Smith 
355639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
356bbf0e169SBarry Smith @*/
357639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
358bbf0e169SBarry Smith {
359639f9d9dSBarry Smith   MatFDColoring c;
360639f9d9dSBarry Smith   MPI_Comm      comm;
361639f9d9dSBarry Smith   int           ierr,M,N;
362639f9d9dSBarry Smith 
363639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
364e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
365639f9d9dSBarry Smith 
366639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
367d4bb536fSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
368639f9d9dSBarry Smith   PLogObjectCreate(c);
369639f9d9dSBarry Smith 
370639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
371639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
372639f9d9dSBarry Smith   } else {
373e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
374639f9d9dSBarry Smith   }
375639f9d9dSBarry Smith 
376639f9d9dSBarry Smith   c->error_rel = 1.e-8;
377ae09f205SBarry Smith   c->umin      = 1.e-6;
378005c665bSBarry Smith   c->freq      = 1;
379005c665bSBarry Smith 
380005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
381639f9d9dSBarry Smith 
382639f9d9dSBarry Smith   *color = c;
383639f9d9dSBarry Smith 
384bbf0e169SBarry Smith   return 0;
385bbf0e169SBarry Smith }
386bbf0e169SBarry Smith 
3875615d1e5SSatish Balay #undef __FUNC__
388d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
389bbf0e169SBarry Smith /*@C
390639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
391639f9d9dSBarry Smith     via MatFDColoringCreate().
392bbf0e169SBarry Smith 
393b4fc646aSLois Curfman McInnes     Input Parameter:
394639f9d9dSBarry Smith .   c - coloring context
395bbf0e169SBarry Smith 
396639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
397bbf0e169SBarry Smith @*/
398639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
399bbf0e169SBarry Smith {
400263760aaSBarry Smith   int i,ierr;
401bbf0e169SBarry Smith 
402d4bb536fSBarry Smith   if (--c->refct > 0) return 0;
403d4bb536fSBarry Smith 
404639f9d9dSBarry Smith 
405639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
406639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
407639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
408639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
409bbf0e169SBarry Smith   }
410639f9d9dSBarry Smith   PetscFree(c->ncolumns);
411639f9d9dSBarry Smith   PetscFree(c->columns);
412639f9d9dSBarry Smith   PetscFree(c->nrows);
413639f9d9dSBarry Smith   PetscFree(c->rows);
414639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
415639f9d9dSBarry Smith   PetscFree(c->scale);
416005c665bSBarry Smith   if (c->w1) {
417005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
418005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
419005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
420005c665bSBarry Smith   }
421639f9d9dSBarry Smith   PLogObjectDestroy(c);
422639f9d9dSBarry Smith   PetscHeaderDestroy(c);
423bbf0e169SBarry Smith   return 0;
424bbf0e169SBarry Smith }
42543a90d84SBarry Smith 
426005c665bSBarry Smith #include "snes.h"
427005c665bSBarry Smith 
4285615d1e5SSatish Balay #undef __FUNC__
4295615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
43043a90d84SBarry Smith /*@
431e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
432e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
43343a90d84SBarry Smith 
43443a90d84SBarry Smith     Input Parameters:
43543a90d84SBarry Smith .   mat - location to store Jacobian
43643a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
43743a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
438005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
43943a90d84SBarry Smith 
44043a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
44143a90d84SBarry Smith 
44243a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
44343a90d84SBarry Smith @*/
444005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
44543a90d84SBarry Smith {
446e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
44743a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
44843a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
44943a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
450005c665bSBarry Smith   Vec           w1,w2,w3;
451dd9fa9a5SLois Curfman McInnes   int           (*f)(void *,Vec,Vec,void *) = coloring->f;
452005c665bSBarry Smith   void          *fctx = coloring->fctx;
453005c665bSBarry Smith 
454e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
455e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
456e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
457e0907662SLois Curfman McInnes 
458d4bb536fSBarry Smith   /*
459005c665bSBarry Smith   ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr);
460005c665bSBarry Smith   if ((freq > 1) && ((it % freq) != 1)) {
461005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq);
462005c665bSBarry Smith     *flag = SAME_PRECONDITIONER;
463005c665bSBarry Smith     return 0;
464005c665bSBarry Smith   } else {
465005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq);
466005c665bSBarry Smith     *flag = SAME_NONZERO_PATTERN;
467d4bb536fSBarry Smith   }*/
468005c665bSBarry Smith 
469005c665bSBarry Smith   if (!coloring->w1) {
470005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
471005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
472005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
473005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
474005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
475005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
476005c665bSBarry Smith   }
477005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
47843a90d84SBarry Smith 
479e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
480e0907662SLois Curfman McInnes   if (fg) {
481e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
482e0907662SLois Curfman McInnes   } else {
48343a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
484e0907662SLois Curfman McInnes   }
48543a90d84SBarry Smith 
48643a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
48743a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
48843a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
48943a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
49043a90d84SBarry Smith 
49143a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
49243a90d84SBarry Smith   /*
49343a90d84SBarry Smith       Loop over each color
49443a90d84SBarry Smith   */
49543a90d84SBarry Smith 
49643a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
49743a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
49843a90d84SBarry Smith     /*
49943a90d84SBarry Smith        Loop over each column associated with color adding the
50043a90d84SBarry Smith        perturbation to the vector w3.
50143a90d84SBarry Smith     */
50243a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
50343a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
50443a90d84SBarry Smith       dx  = xx[col-start];
505ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
50643a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
507ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
508ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
50943a90d84SBarry Smith #else
510ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
511ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
51243a90d84SBarry Smith #endif
51343a90d84SBarry Smith       dx          *= epsilon;
51443a90d84SBarry Smith       wscale[col] = 1.0/dx;
51543a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
51643a90d84SBarry Smith     }
51743a90d84SBarry Smith     VecRestoreArray(x1,&xx);
51843a90d84SBarry Smith     /*
519e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
52043a90d84SBarry Smith     */
52143a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
52243a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
52343a90d84SBarry Smith     /* Communicate scale to all processors */
52443a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
52543a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
52643a90d84SBarry Smith #else
52743a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
52843a90d84SBarry Smith #endif
52943a90d84SBarry Smith     /*
530e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
53143a90d84SBarry Smith     */
53243a90d84SBarry Smith     VecGetArray(w2,&y);
53343a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
53443a90d84SBarry Smith       row    = coloring->rows[k][l];
53543a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
53643a90d84SBarry Smith       y[row] *= scale[col];
53743a90d84SBarry Smith       srow   = row + start;
53843a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
53943a90d84SBarry Smith     }
54043a90d84SBarry Smith     VecRestoreArray(w2,&y);
54143a90d84SBarry Smith   }
54243a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
54343a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
54443a90d84SBarry Smith   return 0;
54543a90d84SBarry Smith }
546