xref: /petsc/src/mat/matfd/fdmatrix.c (revision 263760aa02c8d9eb263ebefc569e51f63f8a9ae2)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*263760aaSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.20 1997/10/10 04:03:21 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__
193005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
194005c665bSBarry Smith /*@
195005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
196005c665bSBarry Smith 
197005c665bSBarry Smith    Input Parameters:
198b4fc646aSLois Curfman McInnes .  coloring - the coloring context
199005c665bSBarry Smith .  f - the function
200005c665bSBarry Smith .  fctx - the function context
201005c665bSBarry Smith 
202b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
203005c665bSBarry Smith @*/
204005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx)
205005c665bSBarry Smith {
206005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
207005c665bSBarry Smith 
208005c665bSBarry Smith   matfd->f    = f;
209005c665bSBarry Smith   matfd->fctx = fctx;
210005c665bSBarry Smith 
211005c665bSBarry Smith   return 0;
212005c665bSBarry Smith }
213005c665bSBarry Smith 
214005c665bSBarry Smith #undef __FUNC__
215d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
216639f9d9dSBarry Smith /*@
217b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
218639f9d9dSBarry Smith    the options database.
219639f9d9dSBarry Smith 
220b4fc646aSLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
221639f9d9dSBarry Smith $       J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
222639f9d9dSBarry Smith $        h = error_rel*u[i]                    if  u[i] > umin
223ae09f205SBarry Smith $          = error_rel*umin                      else
224639f9d9dSBarry Smith $
225639f9d9dSBarry Smith $   dx_{i} = (0, ... 1, .... 0)
226639f9d9dSBarry Smith 
227639f9d9dSBarry Smith    Input Parameters:
228b4fc646aSLois Curfman McInnes .  coloring - the coloring context
229639f9d9dSBarry Smith 
230b4fc646aSLois Curfman McInnes    Options Database Keys:
231b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_error <err>, where <err> is the square root
232b4fc646aSLois Curfman McInnes $           of relative error in the function
233b4fc646aSLois Curfman McInnes $  -mat_fd_coloring_umin  <umin>, where umin is described above
234e0907662SLois Curfman McInnes $  -mat_fd_coloring_freq <freq> where <freq> is the frequency of
235e0907662SLois Curfman McInnes $           computing a new Jacobian
236ca71c51bSBarry Smith $  -mat_fd_coloring_view
237ca71c51bSBarry Smith $  -mat_fd_coloring_view_info
238ca71c51bSBarry Smith $  -mat_fd_coloring_view_draw
239639f9d9dSBarry Smith 
240b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
241639f9d9dSBarry Smith @*/
242639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
243639f9d9dSBarry Smith {
244005c665bSBarry Smith   int    ierr,flag,freq = 1;
245639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
246639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
247639f9d9dSBarry Smith 
248639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
249639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
250639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
251005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
252005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
253005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
254639f9d9dSBarry Smith   if (flag) {
255639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
256639f9d9dSBarry Smith   }
257639f9d9dSBarry Smith   return 0;
258639f9d9dSBarry Smith }
259639f9d9dSBarry Smith 
2605615d1e5SSatish Balay #undef __FUNC__
261d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
262639f9d9dSBarry Smith /*@
263639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
264639f9d9dSBarry Smith     using coloring.
265639f9d9dSBarry Smith 
266639f9d9dSBarry Smith     Input Parameter:
267639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
268639f9d9dSBarry Smith 
269639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
270639f9d9dSBarry Smith @*/
271639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
272639f9d9dSBarry Smith {
273639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
274639f9d9dSBarry Smith 
275e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
276e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
277e0907662SLois Curfman McInnes   PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
278005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view\n");
279005c665bSBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n");
280ae09f205SBarry Smith   PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n");
281005c665bSBarry Smith   return 0;
282005c665bSBarry Smith }
283005c665bSBarry Smith 
284005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
285005c665bSBarry Smith {
286005c665bSBarry Smith   int ierr,flg;
287005c665bSBarry Smith 
288005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
289005c665bSBarry Smith   if (flg) {
290f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
291005c665bSBarry Smith   }
292ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
293ae09f205SBarry Smith   if (flg) {
294f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
295f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
296f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
297ae09f205SBarry Smith   }
298005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
299005c665bSBarry Smith   if (flg) {
300005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
301005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
302005c665bSBarry Smith   }
303bbf0e169SBarry Smith   return 0;
304bbf0e169SBarry Smith }
305bbf0e169SBarry Smith 
3065615d1e5SSatish Balay #undef __FUNC__
3075615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
308bbf0e169SBarry Smith /*@C
309639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
310639f9d9dSBarry Smith    computation of Jacobians.
311bbf0e169SBarry Smith 
312639f9d9dSBarry Smith    Input Parameters:
313639f9d9dSBarry Smith .  mat - the matrix containing the nonzero structure of the Jacobian
314639f9d9dSBarry Smith .  iscoloring - the coloring of the matrix
315bbf0e169SBarry Smith 
316bbf0e169SBarry Smith     Output Parameter:
317639f9d9dSBarry Smith .   color - the new coloring context
318bbf0e169SBarry Smith 
319b4fc646aSLois Curfman McInnes     Options Database Keys:
320b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view
321b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_draw
322b4fc646aSLois Curfman McInnes $    -mat_fd_coloring_view_info
323639f9d9dSBarry Smith 
324639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
325bbf0e169SBarry Smith @*/
326639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
327bbf0e169SBarry Smith {
328639f9d9dSBarry Smith   MatFDColoring c;
329639f9d9dSBarry Smith   MPI_Comm      comm;
330639f9d9dSBarry Smith   int           ierr,M,N;
331639f9d9dSBarry Smith 
332639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
333e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
334639f9d9dSBarry Smith 
335639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
336d4bb536fSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView);
337639f9d9dSBarry Smith   PLogObjectCreate(c);
338639f9d9dSBarry Smith 
339639f9d9dSBarry Smith   if (mat->ops.fdcoloringcreate) {
340639f9d9dSBarry Smith     ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
341639f9d9dSBarry Smith   } else {
342e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
343639f9d9dSBarry Smith   }
344639f9d9dSBarry Smith 
345639f9d9dSBarry Smith   c->error_rel = 1.e-8;
346ae09f205SBarry Smith   c->umin      = 1.e-6;
347005c665bSBarry Smith   c->freq      = 1;
348005c665bSBarry Smith 
349005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
350639f9d9dSBarry Smith 
351639f9d9dSBarry Smith   *color = c;
352639f9d9dSBarry Smith 
353bbf0e169SBarry Smith   return 0;
354bbf0e169SBarry Smith }
355bbf0e169SBarry Smith 
3565615d1e5SSatish Balay #undef __FUNC__
357d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
358bbf0e169SBarry Smith /*@C
359639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
360639f9d9dSBarry Smith     via MatFDColoringCreate().
361bbf0e169SBarry Smith 
362b4fc646aSLois Curfman McInnes     Input Parameter:
363639f9d9dSBarry Smith .   c - coloring context
364bbf0e169SBarry Smith 
365639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
366bbf0e169SBarry Smith @*/
367639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
368bbf0e169SBarry Smith {
369*263760aaSBarry Smith   int i,ierr;
370bbf0e169SBarry Smith 
371d4bb536fSBarry Smith   if (--c->refct > 0) return 0;
372d4bb536fSBarry Smith 
373639f9d9dSBarry Smith 
374639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
375639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
376639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
377639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
378bbf0e169SBarry Smith   }
379639f9d9dSBarry Smith   PetscFree(c->ncolumns);
380639f9d9dSBarry Smith   PetscFree(c->columns);
381639f9d9dSBarry Smith   PetscFree(c->nrows);
382639f9d9dSBarry Smith   PetscFree(c->rows);
383639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
384639f9d9dSBarry Smith   PetscFree(c->scale);
385005c665bSBarry Smith   if (c->w1) {
386005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
387005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
388005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
389005c665bSBarry Smith   }
390639f9d9dSBarry Smith   PLogObjectDestroy(c);
391639f9d9dSBarry Smith   PetscHeaderDestroy(c);
392bbf0e169SBarry Smith   return 0;
393bbf0e169SBarry Smith }
39443a90d84SBarry Smith 
395005c665bSBarry Smith #include "snes.h"
396005c665bSBarry Smith 
3975615d1e5SSatish Balay #undef __FUNC__
3985615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
39943a90d84SBarry Smith /*@
400e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
401e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
40243a90d84SBarry Smith 
40343a90d84SBarry Smith     Input Parameters:
40443a90d84SBarry Smith .   mat - location to store Jacobian
40543a90d84SBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
40643a90d84SBarry Smith .   x1 - location at which Jacobian is to be computed
407005c665bSBarry Smith .   sctx - optional context required by function (actually a SNES context)
40843a90d84SBarry Smith 
40943a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
41043a90d84SBarry Smith 
41143a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
41243a90d84SBarry Smith @*/
413005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
41443a90d84SBarry Smith {
415e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
41643a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
41743a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
41843a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
419005c665bSBarry Smith   Vec           w1,w2,w3;
420dd9fa9a5SLois Curfman McInnes   int           (*f)(void *,Vec,Vec,void *) = coloring->f;
421005c665bSBarry Smith   void          *fctx = coloring->fctx;
422005c665bSBarry Smith 
423e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
424e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
425e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
426e0907662SLois Curfman McInnes 
427d4bb536fSBarry Smith   /*
428005c665bSBarry Smith   ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr);
429005c665bSBarry Smith   if ((freq > 1) && ((it % freq) != 1)) {
430005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq);
431005c665bSBarry Smith     *flag = SAME_PRECONDITIONER;
432005c665bSBarry Smith     return 0;
433005c665bSBarry Smith   } else {
434005c665bSBarry Smith     PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq);
435005c665bSBarry Smith     *flag = SAME_NONZERO_PATTERN;
436d4bb536fSBarry Smith   }*/
437005c665bSBarry Smith 
438005c665bSBarry Smith   if (!coloring->w1) {
439005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
440005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
441005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
442005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
443005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
444005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
445005c665bSBarry Smith   }
446005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
44743a90d84SBarry Smith 
448e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
449e0907662SLois Curfman McInnes   if (fg) {
450e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
451e0907662SLois Curfman McInnes   } else {
45243a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
453e0907662SLois Curfman McInnes   }
45443a90d84SBarry Smith 
45543a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
45643a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
45743a90d84SBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
45843a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
45943a90d84SBarry Smith 
46043a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
46143a90d84SBarry Smith   /*
46243a90d84SBarry Smith       Loop over each color
46343a90d84SBarry Smith   */
46443a90d84SBarry Smith 
46543a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
46643a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
46743a90d84SBarry Smith     /*
46843a90d84SBarry Smith        Loop over each column associated with color adding the
46943a90d84SBarry Smith        perturbation to the vector w3.
47043a90d84SBarry Smith     */
47143a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
47243a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
47343a90d84SBarry Smith       dx  = xx[col-start];
474ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
47543a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
476ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
477ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
47843a90d84SBarry Smith #else
479ae09f205SBarry Smith       if (abs(dx) < umin && real(dx) >= 0.0)     dx = umin;
480ae09f205SBarry Smith       else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin;
48143a90d84SBarry Smith #endif
48243a90d84SBarry Smith       dx          *= epsilon;
48343a90d84SBarry Smith       wscale[col] = 1.0/dx;
48443a90d84SBarry Smith       VecSetValues(w3,1,&col,&dx,ADD_VALUES);
48543a90d84SBarry Smith     }
48643a90d84SBarry Smith     VecRestoreArray(x1,&xx);
48743a90d84SBarry Smith     /*
488e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
48943a90d84SBarry Smith     */
49043a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
49143a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
49243a90d84SBarry Smith     /* Communicate scale to all processors */
49343a90d84SBarry Smith #if !defined(PETSC_COMPLEX)
49443a90d84SBarry Smith     MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);
49543a90d84SBarry Smith #else
49643a90d84SBarry Smith     MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);
49743a90d84SBarry Smith #endif
49843a90d84SBarry Smith     /*
499e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
50043a90d84SBarry Smith     */
50143a90d84SBarry Smith     VecGetArray(w2,&y);
50243a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
50343a90d84SBarry Smith       row    = coloring->rows[k][l];
50443a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
50543a90d84SBarry Smith       y[row] *= scale[col];
50643a90d84SBarry Smith       srow   = row + start;
50743a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
50843a90d84SBarry Smith     }
50943a90d84SBarry Smith     VecRestoreArray(w2,&y);
51043a90d84SBarry Smith   }
51143a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
51243a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
51343a90d84SBarry Smith   return 0;
51443a90d84SBarry Smith }
515