xref: /petsc/src/mat/matfd/fdmatrix.c (revision f1af5d2ffeae1f5fc391a89939f4818e47770ae3)
1*f1af5d2fSBarry Smith /*$Id: fdmatrix.c,v 1.54 1999/10/24 14:02:53 bsmith Exp bsmith $*/
2bbf0e169SBarry Smith 
3bbf0e169SBarry Smith /*
4639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
5639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
6bbf0e169SBarry Smith */
7bbf0e169SBarry Smith 
8bbf0e169SBarry Smith #include "petsc.h"
9bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
10bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
11bbf0e169SBarry Smith 
125615d1e5SSatish Balay #undef __FUNC__
13d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw"
14005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
15005c665bSBarry Smith {
16005c665bSBarry Smith   int         ierr,i,j,pause;
17005c665bSBarry Smith   PetscTruth  isnull;
18005c665bSBarry Smith   Draw        draw;
19d4bb536fSBarry Smith   double      xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0;
20005c665bSBarry Smith   DrawButton  button;
21005c665bSBarry Smith 
223a40ed3dSBarry Smith   PetscFunctionBegin;
2377ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
243a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
252bdab257SBarry Smith   ierr = DrawSynchronizedClear(draw);CHKERRQ(ierr);
26005c665bSBarry Smith 
27005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
28005c665bSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
29005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
30005c665bSBarry Smith 
31005c665bSBarry Smith   /* loop over colors  */
32005c665bSBarry Smith   for (i=0; i<fd->ncolors; i++ ) {
33005c665bSBarry Smith     for ( j=0; j<fd->nrows[i]; j++ ) {
34005c665bSBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
35005c665bSBarry Smith       x = fd->columnsforrow[i][j];
365cd90555SBarry Smith       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
37005c665bSBarry Smith     }
38005c665bSBarry Smith   }
392bdab257SBarry Smith   ierr = DrawSynchronizedFlush(draw);CHKERRQ(ierr);
40005c665bSBarry Smith   ierr = DrawGetPause(draw,&pause);CHKERRQ(ierr);
413a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
425cd90555SBarry Smith   ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr);
435cd90555SBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr);
44005c665bSBarry Smith   while (button != BUTTON_RIGHT) {
452bdab257SBarry Smith     ierr = DrawSynchronizedClear(draw);CHKERRQ(ierr);
46005c665bSBarry Smith     if (button == BUTTON_LEFT) scale = .5;
47005c665bSBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
48005c665bSBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
49005c665bSBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
50005c665bSBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
51005c665bSBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
52005c665bSBarry Smith     w *= scale; h *= scale;
53005c665bSBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
54005c665bSBarry Smith     /* loop over colors  */
55005c665bSBarry Smith     for (i=0; i<fd->ncolors; i++ ) {
56005c665bSBarry Smith       for ( j=0; j<fd->nrows[i]; j++ ) {
57005c665bSBarry Smith         y = fd->M - fd->rows[i][j] - fd->rstart;
58005c665bSBarry Smith         x = fd->columnsforrow[i][j];
595cd90555SBarry Smith         ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
60005c665bSBarry Smith       }
61005c665bSBarry Smith     }
625cd90555SBarry Smith     ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr);
635cd90555SBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr);
64005c665bSBarry Smith   }
65005c665bSBarry Smith 
663a40ed3dSBarry Smith   PetscFunctionReturn(0);
67005c665bSBarry Smith }
68005c665bSBarry Smith 
69005c665bSBarry Smith #undef __FUNC__
70d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView"
71bbf0e169SBarry Smith /*@C
72639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
73bbf0e169SBarry Smith 
74fee21e36SBarry Smith    Collective on MatFDColoring unless Viewer is VIEWER_STDOUT_SELF
75fee21e36SBarry Smith 
76ef5ee4d1SLois Curfman McInnes    Input  Parameters:
77ef5ee4d1SLois Curfman McInnes +  c - the coloring context
78ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
79ef5ee4d1SLois Curfman McInnes 
8015091d37SBarry Smith    Level: intermediate
8115091d37SBarry Smith 
82b4fc646aSLois Curfman McInnes    Notes:
83b4fc646aSLois Curfman McInnes    The available visualization contexts include
84ef5ee4d1SLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
85ef5ee4d1SLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
86ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
87ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
88ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
89c655490fSBarry Smith -     VIEWER_DRAW_WORLD - graphical display of nonzero structure
90bbf0e169SBarry Smith 
91639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
92005c665bSBarry Smith 
93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
94bbf0e169SBarry Smith @*/
95b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
96bbf0e169SBarry Smith {
97639f9d9dSBarry Smith   int        i,j,format,ierr;
986831982aSBarry Smith   PetscTruth isdraw,isascii;
99bbf0e169SBarry Smith 
1003a40ed3dSBarry Smith   PetscFunctionBegin;
101b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
1020f5bd95cSBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_SELF;
1030f5bd95cSBarry Smith   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
1046831982aSBarry Smith   PetscCheckSameComm(c,viewer);
105bbf0e169SBarry Smith 
1066831982aSBarry Smith   ierr  = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
1076831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
1080f5bd95cSBarry Smith   if (isdraw) {
109b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1100f5bd95cSBarry Smith   } else if (isascii) {
111d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
112d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
113d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
114d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
115ae09f205SBarry Smith 
116ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
117ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
118b4fc646aSLois Curfman McInnes       for ( i=0; i<c->ncolors; i++ ) {
119d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
120d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
121b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
122d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
123639f9d9dSBarry Smith         }
124d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
125b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
126d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
127b4fc646aSLois Curfman McInnes         }
128bbf0e169SBarry Smith       }
129bbf0e169SBarry Smith     }
1306831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
1315cd90555SBarry Smith   } else {
1320f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
133bbf0e169SBarry Smith   }
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135639f9d9dSBarry Smith }
136639f9d9dSBarry Smith 
1375615d1e5SSatish Balay #undef __FUNC__
1385615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters"
139639f9d9dSBarry Smith /*@
140b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
141b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
142639f9d9dSBarry Smith 
143ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
144ef5ee4d1SLois Curfman McInnes 
145ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
146ef5ee4d1SLois Curfman McInnes .vb
14765f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
148f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
149f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
150ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
151ef5ee4d1SLois Curfman McInnes .ve
152639f9d9dSBarry Smith 
153639f9d9dSBarry Smith    Input Parameters:
154ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
155639f9d9dSBarry Smith .  error_rel - relative error
156f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
157fee21e36SBarry Smith 
15815091d37SBarry Smith    Level: advanced
15915091d37SBarry Smith 
160b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
161b4fc646aSLois Curfman McInnes 
162b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
163639f9d9dSBarry Smith @*/
164639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
165639f9d9dSBarry Smith {
1663a40ed3dSBarry Smith   PetscFunctionBegin;
167639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
168639f9d9dSBarry Smith 
169639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
170639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1713a40ed3dSBarry Smith   PetscFunctionReturn(0);
172639f9d9dSBarry Smith }
173639f9d9dSBarry Smith 
1745615d1e5SSatish Balay #undef __FUNC__
175005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
176005c665bSBarry Smith /*@
177e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
178e0907662SLois Curfman McInnes    matrices.
179005c665bSBarry Smith 
180fee21e36SBarry Smith    Collective on MatFDColoring
181fee21e36SBarry Smith 
182ef5ee4d1SLois Curfman McInnes    Input Parameters:
183ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
184ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
185ef5ee4d1SLois Curfman McInnes 
18615091d37SBarry Smith    Options Database Keys:
18715091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
18815091d37SBarry Smith 
18915091d37SBarry Smith    Level: advanced
19015091d37SBarry Smith 
191e0907662SLois Curfman McInnes    Notes:
192e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
193e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
194e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
195e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
196e0907662SLois Curfman McInnes 
197b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
198ef5ee4d1SLois Curfman McInnes 
199ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency()
200005c665bSBarry Smith @*/
201005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
202005c665bSBarry Smith {
2033a40ed3dSBarry Smith   PetscFunctionBegin;
204005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
205005c665bSBarry Smith 
206005c665bSBarry Smith   matfd->freq = freq;
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
208005c665bSBarry Smith }
209005c665bSBarry Smith 
210005c665bSBarry Smith #undef __FUNC__
211ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
212ff0cfa39SBarry Smith /*@
213ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
214ff0cfa39SBarry Smith    matrices.
215ff0cfa39SBarry Smith 
216ef5ee4d1SLois Curfman McInnes    Not Collective
217ef5ee4d1SLois Curfman McInnes 
218ff0cfa39SBarry Smith    Input Parameters:
219ff0cfa39SBarry Smith .  coloring - the coloring context
220ff0cfa39SBarry Smith 
221ff0cfa39SBarry Smith    Output Parameters:
222ff0cfa39SBarry Smith .  freq - frequency (default is 1)
223ff0cfa39SBarry Smith 
22415091d37SBarry Smith    Options Database Keys:
22515091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
22615091d37SBarry Smith 
22715091d37SBarry Smith    Level: advanced
22815091d37SBarry Smith 
229ff0cfa39SBarry Smith    Notes:
230ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
231ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
232ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
233ff0cfa39SBarry Smith    <freq> nonlinear iterations.
234ff0cfa39SBarry Smith 
235ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
236ef5ee4d1SLois Curfman McInnes 
237ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
238ff0cfa39SBarry Smith @*/
239ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
240ff0cfa39SBarry Smith {
2413a40ed3dSBarry Smith   PetscFunctionBegin;
242ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
243ff0cfa39SBarry Smith 
244ff0cfa39SBarry Smith   *freq = matfd->freq;
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
246ff0cfa39SBarry Smith }
247ff0cfa39SBarry Smith 
248ff0cfa39SBarry Smith #undef __FUNC__
249005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
250d64ed03dSBarry Smith /*@C
251005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
252005c665bSBarry Smith 
253fee21e36SBarry Smith    Collective on MatFDColoring
254fee21e36SBarry Smith 
255ef5ee4d1SLois Curfman McInnes    Input Parameters:
256ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
257ef5ee4d1SLois Curfman McInnes .  f - the function
258ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
259ef5ee4d1SLois Curfman McInnes 
26015091d37SBarry Smith    Level: intermediate
26115091d37SBarry Smith 
262f881d145SBarry Smith    Notes:
263f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
264f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
265f881d145SBarry Smith   with the TS solvers.
266f881d145SBarry Smith 
267b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
268005c665bSBarry Smith @*/
269840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
270005c665bSBarry Smith {
2713a40ed3dSBarry Smith   PetscFunctionBegin;
272005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
273005c665bSBarry Smith 
274005c665bSBarry Smith   matfd->f    = f;
275005c665bSBarry Smith   matfd->fctx = fctx;
276005c665bSBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionReturn(0);
278005c665bSBarry Smith }
279005c665bSBarry Smith 
280005c665bSBarry Smith #undef __FUNC__
281d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
282639f9d9dSBarry Smith /*@
283b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
284639f9d9dSBarry Smith    the options database.
285639f9d9dSBarry Smith 
286fee21e36SBarry Smith    Collective on MatFDColoring
287fee21e36SBarry Smith 
28865f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
289ef5ee4d1SLois Curfman McInnes .vb
29065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
291f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
292f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
293ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
294ef5ee4d1SLois Curfman McInnes .ve
295ef5ee4d1SLois Curfman McInnes 
296ef5ee4d1SLois Curfman McInnes    Input Parameter:
297ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
298ef5ee4d1SLois Curfman McInnes 
299b4fc646aSLois Curfman McInnes    Options Database Keys:
300ef5ee4d1SLois Curfman McInnes +  -mat_fd_coloring_error <err> - Sets <err> (square root
301ef5ee4d1SLois Curfman McInnes            of relative error in the function)
302f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
303ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
304ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
305ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
306ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
307639f9d9dSBarry Smith 
30815091d37SBarry Smith     Level: intermediate
30915091d37SBarry Smith 
310b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
311639f9d9dSBarry Smith @*/
312639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
313639f9d9dSBarry Smith {
314*f1af5d2fSBarry Smith   int        ierr,freq = 1;
315639f9d9dSBarry Smith   double     error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
316*f1af5d2fSBarry Smith   PetscTruth flag;
3173a40ed3dSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
319639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
320639f9d9dSBarry Smith 
321*f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,PETSC_NULL);CHKERRQ(ierr);
322*f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,PETSC_NULL);CHKERRQ(ierr);
323639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin);CHKERRQ(ierr);
324*f1af5d2fSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,PETSC_NULL);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 
364433994e6SBarry Smith #undef __FUNC__
365433994e6SBarry Smith #define __FUNC__ "MatFDColoringView_Private"
366005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
367005c665bSBarry Smith {
368*f1af5d2fSBarry Smith   int        ierr;
369*f1af5d2fSBarry Smith   PetscTruth flg;
370005c665bSBarry Smith 
3713a40ed3dSBarry Smith   PetscFunctionBegin;
372005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
373005c665bSBarry Smith   if (flg) {
374f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
375005c665bSBarry Smith   }
376ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
377ae09f205SBarry Smith   if (flg) {
378f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
379f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
380f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
381ae09f205SBarry Smith   }
382005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
383005c665bSBarry Smith   if (flg) {
384c655490fSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
385c655490fSBarry Smith     ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
386005c665bSBarry Smith   }
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
388bbf0e169SBarry Smith }
389bbf0e169SBarry Smith 
3905615d1e5SSatish Balay #undef __FUNC__
3915615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
392bbf0e169SBarry Smith /*@C
393639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
394639f9d9dSBarry Smith    computation of Jacobians.
395bbf0e169SBarry Smith 
396ef5ee4d1SLois Curfman McInnes    Collective on Mat
397ef5ee4d1SLois Curfman McInnes 
398639f9d9dSBarry Smith    Input Parameters:
399ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
400ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
401bbf0e169SBarry Smith 
402bbf0e169SBarry Smith     Output Parameter:
403639f9d9dSBarry Smith .   color - the new coloring context
404bbf0e169SBarry Smith 
405b4fc646aSLois Curfman McInnes     Options Database Keys:
406ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
407ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
408ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
409639f9d9dSBarry Smith 
41015091d37SBarry Smith     Level: intermediate
41115091d37SBarry Smith 
4126831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
4136831982aSBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
414bbf0e169SBarry Smith @*/
415639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
416bbf0e169SBarry Smith {
417639f9d9dSBarry Smith   MatFDColoring c;
418639f9d9dSBarry Smith   MPI_Comm      comm;
419639f9d9dSBarry Smith   int           ierr,M,N;
420639f9d9dSBarry Smith 
4213a40ed3dSBarry Smith   PetscFunctionBegin;
422639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
423e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
424639f9d9dSBarry Smith 
425f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4263f1db9ecSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,
4273f1db9ecSBarry Smith                     MatFDColoringDestroy,MatFDColoringView);
428639f9d9dSBarry Smith   PLogObjectCreate(c);
429639f9d9dSBarry Smith 
430f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
431f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
432639f9d9dSBarry Smith   } else {
433e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
434639f9d9dSBarry Smith   }
435639f9d9dSBarry Smith 
436639f9d9dSBarry Smith   c->error_rel = 1.e-8;
437ae09f205SBarry Smith   c->umin      = 1.e-6;
438005c665bSBarry Smith   c->freq      = 1;
439005c665bSBarry Smith 
440005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
441639f9d9dSBarry Smith 
442639f9d9dSBarry Smith   *color = c;
443639f9d9dSBarry Smith 
4443a40ed3dSBarry Smith   PetscFunctionReturn(0);
445bbf0e169SBarry Smith }
446bbf0e169SBarry Smith 
4475615d1e5SSatish Balay #undef __FUNC__
448d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
449bbf0e169SBarry Smith /*@C
450639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
451639f9d9dSBarry Smith     via MatFDColoringCreate().
452bbf0e169SBarry Smith 
453ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
454ef5ee4d1SLois Curfman McInnes 
455b4fc646aSLois Curfman McInnes     Input Parameter:
456639f9d9dSBarry Smith .   c - coloring context
457bbf0e169SBarry Smith 
45815091d37SBarry Smith     Level: intermediate
45915091d37SBarry Smith 
460639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
461bbf0e169SBarry Smith @*/
462639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
463bbf0e169SBarry Smith {
464263760aaSBarry Smith   int i,ierr;
465bbf0e169SBarry Smith 
4663a40ed3dSBarry Smith   PetscFunctionBegin;
4673a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
468d4bb536fSBarry Smith 
469639f9d9dSBarry Smith 
470639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
471606d414cSSatish Balay     if (c->columns[i])       {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
472606d414cSSatish Balay     if (c->rows[i])          {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
473606d414cSSatish Balay     if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
474bbf0e169SBarry Smith   }
475606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
476606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
477606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
478606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
479606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
480606d414cSSatish Balay   ierr = PetscFree(c->scale);CHKERRQ(ierr);
481005c665bSBarry Smith   if (c->w1) {
482005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
483005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
484005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
485005c665bSBarry Smith   }
486639f9d9dSBarry Smith   PLogObjectDestroy(c);
487639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4883a40ed3dSBarry Smith   PetscFunctionReturn(0);
489bbf0e169SBarry Smith }
49043a90d84SBarry Smith 
491005c665bSBarry Smith 
4925615d1e5SSatish Balay #undef __FUNC__
4935615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
49443a90d84SBarry Smith /*@
495e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
496e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49743a90d84SBarry Smith 
498fee21e36SBarry Smith     Collective on MatFDColoring
499fee21e36SBarry Smith 
500ef5ee4d1SLois Curfman McInnes     Input Parameters:
501ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
502ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
503ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
504ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
505ef5ee4d1SLois Curfman McInnes 
5068bba8e72SBarry Smith    Options Database Keys:
507ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
5088bba8e72SBarry Smith 
50915091d37SBarry Smith    Level: intermediate
51015091d37SBarry Smith 
51143a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
51243a90d84SBarry Smith 
51343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
51443a90d84SBarry Smith @*/
515005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51643a90d84SBarry Smith {
517*f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
51843a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
51943a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
52043a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
521005c665bSBarry Smith   Vec           w1,w2,w3;
522840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
523005c665bSBarry Smith   void          *fctx = coloring->fctx;
524*f1af5d2fSBarry Smith   PetscTruth    flg;
525005c665bSBarry Smith 
5263a40ed3dSBarry Smith   PetscFunctionBegin;
527e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
528e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
529e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
530e0907662SLois Curfman McInnes 
531005c665bSBarry Smith 
532005c665bSBarry Smith   if (!coloring->w1) {
533005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
534005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
535005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
536005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
537005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
538005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
539005c665bSBarry Smith   }
540005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
54143a90d84SBarry Smith 
542*f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
543*f1af5d2fSBarry Smith   if (flg) {
544e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
545e0907662SLois Curfman McInnes   } else {
54643a90d84SBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
547e0907662SLois Curfman McInnes   }
54843a90d84SBarry Smith 
54943a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
55043a90d84SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
55143a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
55243a90d84SBarry Smith 
553549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
55443a90d84SBarry Smith   /*
55543a90d84SBarry Smith       Loop over each color
55643a90d84SBarry Smith   */
55743a90d84SBarry Smith 
5583b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
55943a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
56043a90d84SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
56143a90d84SBarry Smith     /*
56243a90d84SBarry Smith        Loop over each column associated with color adding the
56343a90d84SBarry Smith        perturbation to the vector w3.
56443a90d84SBarry Smith     */
56543a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
56643a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
56743a90d84SBarry Smith       dx  = xx[col-start];
568ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
569aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
570ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
571ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
57243a90d84SBarry Smith #else
573e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
574e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
57543a90d84SBarry Smith #endif
57643a90d84SBarry Smith       dx          *= epsilon;
57743a90d84SBarry Smith       wscale[col] = 1.0/dx;
5783b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
57943a90d84SBarry Smith     }
5803b28642cSBarry Smith 
58143a90d84SBarry Smith     /*
582e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
58343a90d84SBarry Smith     */
58443a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
58543a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
58643a90d84SBarry Smith     /* Communicate scale to all processors */
5876831982aSBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
58843a90d84SBarry Smith     /*
589e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
59043a90d84SBarry Smith     */
5913b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
59243a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
59343a90d84SBarry Smith       row    = coloring->rows[k][l];
59443a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
59543a90d84SBarry Smith       y[row] *= scale[col];
59643a90d84SBarry Smith       srow   = row + start;
59743a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
59843a90d84SBarry Smith     }
5993b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
60043a90d84SBarry Smith   }
6013b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
60243a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60343a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6043a40ed3dSBarry Smith   PetscFunctionReturn(0);
60543a90d84SBarry Smith }
606840b8ebdSBarry Smith 
607840b8ebdSBarry Smith #undef __FUNC__
608840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
609840b8ebdSBarry Smith /*@
610840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
611840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
612840b8ebdSBarry Smith 
613fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
614fee21e36SBarry Smith 
615ef5ee4d1SLois Curfman McInnes     Input Parameters:
6163b28642cSBarry Smith +   mat - location to store Jacobian
617ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
618ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
619ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
620ef5ee4d1SLois Curfman McInnes 
621840b8ebdSBarry Smith    Options Database Keys:
622ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
623840b8ebdSBarry Smith 
62415091d37SBarry Smith    Level: intermediate
62515091d37SBarry Smith 
626840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
627840b8ebdSBarry Smith 
628840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
629840b8ebdSBarry Smith @*/
630840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
631840b8ebdSBarry Smith {
632*f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
633*f1af5d2fSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
634840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
635840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
636840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
637840b8ebdSBarry Smith   Vec           w1,w2,w3;
638840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
639*f1af5d2fSBarry Smith   PetscTruth    flg;
640840b8ebdSBarry Smith 
6413a40ed3dSBarry Smith   PetscFunctionBegin;
642840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
643840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
644840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
645840b8ebdSBarry Smith 
646840b8ebdSBarry Smith   if (!coloring->w1) {
647840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
648840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
649840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
650840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
651840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
652840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
653840b8ebdSBarry Smith   }
654840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
655840b8ebdSBarry Smith 
656*f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
657*f1af5d2fSBarry Smith   if (flg) {
658840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
659840b8ebdSBarry Smith   } else {
660840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
661840b8ebdSBarry Smith   }
662840b8ebdSBarry Smith 
663840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
664840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
665840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
666840b8ebdSBarry Smith 
667549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
668840b8ebdSBarry Smith   /*
669840b8ebdSBarry Smith       Loop over each color
670840b8ebdSBarry Smith   */
671840b8ebdSBarry Smith 
6723b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
673840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
674840b8ebdSBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
675840b8ebdSBarry Smith     /*
676840b8ebdSBarry Smith        Loop over each column associated with color adding the
677840b8ebdSBarry Smith        perturbation to the vector w3.
678840b8ebdSBarry Smith     */
679840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
680840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
681840b8ebdSBarry Smith       dx  = xx[col-start];
682840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
683aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
684840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
685840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
686840b8ebdSBarry Smith #else
687e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
688e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
689840b8ebdSBarry Smith #endif
690840b8ebdSBarry Smith       dx          *= epsilon;
691840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
6923b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
693840b8ebdSBarry Smith     }
694840b8ebdSBarry Smith     /*
695840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
696840b8ebdSBarry Smith     */
697840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
698840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
699840b8ebdSBarry Smith     /* Communicate scale to all processors */
7006831982aSBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
701840b8ebdSBarry Smith     /*
702840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
703840b8ebdSBarry Smith     */
7043b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
705840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
706840b8ebdSBarry Smith       row    = coloring->rows[k][l];
707840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
708840b8ebdSBarry Smith       y[row] *= scale[col];
709840b8ebdSBarry Smith       srow   = row + start;
710840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
711840b8ebdSBarry Smith     }
7123b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
713840b8ebdSBarry Smith   }
7143b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
715840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
716840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7173a40ed3dSBarry Smith   PetscFunctionReturn(0);
718840b8ebdSBarry Smith }
7193b28642cSBarry Smith 
7203b28642cSBarry Smith 
7213b28642cSBarry Smith 
722