xref: /petsc/src/mat/matfd/fdmatrix.c (revision 433994e6f6ce70ad8bdf6278bd516e26e00757e0)
1840b8ebdSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*433994e6SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.48 1999/06/30 23:52:39 balay 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 
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 
253a40ed3dSBarry Smith   PetscFunctionBegin;
2677ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
273a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(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];
395cd90555SBarry Smith       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
40005c665bSBarry Smith     }
41005c665bSBarry Smith   }
422bdab257SBarry Smith   ierr = DrawSynchronizedFlush(draw);CHKERRQ(ierr);
43005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause);CHKERRQ(ierr);
443a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
455cd90555SBarry Smith   ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr);
465cd90555SBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr);
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];
625cd90555SBarry Smith         ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
63005c665bSBarry Smith       }
64005c665bSBarry Smith     }
655cd90555SBarry Smith     ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr);
665cd90555SBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr);
67005c665bSBarry Smith   }
68005c665bSBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionReturn(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 
77fee21e36SBarry Smith    Collective on MatFDColoring unless Viewer is VIEWER_STDOUT_SELF
78fee21e36SBarry Smith 
79ef5ee4d1SLois Curfman McInnes    Input  Parameters:
80ef5ee4d1SLois Curfman McInnes +  c - the coloring context
81ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
82ef5ee4d1SLois Curfman McInnes 
8315091d37SBarry Smith    Level: intermediate
8415091d37SBarry Smith 
85b4fc646aSLois Curfman McInnes    Notes:
86b4fc646aSLois Curfman McInnes    The available visualization contexts include
87ef5ee4d1SLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
88ef5ee4d1SLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
89ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
90ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
91ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
92c655490fSBarry Smith -     VIEWER_DRAW_WORLD - graphical display of nonzero structure
93bbf0e169SBarry Smith 
94639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
95005c665bSBarry Smith 
96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
97bbf0e169SBarry Smith @*/
98b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
99bbf0e169SBarry Smith {
100005c665bSBarry Smith   ViewerType vtype;
101639f9d9dSBarry Smith   int        i,j,format,ierr;
102bbf0e169SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
104b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
105b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
106b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
107bbf0e169SBarry Smith 
108005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
1093f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,DRAW_VIEWER)) {
110b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1113a40ed3dSBarry Smith     PetscFunctionReturn(0);
1123f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
113d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
114d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
115d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
116d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
117ae09f205SBarry Smith 
118ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
119ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
120b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
121d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
122d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
123b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
124d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
125639f9d9dSBarry Smith         }
126d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
127b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
128d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
129b4fc646aSLois Curfman McInnes         }
130bbf0e169SBarry Smith       }
131bbf0e169SBarry Smith     }
1325cd90555SBarry Smith   } else {
1335cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
134bbf0e169SBarry Smith   }
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
136639f9d9dSBarry Smith }
137639f9d9dSBarry Smith 
1385615d1e5SSatish Balay #undef __FUNC__
1395615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
140639f9d9dSBarry Smith /*@
141b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
142b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
143639f9d9dSBarry Smith 
144ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
145ef5ee4d1SLois Curfman McInnes 
146ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
147ef5ee4d1SLois Curfman McInnes .vb
14865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
149f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
150f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
151ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
152ef5ee4d1SLois Curfman McInnes .ve
153639f9d9dSBarry Smith 
154639f9d9dSBarry Smith    Input Parameters:
155ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
156639f9d9dSBarry Smith .  error_rel - relative error
157f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
158fee21e36SBarry Smith 
15915091d37SBarry Smith    Level: advanced
16015091d37SBarry Smith 
161b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
162b4fc646aSLois Curfman McInnes 
163b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
164639f9d9dSBarry Smith @*/
165639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
166639f9d9dSBarry Smith {
1673a40ed3dSBarry Smith   PetscFunctionBegin;
168639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
169639f9d9dSBarry Smith 
170639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
171639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
173639f9d9dSBarry Smith }
174639f9d9dSBarry Smith 
1755615d1e5SSatish Balay #undef __FUNC__
176005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
177005c665bSBarry Smith /*@
178e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
179e0907662SLois Curfman McInnes    matrices.
180005c665bSBarry Smith 
181fee21e36SBarry Smith    Collective on MatFDColoring
182fee21e36SBarry Smith 
183ef5ee4d1SLois Curfman McInnes    Input Parameters:
184ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
185ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
186ef5ee4d1SLois Curfman McInnes 
18715091d37SBarry Smith    Options Database Keys:
18815091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
18915091d37SBarry Smith 
19015091d37SBarry Smith    Level: advanced
19115091d37SBarry Smith 
192e0907662SLois Curfman McInnes    Notes:
193e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
194e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
195e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
196e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
197e0907662SLois Curfman McInnes 
198b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
199ef5ee4d1SLois Curfman McInnes 
200ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency()
201005c665bSBarry Smith @*/
202005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
203005c665bSBarry Smith {
2043a40ed3dSBarry Smith   PetscFunctionBegin;
205005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
206005c665bSBarry Smith 
207005c665bSBarry Smith   matfd->freq = freq;
2083a40ed3dSBarry Smith   PetscFunctionReturn(0);
209005c665bSBarry Smith }
210005c665bSBarry Smith 
211005c665bSBarry Smith #undef __FUNC__
212ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
213ff0cfa39SBarry Smith /*@
214ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
215ff0cfa39SBarry Smith    matrices.
216ff0cfa39SBarry Smith 
217ef5ee4d1SLois Curfman McInnes    Not Collective
218ef5ee4d1SLois Curfman McInnes 
219ff0cfa39SBarry Smith    Input Parameters:
220ff0cfa39SBarry Smith .  coloring - the coloring context
221ff0cfa39SBarry Smith 
222ff0cfa39SBarry Smith    Output Parameters:
223ff0cfa39SBarry Smith .  freq - frequency (default is 1)
224ff0cfa39SBarry Smith 
22515091d37SBarry Smith    Options Database Keys:
22615091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
22715091d37SBarry Smith 
22815091d37SBarry Smith    Level: advanced
22915091d37SBarry Smith 
230ff0cfa39SBarry Smith    Notes:
231ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
232ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
233ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
234ff0cfa39SBarry Smith    <freq> nonlinear iterations.
235ff0cfa39SBarry Smith 
236ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
237ef5ee4d1SLois Curfman McInnes 
238ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
239ff0cfa39SBarry Smith @*/
240ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
241ff0cfa39SBarry Smith {
2423a40ed3dSBarry Smith   PetscFunctionBegin;
243ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
244ff0cfa39SBarry Smith 
245ff0cfa39SBarry Smith   *freq = matfd->freq;
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
247ff0cfa39SBarry Smith }
248ff0cfa39SBarry Smith 
249ff0cfa39SBarry Smith #undef __FUNC__
250005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
251d64ed03dSBarry Smith /*@C
252005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
253005c665bSBarry Smith 
254fee21e36SBarry Smith    Collective on MatFDColoring
255fee21e36SBarry Smith 
256ef5ee4d1SLois Curfman McInnes    Input Parameters:
257ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
258ef5ee4d1SLois Curfman McInnes .  f - the function
259ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
260ef5ee4d1SLois Curfman McInnes 
26115091d37SBarry Smith    Level: intermediate
26215091d37SBarry Smith 
263f881d145SBarry Smith    Notes:
264f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
265f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
266f881d145SBarry Smith   with the TS solvers.
267f881d145SBarry Smith 
268b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
269005c665bSBarry Smith @*/
270840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
271005c665bSBarry Smith {
2723a40ed3dSBarry Smith   PetscFunctionBegin;
273005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
274005c665bSBarry Smith 
275005c665bSBarry Smith   matfd->f    = f;
276005c665bSBarry Smith   matfd->fctx = fctx;
277005c665bSBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionReturn(0);
279005c665bSBarry Smith }
280005c665bSBarry Smith 
281005c665bSBarry Smith #undef __FUNC__
282d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
283639f9d9dSBarry Smith /*@
284b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
285639f9d9dSBarry Smith    the options database.
286639f9d9dSBarry Smith 
287fee21e36SBarry Smith    Collective on MatFDColoring
288fee21e36SBarry Smith 
28965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
290ef5ee4d1SLois Curfman McInnes .vb
29165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
292f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
293f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
294ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
295ef5ee4d1SLois Curfman McInnes .ve
296ef5ee4d1SLois Curfman McInnes 
297ef5ee4d1SLois Curfman McInnes    Input Parameter:
298ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
299ef5ee4d1SLois Curfman McInnes 
300b4fc646aSLois Curfman McInnes    Options Database Keys:
301ef5ee4d1SLois Curfman McInnes +  -mat_fd_coloring_error <err> - Sets <err> (square root
302ef5ee4d1SLois Curfman McInnes            of relative error in the function)
303f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
304ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
305ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
306ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
307ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
308639f9d9dSBarry Smith 
30915091d37SBarry Smith     Level: intermediate
31015091d37SBarry Smith 
311b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
312639f9d9dSBarry Smith @*/
313639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
314639f9d9dSBarry Smith {
315005c665bSBarry Smith   int    ierr,flag,freq = 1;
316639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
3173a40ed3dSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
319639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
320639f9d9dSBarry Smith 
321639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
322639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
323639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin);CHKERRQ(ierr);
324005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
325005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
326005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag);CHKERRQ(ierr);
327639f9d9dSBarry Smith   if (flag) {
328639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd);CHKERRQ(ierr);
329639f9d9dSBarry Smith   }
3303a40ed3dSBarry Smith   PetscFunctionReturn(0);
331639f9d9dSBarry Smith }
332639f9d9dSBarry Smith 
3335615d1e5SSatish Balay #undef __FUNC__
334d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
335639f9d9dSBarry Smith /*@
336639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
337639f9d9dSBarry Smith     using coloring.
338639f9d9dSBarry Smith 
339ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
340ef5ee4d1SLois Curfman McInnes 
341639f9d9dSBarry Smith     Input Parameter:
342639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
343639f9d9dSBarry Smith 
34415091d37SBarry Smith     Level: intermediate
34515091d37SBarry Smith 
346639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
347639f9d9dSBarry Smith @*/
348639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
349639f9d9dSBarry Smith {
350d132466eSBarry Smith   int ierr;
351d132466eSBarry Smith 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
353639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
354639f9d9dSBarry Smith 
355d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);CHKERRQ(ierr);
356d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);CHKERRQ(ierr);
357d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);CHKERRQ(ierr);
358d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");CHKERRQ(ierr);
359d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");CHKERRQ(ierr);
360d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");CHKERRQ(ierr);
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
362005c665bSBarry Smith }
363005c665bSBarry Smith 
364*433994e6SBarry Smith #undef __FUNC__
365*433994e6SBarry Smith #define __FUNC__ "MatFDColoringView_Private"
366005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
367005c665bSBarry Smith {
368005c665bSBarry Smith   int ierr,flg;
369005c665bSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
372005c665bSBarry Smith   if (flg) {
373f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
374005c665bSBarry Smith   }
375ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
376ae09f205SBarry Smith   if (flg) {
377f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
378f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
379f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
380ae09f205SBarry Smith   }
381005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
382005c665bSBarry Smith   if (flg) {
383c655490fSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
384c655490fSBarry Smith     ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
385005c665bSBarry Smith   }
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
387bbf0e169SBarry Smith }
388bbf0e169SBarry Smith 
3895615d1e5SSatish Balay #undef __FUNC__
3905615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
391bbf0e169SBarry Smith /*@C
392639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
393639f9d9dSBarry Smith    computation of Jacobians.
394bbf0e169SBarry Smith 
395ef5ee4d1SLois Curfman McInnes    Collective on Mat
396ef5ee4d1SLois Curfman McInnes 
397639f9d9dSBarry Smith    Input Parameters:
398ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
399ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
400bbf0e169SBarry Smith 
401bbf0e169SBarry Smith     Output Parameter:
402639f9d9dSBarry Smith .   color - the new coloring context
403bbf0e169SBarry Smith 
404b4fc646aSLois Curfman McInnes     Options Database Keys:
405ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
406ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
407ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
408639f9d9dSBarry Smith 
40915091d37SBarry Smith     Level: intermediate
41015091d37SBarry Smith 
411639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
412bbf0e169SBarry Smith @*/
413639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
414bbf0e169SBarry Smith {
415639f9d9dSBarry Smith   MatFDColoring c;
416639f9d9dSBarry Smith   MPI_Comm      comm;
417639f9d9dSBarry Smith   int           ierr,M,N;
418639f9d9dSBarry Smith 
4193a40ed3dSBarry Smith   PetscFunctionBegin;
420639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
421e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
422639f9d9dSBarry Smith 
423f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4243f1db9ecSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,
4253f1db9ecSBarry Smith                     MatFDColoringDestroy,MatFDColoringView);
426639f9d9dSBarry Smith   PLogObjectCreate(c);
427639f9d9dSBarry Smith 
428f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
429f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
430639f9d9dSBarry Smith   } else {
431e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
432639f9d9dSBarry Smith   }
433639f9d9dSBarry Smith 
434639f9d9dSBarry Smith   c->error_rel = 1.e-8;
435ae09f205SBarry Smith   c->umin      = 1.e-6;
436005c665bSBarry Smith   c->freq      = 1;
437005c665bSBarry Smith 
438005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
439639f9d9dSBarry Smith 
440639f9d9dSBarry Smith   *color = c;
441639f9d9dSBarry Smith 
4423a40ed3dSBarry Smith   PetscFunctionReturn(0);
443bbf0e169SBarry Smith }
444bbf0e169SBarry Smith 
4455615d1e5SSatish Balay #undef __FUNC__
446d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
447bbf0e169SBarry Smith /*@C
448639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
449639f9d9dSBarry Smith     via MatFDColoringCreate().
450bbf0e169SBarry Smith 
451ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
452ef5ee4d1SLois Curfman McInnes 
453b4fc646aSLois Curfman McInnes     Input Parameter:
454639f9d9dSBarry Smith .   c - coloring context
455bbf0e169SBarry Smith 
45615091d37SBarry Smith     Level: intermediate
45715091d37SBarry Smith 
458639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
459bbf0e169SBarry Smith @*/
460639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
461bbf0e169SBarry Smith {
462263760aaSBarry Smith   int i,ierr;
463bbf0e169SBarry Smith 
4643a40ed3dSBarry Smith   PetscFunctionBegin;
4653a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
466d4bb536fSBarry Smith 
467639f9d9dSBarry Smith 
468639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
469606d414cSSatish Balay     if (c->columns[i])       {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
470606d414cSSatish Balay     if (c->rows[i])          {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
471606d414cSSatish Balay     if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
472bbf0e169SBarry Smith   }
473606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
474606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
475606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
476606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
477606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
478606d414cSSatish Balay   ierr = PetscFree(c->scale);CHKERRQ(ierr);
479005c665bSBarry Smith   if (c->w1) {
480005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
481005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
482005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
483005c665bSBarry Smith   }
484639f9d9dSBarry Smith   PLogObjectDestroy(c);
485639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4863a40ed3dSBarry Smith   PetscFunctionReturn(0);
487bbf0e169SBarry Smith }
48843a90d84SBarry Smith 
489005c665bSBarry Smith #include "snes.h"
490005c665bSBarry Smith 
4915615d1e5SSatish Balay #undef __FUNC__
4925615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
49343a90d84SBarry Smith /*@
494e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
495e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49643a90d84SBarry Smith 
497fee21e36SBarry Smith     Collective on MatFDColoring
498fee21e36SBarry Smith 
499ef5ee4d1SLois Curfman McInnes     Input Parameters:
500ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
501ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
502ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
503ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
504ef5ee4d1SLois Curfman McInnes 
5058bba8e72SBarry Smith    Options Database Keys:
506ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
5078bba8e72SBarry Smith 
50815091d37SBarry Smith    Level: intermediate
50915091d37SBarry Smith 
51043a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
51143a90d84SBarry Smith 
51243a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
51343a90d84SBarry Smith @*/
514005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51543a90d84SBarry Smith {
516e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
51743a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
51843a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
51943a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
520005c665bSBarry Smith   Vec           w1,w2,w3;
521840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
522005c665bSBarry Smith   void          *fctx = coloring->fctx;
523005c665bSBarry Smith 
5243a40ed3dSBarry Smith   PetscFunctionBegin;
525e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
526e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
527e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
528e0907662SLois Curfman McInnes 
529005c665bSBarry Smith 
530005c665bSBarry Smith   if (!coloring->w1) {
531005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
532005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
533005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
534005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
535005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
536005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
537005c665bSBarry Smith   }
538005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
53943a90d84SBarry Smith 
540e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg);CHKERRQ(ierr);
541e0907662SLois Curfman McInnes   if (fg) {
542e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
543e0907662SLois Curfman McInnes   } else {
54443a90d84SBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
545e0907662SLois Curfman McInnes   }
54643a90d84SBarry Smith 
54743a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
54843a90d84SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
54943a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
55043a90d84SBarry Smith 
551549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
55243a90d84SBarry Smith   /*
55343a90d84SBarry Smith       Loop over each color
55443a90d84SBarry Smith   */
55543a90d84SBarry Smith 
5563b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
55743a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
55843a90d84SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
55943a90d84SBarry Smith     /*
56043a90d84SBarry Smith        Loop over each column associated with color adding the
56143a90d84SBarry Smith        perturbation to the vector w3.
56243a90d84SBarry Smith     */
56343a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
56443a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
56543a90d84SBarry Smith       dx  = xx[col-start];
566ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
567aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
568ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
569ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
57043a90d84SBarry Smith #else
571e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
572e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
57343a90d84SBarry Smith #endif
57443a90d84SBarry Smith       dx          *= epsilon;
57543a90d84SBarry Smith       wscale[col] = 1.0/dx;
5763b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
57743a90d84SBarry Smith     }
5783b28642cSBarry Smith 
57943a90d84SBarry Smith     /*
580e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
58143a90d84SBarry Smith     */
58243a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
58343a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
58443a90d84SBarry Smith     /* Communicate scale to all processors */
585aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
586ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
58743a90d84SBarry Smith #else
588ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
58943a90d84SBarry Smith #endif
59043a90d84SBarry Smith     /*
591e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
59243a90d84SBarry Smith     */
5933b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
59443a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
59543a90d84SBarry Smith       row    = coloring->rows[k][l];
59643a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
59743a90d84SBarry Smith       y[row] *= scale[col];
59843a90d84SBarry Smith       srow   = row + start;
59943a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
60043a90d84SBarry Smith     }
6013b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
60243a90d84SBarry Smith   }
6033b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
60443a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60543a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6063a40ed3dSBarry Smith   PetscFunctionReturn(0);
60743a90d84SBarry Smith }
608840b8ebdSBarry Smith 
609840b8ebdSBarry Smith #include "ts.h"
610840b8ebdSBarry Smith 
611840b8ebdSBarry Smith #undef __FUNC__
612840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
613840b8ebdSBarry Smith /*@
614840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
615840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
616840b8ebdSBarry Smith 
617fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
618fee21e36SBarry Smith 
619ef5ee4d1SLois Curfman McInnes     Input Parameters:
6203b28642cSBarry Smith +   mat - location to store Jacobian
621ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
622ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
623ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
624ef5ee4d1SLois Curfman McInnes 
625840b8ebdSBarry Smith    Options Database Keys:
626ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
627840b8ebdSBarry Smith 
62815091d37SBarry Smith    Level: intermediate
62915091d37SBarry Smith 
630840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
631840b8ebdSBarry Smith 
632840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
633840b8ebdSBarry Smith @*/
634840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
635840b8ebdSBarry Smith {
636840b8ebdSBarry Smith   int           k,fg,ierr,N,start,end,l,row,col,srow;
637840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
638840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
639840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
640840b8ebdSBarry Smith   Vec           w1,w2,w3;
641840b8ebdSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
642840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
643840b8ebdSBarry Smith 
6443a40ed3dSBarry Smith   PetscFunctionBegin;
645840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
646840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
647840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
648840b8ebdSBarry Smith 
649840b8ebdSBarry Smith   if (!coloring->w1) {
650840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
651840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
652840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
653840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
654840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
655840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
656840b8ebdSBarry Smith   }
657840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
658840b8ebdSBarry Smith 
659840b8ebdSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg);CHKERRQ(ierr);
660840b8ebdSBarry Smith   if (fg) {
661840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
662840b8ebdSBarry Smith   } else {
663840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
664840b8ebdSBarry Smith   }
665840b8ebdSBarry Smith 
666840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
667840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
668840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
669840b8ebdSBarry Smith 
670549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
671840b8ebdSBarry Smith   /*
672840b8ebdSBarry Smith       Loop over each color
673840b8ebdSBarry Smith   */
674840b8ebdSBarry Smith 
6753b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
676840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
677840b8ebdSBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
678840b8ebdSBarry Smith     /*
679840b8ebdSBarry Smith        Loop over each column associated with color adding the
680840b8ebdSBarry Smith        perturbation to the vector w3.
681840b8ebdSBarry Smith     */
682840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
683840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
684840b8ebdSBarry Smith       dx  = xx[col-start];
685840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
686aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
687840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
688840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
689840b8ebdSBarry Smith #else
690e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
691e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
692840b8ebdSBarry Smith #endif
693840b8ebdSBarry Smith       dx          *= epsilon;
694840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
6953b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
696840b8ebdSBarry Smith     }
697840b8ebdSBarry Smith     /*
698840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
699840b8ebdSBarry Smith     */
700840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
701840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
702840b8ebdSBarry Smith     /* Communicate scale to all processors */
703aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
704ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
705840b8ebdSBarry Smith #else
706ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
707840b8ebdSBarry Smith #endif
708840b8ebdSBarry Smith     /*
709840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
710840b8ebdSBarry Smith     */
7113b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
712840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
713840b8ebdSBarry Smith       row    = coloring->rows[k][l];
714840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
715840b8ebdSBarry Smith       y[row] *= scale[col];
716840b8ebdSBarry Smith       srow   = row + start;
717840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
718840b8ebdSBarry Smith     }
7193b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
720840b8ebdSBarry Smith   }
7213b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
722840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
723840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7243a40ed3dSBarry Smith   PetscFunctionReturn(0);
725840b8ebdSBarry Smith }
7263b28642cSBarry Smith 
7273b28642cSBarry Smith 
7283b28642cSBarry Smith 
729