xref: /petsc/src/mat/matfd/fdmatrix.c (revision 3f1db9ec2fd39765c6c3a00831044586630c4cca)
1840b8ebdSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*3f1db9ecSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.38 1998/12/03 03:59:42 bsmith Exp bsmith $";
4bbf0e169SBarry Smith #endif
5bbf0e169SBarry Smith 
6bbf0e169SBarry Smith /*
7639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
8639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
9bbf0e169SBarry Smith */
10bbf0e169SBarry Smith 
11bbf0e169SBarry Smith #include "petsc.h"
12bbf0e169SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
13bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
14bbf0e169SBarry Smith 
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 
83b4fc646aSLois Curfman McInnes    Notes:
84b4fc646aSLois Curfman McInnes    The available visualization contexts include
85ef5ee4d1SLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
86ef5ee4d1SLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
87ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
88ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
89ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
90ef5ee4d1SLois Curfman McInnes -     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
91bbf0e169SBarry Smith 
92639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
93005c665bSBarry Smith 
94b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
95bbf0e169SBarry Smith @*/
96b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
97bbf0e169SBarry Smith {
98005c665bSBarry Smith   ViewerType vtype;
99639f9d9dSBarry Smith   int        i,j,format,ierr;
100b4fc646aSLois Curfman McInnes   FILE       *fd;
101bbf0e169SBarry Smith 
1023a40ed3dSBarry Smith   PetscFunctionBegin;
103b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
104b4fc646aSLois Curfman McInnes   if (viewer) {PetscValidHeader(viewer);}
105b4fc646aSLois Curfman McInnes   else {viewer = VIEWER_STDOUT_SELF;}
106bbf0e169SBarry Smith 
107005c665bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
108*3f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,DRAW_VIEWER)) {
109b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr);
1103a40ed3dSBarry Smith     PetscFunctionReturn(0);
111*3f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
112b4fc646aSLois Curfman McInnes     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
113b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n");
114b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Error tolerance=%g\n",c->error_rel);
115b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Umin=%g\n",c->umin);
116b4fc646aSLois Curfman McInnes     PetscFPrintf(c->comm,fd,"  Number of colors=%d\n",c->ncolors);
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++ ) {
121b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"  Information for color %d\n",i);
122b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of columns %d\n",c->ncolumns[i]);
123b4fc646aSLois Curfman McInnes         for ( j=0; j<c->ncolumns[i]; j++ ) {
124b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d\n",c->columns[i][j]);
125639f9d9dSBarry Smith         }
126b4fc646aSLois Curfman McInnes         PetscFPrintf(c->comm,fd,"    Number of rows %d\n",c->nrows[i]);
127b4fc646aSLois Curfman McInnes         for ( j=0; j<c->nrows[i]; j++ ) {
128b4fc646aSLois Curfman McInnes           PetscFPrintf(c->comm,fd,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);
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
148ef5ee4d1SLois Curfman McInnes        J(u)_{:,i} = [J(u+h*dx_{i}) - J(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 
159b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
160b4fc646aSLois Curfman McInnes 
161b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
162639f9d9dSBarry Smith @*/
163639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin)
164639f9d9dSBarry Smith {
1653a40ed3dSBarry Smith   PetscFunctionBegin;
166639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
167639f9d9dSBarry Smith 
168639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
169639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1703a40ed3dSBarry Smith   PetscFunctionReturn(0);
171639f9d9dSBarry Smith }
172639f9d9dSBarry Smith 
1735615d1e5SSatish Balay #undef __FUNC__
174005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency"
175005c665bSBarry Smith /*@
176e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
177e0907662SLois Curfman McInnes    matrices.
178005c665bSBarry Smith 
179fee21e36SBarry Smith    Collective on MatFDColoring
180fee21e36SBarry Smith 
181ef5ee4d1SLois Curfman McInnes    Input Parameters:
182ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
183ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
184ef5ee4d1SLois Curfman McInnes 
185e0907662SLois Curfman McInnes    Notes:
186e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
187e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
188e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
189e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
190e0907662SLois Curfman McInnes 
191e0907662SLois Curfman McInnes    Options Database Keys:
192ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
193005c665bSBarry Smith 
194b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
195ef5ee4d1SLois Curfman McInnes 
196ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency()
197005c665bSBarry Smith @*/
198005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
199005c665bSBarry Smith {
2003a40ed3dSBarry Smith   PetscFunctionBegin;
201005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
202005c665bSBarry Smith 
203005c665bSBarry Smith   matfd->freq = freq;
2043a40ed3dSBarry Smith   PetscFunctionReturn(0);
205005c665bSBarry Smith }
206005c665bSBarry Smith 
207005c665bSBarry Smith #undef __FUNC__
208ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency"
209ff0cfa39SBarry Smith /*@
210ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
211ff0cfa39SBarry Smith    matrices.
212ff0cfa39SBarry Smith 
213ef5ee4d1SLois Curfman McInnes    Not Collective
214ef5ee4d1SLois Curfman McInnes 
215ff0cfa39SBarry Smith    Input Parameters:
216ff0cfa39SBarry Smith .  coloring - the coloring context
217ff0cfa39SBarry Smith 
218ff0cfa39SBarry Smith    Output Parameters:
219ff0cfa39SBarry Smith .  freq - frequency (default is 1)
220ff0cfa39SBarry Smith 
221ff0cfa39SBarry Smith    Notes:
222ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
223ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
224ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
225ff0cfa39SBarry Smith    <freq> nonlinear iterations.
226ff0cfa39SBarry Smith 
227ff0cfa39SBarry Smith    Options Database Keys:
228ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
229ff0cfa39SBarry Smith 
230ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
231ef5ee4d1SLois Curfman McInnes 
232ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
233ff0cfa39SBarry Smith @*/
234ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
235ff0cfa39SBarry Smith {
2363a40ed3dSBarry Smith   PetscFunctionBegin;
237ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
238ff0cfa39SBarry Smith 
239ff0cfa39SBarry Smith   *freq = matfd->freq;
2403a40ed3dSBarry Smith   PetscFunctionReturn(0);
241ff0cfa39SBarry Smith }
242ff0cfa39SBarry Smith 
243ff0cfa39SBarry Smith #undef __FUNC__
244005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction"
245d64ed03dSBarry Smith /*@C
246005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
247005c665bSBarry Smith 
248fee21e36SBarry Smith    Collective on MatFDColoring
249fee21e36SBarry Smith 
250ef5ee4d1SLois Curfman McInnes    Input Parameters:
251ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
252ef5ee4d1SLois Curfman McInnes .  f - the function
253ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
254ef5ee4d1SLois Curfman McInnes 
255b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
256005c665bSBarry Smith @*/
257840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
258005c665bSBarry Smith {
2593a40ed3dSBarry Smith   PetscFunctionBegin;
260005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
261005c665bSBarry Smith 
262005c665bSBarry Smith   matfd->f    = f;
263005c665bSBarry Smith   matfd->fctx = fctx;
264005c665bSBarry Smith 
2653a40ed3dSBarry Smith   PetscFunctionReturn(0);
266005c665bSBarry Smith }
267005c665bSBarry Smith 
268005c665bSBarry Smith #undef __FUNC__
269d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions"
270639f9d9dSBarry Smith /*@
271b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
272639f9d9dSBarry Smith    the options database.
273639f9d9dSBarry Smith 
274fee21e36SBarry Smith    Collective on MatFDColoring
275fee21e36SBarry Smith 
276ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
277ef5ee4d1SLois Curfman McInnes .vb
278ef5ee4d1SLois Curfman McInnes        J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where
279f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
280f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
281ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
282ef5ee4d1SLois Curfman McInnes .ve
283ef5ee4d1SLois Curfman McInnes 
284ef5ee4d1SLois Curfman McInnes    Input Parameter:
285ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
286ef5ee4d1SLois Curfman McInnes 
287b4fc646aSLois Curfman McInnes    Options Database Keys:
288ef5ee4d1SLois Curfman McInnes +  -mat_fd_coloring_error <err> - Sets <err> (square root
289ef5ee4d1SLois Curfman McInnes            of relative error in the function)
290f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
291ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
292ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
293ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
294ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
295639f9d9dSBarry Smith 
296b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
297639f9d9dSBarry Smith @*/
298639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
299639f9d9dSBarry Smith {
300005c665bSBarry Smith   int    ierr,flag,freq = 1;
301639f9d9dSBarry Smith   double error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
3023a40ed3dSBarry Smith 
3033a40ed3dSBarry Smith   PetscFunctionBegin;
304639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
305639f9d9dSBarry Smith 
306639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr);
307639f9d9dSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr);
308639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr);
309005c665bSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr);
310005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
311005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr);
312639f9d9dSBarry Smith   if (flag) {
313639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr);
314639f9d9dSBarry Smith   }
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
316639f9d9dSBarry Smith }
317639f9d9dSBarry Smith 
3185615d1e5SSatish Balay #undef __FUNC__
319d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp"
320639f9d9dSBarry Smith /*@
321639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
322639f9d9dSBarry Smith     using coloring.
323639f9d9dSBarry Smith 
324ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
325ef5ee4d1SLois Curfman McInnes 
326639f9d9dSBarry Smith     Input Parameter:
327639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
328639f9d9dSBarry Smith 
329639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
330639f9d9dSBarry Smith @*/
331639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
332639f9d9dSBarry Smith {
3333a40ed3dSBarry Smith   PetscFunctionBegin;
334639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
335639f9d9dSBarry Smith 
33676be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);
33776be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);
33876be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);
33976be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");
34076be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");
34176be9ce4SBarry Smith   (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
343005c665bSBarry Smith }
344005c665bSBarry Smith 
345005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
346005c665bSBarry Smith {
347005c665bSBarry Smith   int ierr,flg;
348005c665bSBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
350005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr);
351005c665bSBarry Smith   if (flg) {
352f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
353005c665bSBarry Smith   }
354ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr);
355ae09f205SBarry Smith   if (flg) {
356f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
357f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr);
358f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
359ae09f205SBarry Smith   }
360005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr);
361005c665bSBarry Smith   if (flg) {
362005c665bSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
363005c665bSBarry Smith     ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr);
364005c665bSBarry Smith   }
3653a40ed3dSBarry Smith   PetscFunctionReturn(0);
366bbf0e169SBarry Smith }
367bbf0e169SBarry Smith 
3685615d1e5SSatish Balay #undef __FUNC__
3695615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate"
370bbf0e169SBarry Smith /*@C
371639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
372639f9d9dSBarry Smith    computation of Jacobians.
373bbf0e169SBarry Smith 
374ef5ee4d1SLois Curfman McInnes    Collective on Mat
375ef5ee4d1SLois Curfman McInnes 
376639f9d9dSBarry Smith    Input Parameters:
377ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
378ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
379bbf0e169SBarry Smith 
380bbf0e169SBarry Smith     Output Parameter:
381639f9d9dSBarry Smith .   color - the new coloring context
382bbf0e169SBarry Smith 
383b4fc646aSLois Curfman McInnes     Options Database Keys:
384ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
385ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
386ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
387639f9d9dSBarry Smith 
388639f9d9dSBarry Smith .seealso: MatFDColoringDestroy()
389bbf0e169SBarry Smith @*/
390639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
391bbf0e169SBarry Smith {
392639f9d9dSBarry Smith   MatFDColoring c;
393639f9d9dSBarry Smith   MPI_Comm      comm;
394639f9d9dSBarry Smith   int           ierr,M,N;
395639f9d9dSBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
397639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr);
398e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
399639f9d9dSBarry Smith 
400639f9d9dSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
401*3f1db9ecSBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,
402*3f1db9ecSBarry Smith                     MatFDColoringDestroy,MatFDColoringView);
403639f9d9dSBarry Smith   PLogObjectCreate(c);
404639f9d9dSBarry Smith 
405f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
406f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr);
407639f9d9dSBarry Smith   } else {
408e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
409639f9d9dSBarry Smith   }
410639f9d9dSBarry Smith 
411639f9d9dSBarry Smith   c->error_rel = 1.e-8;
412ae09f205SBarry Smith   c->umin      = 1.e-6;
413005c665bSBarry Smith   c->freq      = 1;
414005c665bSBarry Smith 
415005c665bSBarry Smith   ierr = MatFDColoringView_Private(c); CHKERRQ(ierr);
416639f9d9dSBarry Smith 
417639f9d9dSBarry Smith   *color = c;
418639f9d9dSBarry Smith 
4193a40ed3dSBarry Smith   PetscFunctionReturn(0);
420bbf0e169SBarry Smith }
421bbf0e169SBarry Smith 
4225615d1e5SSatish Balay #undef __FUNC__
423d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy"
424bbf0e169SBarry Smith /*@C
425639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
426639f9d9dSBarry Smith     via MatFDColoringCreate().
427bbf0e169SBarry Smith 
428ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
429ef5ee4d1SLois Curfman McInnes 
430b4fc646aSLois Curfman McInnes     Input Parameter:
431639f9d9dSBarry Smith .   c - coloring context
432bbf0e169SBarry Smith 
433639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
434bbf0e169SBarry Smith @*/
435639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
436bbf0e169SBarry Smith {
437263760aaSBarry Smith   int i,ierr;
438bbf0e169SBarry Smith 
4393a40ed3dSBarry Smith   PetscFunctionBegin;
4403a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
441d4bb536fSBarry Smith 
442639f9d9dSBarry Smith 
443639f9d9dSBarry Smith   for ( i=0; i<c->ncolors; i++ ) {
444639f9d9dSBarry Smith     if (c->columns[i])       PetscFree(c->columns[i]);
445639f9d9dSBarry Smith     if (c->rows[i])          PetscFree(c->rows[i]);
446639f9d9dSBarry Smith     if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]);
447bbf0e169SBarry Smith   }
448639f9d9dSBarry Smith   PetscFree(c->ncolumns);
449639f9d9dSBarry Smith   PetscFree(c->columns);
450639f9d9dSBarry Smith   PetscFree(c->nrows);
451639f9d9dSBarry Smith   PetscFree(c->rows);
452639f9d9dSBarry Smith   PetscFree(c->columnsforrow);
453639f9d9dSBarry Smith   PetscFree(c->scale);
454005c665bSBarry Smith   if (c->w1) {
455005c665bSBarry Smith     ierr = VecDestroy(c->w1); CHKERRQ(ierr);
456005c665bSBarry Smith     ierr = VecDestroy(c->w2); CHKERRQ(ierr);
457005c665bSBarry Smith     ierr = VecDestroy(c->w3); CHKERRQ(ierr);
458005c665bSBarry Smith   }
459639f9d9dSBarry Smith   PLogObjectDestroy(c);
460639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4613a40ed3dSBarry Smith   PetscFunctionReturn(0);
462bbf0e169SBarry Smith }
46343a90d84SBarry Smith 
464005c665bSBarry Smith #include "snes.h"
465005c665bSBarry Smith 
4665615d1e5SSatish Balay #undef __FUNC__
4675615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply"
46843a90d84SBarry Smith /*@
469e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
470e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
47143a90d84SBarry Smith 
472fee21e36SBarry Smith     Collective on MatFDColoring
473fee21e36SBarry Smith 
474ef5ee4d1SLois Curfman McInnes     Input Parameters:
475ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
476ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
477ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
478ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
479ef5ee4d1SLois Curfman McInnes 
4808bba8e72SBarry Smith    Options Database Keys:
481ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
4828bba8e72SBarry Smith 
48343a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
48443a90d84SBarry Smith 
48543a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
48643a90d84SBarry Smith @*/
487005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
48843a90d84SBarry Smith {
489e0907662SLois Curfman McInnes   int           k,fg,ierr,N,start,end,l,row,col,srow;
49043a90d84SBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
49143a90d84SBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
49243a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
493005c665bSBarry Smith   Vec           w1,w2,w3;
494840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f;
495005c665bSBarry Smith   void          *fctx = coloring->fctx;
496005c665bSBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
498e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
499e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
500e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
501e0907662SLois Curfman McInnes 
502005c665bSBarry Smith 
503005c665bSBarry Smith   if (!coloring->w1) {
504005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
505005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
506005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
507005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
508005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
509005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
510005c665bSBarry Smith   }
511005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
51243a90d84SBarry Smith 
513e0907662SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
514e0907662SLois Curfman McInnes   if (fg) {
515e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
516e0907662SLois Curfman McInnes   } else {
51743a90d84SBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
518e0907662SLois Curfman McInnes   }
51943a90d84SBarry Smith 
52043a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
52143a90d84SBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
52243a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr);
52343a90d84SBarry Smith 
52443a90d84SBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
52543a90d84SBarry Smith   /*
52643a90d84SBarry Smith       Loop over each color
52743a90d84SBarry Smith   */
52843a90d84SBarry Smith 
5293b28642cSBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
53043a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
53143a90d84SBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
53243a90d84SBarry Smith     /*
53343a90d84SBarry Smith        Loop over each column associated with color adding the
53443a90d84SBarry Smith        perturbation to the vector w3.
53543a90d84SBarry Smith     */
53643a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
53743a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
53843a90d84SBarry Smith       dx  = xx[col-start];
539ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
5403a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
541ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
542ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
54343a90d84SBarry Smith #else
544e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
545e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
54643a90d84SBarry Smith #endif
54743a90d84SBarry Smith       dx          *= epsilon;
54843a90d84SBarry Smith       wscale[col] = 1.0/dx;
5493b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
55043a90d84SBarry Smith     }
5513b28642cSBarry Smith 
55243a90d84SBarry Smith     /*
553e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
55443a90d84SBarry Smith     */
55543a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr);
55643a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
55743a90d84SBarry Smith     /* Communicate scale to all processors */
5583a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
559ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
56043a90d84SBarry Smith #else
561ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
56243a90d84SBarry Smith #endif
56343a90d84SBarry Smith     /*
564e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
56543a90d84SBarry Smith     */
5663b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
56743a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
56843a90d84SBarry Smith       row    = coloring->rows[k][l];
56943a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
57043a90d84SBarry Smith       y[row] *= scale[col];
57143a90d84SBarry Smith       srow   = row + start;
57243a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
57343a90d84SBarry Smith     }
5743b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y); CHKERRQ(ierr);
57543a90d84SBarry Smith   }
5763b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
57743a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
57843a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5793a40ed3dSBarry Smith   PetscFunctionReturn(0);
58043a90d84SBarry Smith }
581840b8ebdSBarry Smith 
582840b8ebdSBarry Smith #include "ts.h"
583840b8ebdSBarry Smith 
584840b8ebdSBarry Smith #undef __FUNC__
585840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS"
586840b8ebdSBarry Smith /*@
587840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
588840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
589840b8ebdSBarry Smith 
590fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
591fee21e36SBarry Smith 
592ef5ee4d1SLois Curfman McInnes     Input Parameters:
5933b28642cSBarry Smith +   mat - location to store Jacobian
594ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
595ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
596ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
597ef5ee4d1SLois Curfman McInnes 
598840b8ebdSBarry Smith    Options Database Keys:
599ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
600840b8ebdSBarry Smith 
601840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
602840b8ebdSBarry Smith 
603840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
604840b8ebdSBarry Smith @*/
605840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx)
606840b8ebdSBarry Smith {
607840b8ebdSBarry Smith   int           k,fg,ierr,N,start,end,l,row,col,srow;
608840b8ebdSBarry Smith   Scalar        dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
609840b8ebdSBarry Smith   double        epsilon = coloring->error_rel, umin = coloring->umin;
610840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
611840b8ebdSBarry Smith   Vec           w1,w2,w3;
612840b8ebdSBarry Smith   int           (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f;
613840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
614840b8ebdSBarry Smith 
6153a40ed3dSBarry Smith   PetscFunctionBegin;
616840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
617840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
618840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
619840b8ebdSBarry Smith 
620840b8ebdSBarry Smith   if (!coloring->w1) {
621840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr);
622840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
623840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr);
624840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
625840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr);
626840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
627840b8ebdSBarry Smith   }
628840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
629840b8ebdSBarry Smith 
630840b8ebdSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr);
631840b8ebdSBarry Smith   if (fg) {
632840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
633840b8ebdSBarry Smith   } else {
634840b8ebdSBarry Smith     ierr = MatZeroEntries(J); CHKERRQ(ierr);
635840b8ebdSBarry Smith   }
636840b8ebdSBarry Smith 
637840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr);
638840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N); CHKERRQ(ierr);
639840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr);
640840b8ebdSBarry Smith 
641840b8ebdSBarry Smith   PetscMemzero(wscale,N*sizeof(Scalar));
642840b8ebdSBarry Smith   /*
643840b8ebdSBarry Smith       Loop over each color
644840b8ebdSBarry Smith   */
645840b8ebdSBarry Smith 
6463b28642cSBarry Smith   ierr = VecGetArray(x1,&xx); CHKERRQ(ierr);
647840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
648840b8ebdSBarry Smith     ierr = VecCopy(x1,w3); CHKERRQ(ierr);
649840b8ebdSBarry Smith     /*
650840b8ebdSBarry Smith        Loop over each column associated with color adding the
651840b8ebdSBarry Smith        perturbation to the vector w3.
652840b8ebdSBarry Smith     */
653840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
654840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
655840b8ebdSBarry Smith       dx  = xx[col-start];
656840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
6573a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
658840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
659840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
660840b8ebdSBarry Smith #else
661e20fef11SSatish Balay       if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0)     dx = umin;
662e20fef11SSatish Balay       else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
663840b8ebdSBarry Smith #endif
664840b8ebdSBarry Smith       dx          *= epsilon;
665840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
6663b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES); CHKERRQ(ierr);
667840b8ebdSBarry Smith     }
668840b8ebdSBarry Smith     /*
669840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
670840b8ebdSBarry Smith     */
671840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr);
672840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr);
673840b8ebdSBarry Smith     /* Communicate scale to all processors */
6743a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
675ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
676840b8ebdSBarry Smith #else
677ca161407SBarry Smith     ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
678840b8ebdSBarry Smith #endif
679840b8ebdSBarry Smith     /*
680840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
681840b8ebdSBarry Smith     */
6823b28642cSBarry Smith     ierr = VecGetArray(w2,&y); CHKERRQ(ierr);
683840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
684840b8ebdSBarry Smith       row    = coloring->rows[k][l];
685840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
686840b8ebdSBarry Smith       y[row] *= scale[col];
687840b8ebdSBarry Smith       srow   = row + start;
688840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
689840b8ebdSBarry Smith     }
6903b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y); CHKERRQ(ierr);
691840b8ebdSBarry Smith   }
6923b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx); CHKERRQ(ierr);
693840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
694840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
6953a40ed3dSBarry Smith   PetscFunctionReturn(0);
696840b8ebdSBarry Smith }
6973b28642cSBarry Smith 
6983b28642cSBarry Smith 
6993b28642cSBarry Smith 
700