xref: /petsc/src/mat/matfd/fdmatrix.c (revision 54d96fa383ba80e0a45dd9f0fcec1909a0522f96)
1*54d96fa3SBarry Smith /*$Id: fdmatrix.c,v 1.62 2000/05/13 04:18:00 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"
9e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
10bbf0e169SBarry Smith #include "src/vec/vecimpl.h"
11bbf0e169SBarry Smith 
125615d1e5SSatish Balay #undef __FUNC__
139194eea9SBarry Smith #define __FUNC__ /*<a name="MatFDColoringView_Draw_Zoom"></a>*/"MatFDColoringView_Draw_Zoom"
149194eea9SBarry Smith static int MatFDColoringView_Draw_Zoom(Draw draw,void *Aa)
159194eea9SBarry Smith {
169194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
179194eea9SBarry Smith   int           ierr,i,j;
189194eea9SBarry Smith   PetscReal     x,y;
199194eea9SBarry Smith 
209194eea9SBarry Smith   PetscFunctionBegin;
219194eea9SBarry Smith 
229194eea9SBarry Smith   /* loop over colors  */
239194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
249194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
259194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
269194eea9SBarry Smith       x = fd->columnsforrow[i][j];
279194eea9SBarry Smith       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
289194eea9SBarry Smith     }
299194eea9SBarry Smith   }
309194eea9SBarry Smith   PetscFunctionReturn(0);
319194eea9SBarry Smith }
329194eea9SBarry Smith 
339194eea9SBarry Smith #undef __FUNC__
349194eea9SBarry Smith #define __FUNC__ /*<a name="MatFDColoringView_Draw"></a>*/"MatFDColoringView_Draw"
35005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
36005c665bSBarry Smith {
379194eea9SBarry Smith   int         ierr;
38005c665bSBarry Smith   PetscTruth  isnull;
39005c665bSBarry Smith   Draw        draw;
40*54d96fa3SBarry Smith   PetscReal   xr,yr,xl,yl,h,w;
41005c665bSBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
4377ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
443a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
459194eea9SBarry Smith 
469194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
47005c665bSBarry Smith 
48005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
49005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
50005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
519194eea9SBarry Smith   ierr = DrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
529194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54005c665bSBarry Smith }
55005c665bSBarry Smith 
56005c665bSBarry Smith #undef __FUNC__
57b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView"
58bbf0e169SBarry Smith /*@C
59639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
60bbf0e169SBarry Smith 
614c49b128SBarry Smith    Collective on MatFDColoring
62fee21e36SBarry Smith 
63ef5ee4d1SLois Curfman McInnes    Input  Parameters:
64ef5ee4d1SLois Curfman McInnes +  c - the coloring context
65ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
66ef5ee4d1SLois Curfman McInnes 
6715091d37SBarry Smith    Level: intermediate
6815091d37SBarry Smith 
69b4fc646aSLois Curfman McInnes    Notes:
70b4fc646aSLois Curfman McInnes    The available visualization contexts include
71ef5ee4d1SLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
72ef5ee4d1SLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
73ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
74ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
75ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
76c655490fSBarry Smith -     VIEWER_DRAW_WORLD - graphical display of nonzero structure
77bbf0e169SBarry Smith 
789194eea9SBarry Smith    Notes:
799194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
809194eea9SBarry Smith    involves moe than 33 then some seemingly identical colors are displayed making it look
819194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
829194eea9SBarry Smith 
83639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
84005c665bSBarry Smith 
85b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
86bbf0e169SBarry Smith @*/
87b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
88bbf0e169SBarry Smith {
89639f9d9dSBarry Smith   int        i,j,format,ierr;
906831982aSBarry Smith   PetscTruth isdraw,isascii;
91bbf0e169SBarry Smith 
923a40ed3dSBarry Smith   PetscFunctionBegin;
93b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
943eda8832SBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_(c->comm);
950f5bd95cSBarry Smith   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
966831982aSBarry Smith   PetscCheckSameComm(c,viewer);
97bbf0e169SBarry Smith 
986831982aSBarry Smith   ierr  = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
996831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
1000f5bd95cSBarry Smith   if (isdraw) {
101b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1020f5bd95cSBarry Smith   } else if (isascii) {
103d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
104d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
105d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
106d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
107ae09f205SBarry Smith 
108ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
109ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
110b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
111d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
112d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
113b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
114d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
115639f9d9dSBarry Smith         }
116d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
117b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
118d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
119b4fc646aSLois Curfman McInnes         }
120bbf0e169SBarry Smith       }
121bbf0e169SBarry Smith     }
1226831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
1235cd90555SBarry Smith   } else {
1240f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
125bbf0e169SBarry Smith   }
1263a40ed3dSBarry Smith   PetscFunctionReturn(0);
127639f9d9dSBarry Smith }
128639f9d9dSBarry Smith 
1295615d1e5SSatish Balay #undef __FUNC__
130b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetParameters"
131639f9d9dSBarry Smith /*@
132b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
133b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
134639f9d9dSBarry Smith 
135ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
136ef5ee4d1SLois Curfman McInnes 
137ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
138ef5ee4d1SLois Curfman McInnes .vb
13965f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
140f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
141f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
142ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
143ef5ee4d1SLois Curfman McInnes .ve
144639f9d9dSBarry Smith 
145639f9d9dSBarry Smith    Input Parameters:
146ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
147639f9d9dSBarry Smith .  error_rel - relative error
148f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
149fee21e36SBarry Smith 
15015091d37SBarry Smith    Level: advanced
15115091d37SBarry Smith 
152b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
153b4fc646aSLois Curfman McInnes 
154b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
155639f9d9dSBarry Smith @*/
156329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
157639f9d9dSBarry Smith {
1583a40ed3dSBarry Smith   PetscFunctionBegin;
159639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
160639f9d9dSBarry Smith 
161639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
162639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1633a40ed3dSBarry Smith   PetscFunctionReturn(0);
164639f9d9dSBarry Smith }
165639f9d9dSBarry Smith 
1665615d1e5SSatish Balay #undef __FUNC__
167b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFrequency"
168005c665bSBarry Smith /*@
169e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
170e0907662SLois Curfman McInnes    matrices.
171005c665bSBarry Smith 
172fee21e36SBarry Smith    Collective on MatFDColoring
173fee21e36SBarry Smith 
174ef5ee4d1SLois Curfman McInnes    Input Parameters:
175ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
176ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
177ef5ee4d1SLois Curfman McInnes 
17815091d37SBarry Smith    Options Database Keys:
17915091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
18015091d37SBarry Smith 
18115091d37SBarry Smith    Level: advanced
18215091d37SBarry Smith 
183e0907662SLois Curfman McInnes    Notes:
184e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
185e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
186e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
187e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
188e0907662SLois Curfman McInnes 
189b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
190ef5ee4d1SLois Curfman McInnes 
191ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency()
192005c665bSBarry Smith @*/
193005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
194005c665bSBarry Smith {
1953a40ed3dSBarry Smith   PetscFunctionBegin;
196005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
197005c665bSBarry Smith 
198005c665bSBarry Smith   matfd->freq = freq;
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
200005c665bSBarry Smith }
201005c665bSBarry Smith 
202005c665bSBarry Smith #undef __FUNC__
203b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringGetFrequency"
204ff0cfa39SBarry Smith /*@
205ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
206ff0cfa39SBarry Smith    matrices.
207ff0cfa39SBarry Smith 
208ef5ee4d1SLois Curfman McInnes    Not Collective
209ef5ee4d1SLois Curfman McInnes 
210ff0cfa39SBarry Smith    Input Parameters:
211ff0cfa39SBarry Smith .  coloring - the coloring context
212ff0cfa39SBarry Smith 
213ff0cfa39SBarry Smith    Output Parameters:
214ff0cfa39SBarry Smith .  freq - frequency (default is 1)
215ff0cfa39SBarry Smith 
21615091d37SBarry Smith    Options Database Keys:
21715091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
21815091d37SBarry Smith 
21915091d37SBarry Smith    Level: advanced
22015091d37SBarry 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 .keywords: Mat, finite differences, coloring, get, frequency
228ef5ee4d1SLois Curfman McInnes 
229ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
230ff0cfa39SBarry Smith @*/
231ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
232ff0cfa39SBarry Smith {
2333a40ed3dSBarry Smith   PetscFunctionBegin;
234ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
235ff0cfa39SBarry Smith 
236ff0cfa39SBarry Smith   *freq = matfd->freq;
2373a40ed3dSBarry Smith   PetscFunctionReturn(0);
238ff0cfa39SBarry Smith }
239ff0cfa39SBarry Smith 
240ff0cfa39SBarry Smith #undef __FUNC__
241b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFunction"
242d64ed03dSBarry Smith /*@C
243005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
244005c665bSBarry Smith 
245fee21e36SBarry Smith    Collective on MatFDColoring
246fee21e36SBarry Smith 
247ef5ee4d1SLois Curfman McInnes    Input Parameters:
248ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
249ef5ee4d1SLois Curfman McInnes .  f - the function
250ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
251ef5ee4d1SLois Curfman McInnes 
25215091d37SBarry Smith    Level: intermediate
25315091d37SBarry Smith 
254f881d145SBarry Smith    Notes:
255f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
256f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
257f881d145SBarry Smith   with the TS solvers.
258f881d145SBarry Smith 
259b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
260005c665bSBarry Smith @*/
261840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
262005c665bSBarry Smith {
2633a40ed3dSBarry Smith   PetscFunctionBegin;
264005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
265005c665bSBarry Smith 
266005c665bSBarry Smith   matfd->f    = f;
267005c665bSBarry Smith   matfd->fctx = fctx;
268005c665bSBarry Smith 
2693a40ed3dSBarry Smith   PetscFunctionReturn(0);
270005c665bSBarry Smith }
271005c665bSBarry Smith 
272005c665bSBarry Smith #undef __FUNC__
273b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFromOptions"
274639f9d9dSBarry Smith /*@
275b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
276639f9d9dSBarry Smith    the options database.
277639f9d9dSBarry Smith 
278fee21e36SBarry Smith    Collective on MatFDColoring
279fee21e36SBarry Smith 
28065f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
281ef5ee4d1SLois Curfman McInnes .vb
28265f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
283f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
284f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
285ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
286ef5ee4d1SLois Curfman McInnes .ve
287ef5ee4d1SLois Curfman McInnes 
288ef5ee4d1SLois Curfman McInnes    Input Parameter:
289ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
290ef5ee4d1SLois Curfman McInnes 
291b4fc646aSLois Curfman McInnes    Options Database Keys:
292ef5ee4d1SLois Curfman McInnes +  -mat_fd_coloring_error <err> - Sets <err> (square root
293ef5ee4d1SLois Curfman McInnes            of relative error in the function)
294f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
295ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
296ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
297ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
298ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
299639f9d9dSBarry Smith 
30015091d37SBarry Smith     Level: intermediate
30115091d37SBarry Smith 
302b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
303639f9d9dSBarry Smith @*/
304639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
305639f9d9dSBarry Smith {
306f1af5d2fSBarry Smith   int        ierr,freq = 1;
307329f5518SBarry Smith   PetscReal  error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
308f1af5d2fSBarry Smith   PetscTruth flag;
3093a40ed3dSBarry Smith 
3103a40ed3dSBarry Smith   PetscFunctionBegin;
311639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
312639f9d9dSBarry Smith 
313f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,PETSC_NULL);CHKERRQ(ierr);
314f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,PETSC_NULL);CHKERRQ(ierr);
315639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin);CHKERRQ(ierr);
316f1af5d2fSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,PETSC_NULL);CHKERRQ(ierr);
317005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
318005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag);CHKERRQ(ierr);
319639f9d9dSBarry Smith   if (flag) {
320639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd);CHKERRQ(ierr);
321639f9d9dSBarry Smith   }
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
323639f9d9dSBarry Smith }
324639f9d9dSBarry Smith 
3255615d1e5SSatish Balay #undef __FUNC__
326b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringPrintHelp"
327639f9d9dSBarry Smith /*@
328639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
329639f9d9dSBarry Smith     using coloring.
330639f9d9dSBarry Smith 
331ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
332ef5ee4d1SLois Curfman McInnes 
333639f9d9dSBarry Smith     Input Parameter:
334639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
335639f9d9dSBarry Smith 
33615091d37SBarry Smith     Level: intermediate
33715091d37SBarry Smith 
338639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
339639f9d9dSBarry Smith @*/
340639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
341639f9d9dSBarry Smith {
342d132466eSBarry Smith   int ierr;
343d132466eSBarry Smith 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
345639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
346639f9d9dSBarry Smith 
347d132466eSBarry 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);
348d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);CHKERRQ(ierr);
349d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);CHKERRQ(ierr);
350d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");CHKERRQ(ierr);
351d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");CHKERRQ(ierr);
352d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");CHKERRQ(ierr);
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
354005c665bSBarry Smith }
355005c665bSBarry Smith 
356433994e6SBarry Smith #undef __FUNC__
357b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView_Private"
358005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
359005c665bSBarry Smith {
360f1af5d2fSBarry Smith   int        ierr;
361f1af5d2fSBarry Smith   PetscTruth flg;
362005c665bSBarry Smith 
3633a40ed3dSBarry Smith   PetscFunctionBegin;
364005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
365005c665bSBarry Smith   if (flg) {
366f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
367005c665bSBarry Smith   }
368ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
369ae09f205SBarry Smith   if (flg) {
370f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
371f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
372f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
373ae09f205SBarry Smith   }
374005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
375005c665bSBarry Smith   if (flg) {
376c655490fSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
377c655490fSBarry Smith     ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
378005c665bSBarry Smith   }
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
380bbf0e169SBarry Smith }
381bbf0e169SBarry Smith 
3825615d1e5SSatish Balay #undef __FUNC__
383b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringCreate"
384bbf0e169SBarry Smith /*@C
385639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
386639f9d9dSBarry Smith    computation of Jacobians.
387bbf0e169SBarry Smith 
388ef5ee4d1SLois Curfman McInnes    Collective on Mat
389ef5ee4d1SLois Curfman McInnes 
390639f9d9dSBarry Smith    Input Parameters:
391ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
392ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
393bbf0e169SBarry Smith 
394bbf0e169SBarry Smith     Output Parameter:
395639f9d9dSBarry Smith .   color - the new coloring context
396bbf0e169SBarry Smith 
397b4fc646aSLois Curfman McInnes     Options Database Keys:
398ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
399ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
400ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
401639f9d9dSBarry Smith 
40215091d37SBarry Smith     Level: intermediate
40315091d37SBarry Smith 
4046831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
4056831982aSBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
406bbf0e169SBarry Smith @*/
407639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
408bbf0e169SBarry Smith {
409639f9d9dSBarry Smith   MatFDColoring c;
410639f9d9dSBarry Smith   MPI_Comm      comm;
411639f9d9dSBarry Smith   int           ierr,M,N;
412639f9d9dSBarry Smith 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
414639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
415e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
416639f9d9dSBarry Smith 
417f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4189194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
419639f9d9dSBarry Smith   PLogObjectCreate(c);
420639f9d9dSBarry Smith 
421f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
422f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
423639f9d9dSBarry Smith   } else {
424e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
425639f9d9dSBarry Smith   }
426639f9d9dSBarry Smith 
427639f9d9dSBarry Smith   c->error_rel = 1.e-8;
428ae09f205SBarry Smith   c->umin      = 1.e-6;
429005c665bSBarry Smith   c->freq      = 1;
430005c665bSBarry Smith 
431005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
432639f9d9dSBarry Smith 
433639f9d9dSBarry Smith   *color = c;
434639f9d9dSBarry Smith 
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
436bbf0e169SBarry Smith }
437bbf0e169SBarry Smith 
4385615d1e5SSatish Balay #undef __FUNC__
439b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringDestroy"
440bbf0e169SBarry Smith /*@C
441639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
442639f9d9dSBarry Smith     via MatFDColoringCreate().
443bbf0e169SBarry Smith 
444ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
445ef5ee4d1SLois Curfman McInnes 
446b4fc646aSLois Curfman McInnes     Input Parameter:
447639f9d9dSBarry Smith .   c - coloring context
448bbf0e169SBarry Smith 
44915091d37SBarry Smith     Level: intermediate
45015091d37SBarry Smith 
451639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
452bbf0e169SBarry Smith @*/
453639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
454bbf0e169SBarry Smith {
455263760aaSBarry Smith   int i,ierr;
456bbf0e169SBarry Smith 
4573a40ed3dSBarry Smith   PetscFunctionBegin;
4583a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
459d4bb536fSBarry Smith 
460639f9d9dSBarry Smith 
461639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
462606d414cSSatish Balay     if (c->columns[i])       {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
463606d414cSSatish Balay     if (c->rows[i])          {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
464606d414cSSatish Balay     if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
465bbf0e169SBarry Smith   }
466606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
467606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
468606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
469606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
470606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
471606d414cSSatish Balay   ierr = PetscFree(c->scale);CHKERRQ(ierr);
472005c665bSBarry Smith   if (c->w1) {
473005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
474005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
475005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
476005c665bSBarry Smith   }
477639f9d9dSBarry Smith   PLogObjectDestroy(c);
478639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
480bbf0e169SBarry Smith }
48143a90d84SBarry Smith 
482005c665bSBarry Smith 
4835615d1e5SSatish Balay #undef __FUNC__
484b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApply"
48543a90d84SBarry Smith /*@
486e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
487e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48843a90d84SBarry Smith 
489fee21e36SBarry Smith     Collective on MatFDColoring
490fee21e36SBarry Smith 
491ef5ee4d1SLois Curfman McInnes     Input Parameters:
492ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
493ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
494ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
495ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
496ef5ee4d1SLois Curfman McInnes 
4978bba8e72SBarry Smith    Options Database Keys:
498ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
4998bba8e72SBarry Smith 
50015091d37SBarry Smith    Level: intermediate
50115091d37SBarry Smith 
50243a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
50343a90d84SBarry Smith 
50443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50543a90d84SBarry Smith @*/
506005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50743a90d84SBarry Smith {
508f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
5099194eea9SBarry Smith   Scalar        dx,mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale,*w3_array;
510329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
51143a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
512005c665bSBarry Smith   Vec           w1,w2,w3;
513840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void*)= (int (*)(void *,Vec,Vec,void *))coloring->f;
514005c665bSBarry Smith   void          *fctx = coloring->fctx;
515f1af5d2fSBarry Smith   PetscTruth    flg;
516005c665bSBarry Smith 
5173a40ed3dSBarry Smith   PetscFunctionBegin;
518e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
519e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
520e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
521e0907662SLois Curfman McInnes 
522005c665bSBarry Smith 
523005c665bSBarry Smith   if (!coloring->w1) {
524005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
525005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
526005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
527005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
528005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
529005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
530005c665bSBarry Smith   }
531005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
53243a90d84SBarry Smith 
5334c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
534f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
535f1af5d2fSBarry Smith   if (flg) {
536e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
537e0907662SLois Curfman McInnes   } else {
53843a90d84SBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
539e0907662SLois Curfman McInnes   }
54043a90d84SBarry Smith 
54143a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
54243a90d84SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
54343a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
54443a90d84SBarry Smith 
545549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
54643a90d84SBarry Smith   /*
54743a90d84SBarry Smith       Loop over each color
54843a90d84SBarry Smith   */
54943a90d84SBarry Smith 
5503b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
55143a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
5529194eea9SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
55343a90d84SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
55443a90d84SBarry Smith     /*
55543a90d84SBarry Smith        Loop over each column associated with color adding the
55643a90d84SBarry Smith        perturbation to the vector w3.
55743a90d84SBarry Smith     */
55843a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
55943a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
56043a90d84SBarry Smith       dx  = xx[col-start];
561ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
562aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
563ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
564ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
56543a90d84SBarry Smith #else
566329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
567329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
56843a90d84SBarry Smith #endif
56943a90d84SBarry Smith       dx                    *= epsilon;
57043a90d84SBarry Smith       wscale[col]            = 1.0/dx;
5719194eea9SBarry Smith       w3_array[col - start] += dx;
57243a90d84SBarry Smith     }
5739194eea9SBarry Smith     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
5743b28642cSBarry Smith 
57543a90d84SBarry Smith     /*
576e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
57743a90d84SBarry Smith     */
57843a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
57943a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
58043a90d84SBarry Smith     /* Communicate scale to all processors */
5816831982aSBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
58243a90d84SBarry Smith     /*
583e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
58443a90d84SBarry Smith     */
5853b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
58643a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
58743a90d84SBarry Smith       row    = coloring->rows[k][l];
58843a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
58943a90d84SBarry Smith       y[row] *= scale[col];
59043a90d84SBarry Smith       srow   = row + start;
59143a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
59243a90d84SBarry Smith     }
5933b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
59443a90d84SBarry Smith   }
5953b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
59643a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59743a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5983a40ed3dSBarry Smith   PetscFunctionReturn(0);
59943a90d84SBarry Smith }
600840b8ebdSBarry Smith 
601840b8ebdSBarry Smith #undef __FUNC__
602b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApplyTS"
603840b8ebdSBarry Smith /*@
604840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
605840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
606840b8ebdSBarry Smith 
607fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
608fee21e36SBarry Smith 
609ef5ee4d1SLois Curfman McInnes     Input Parameters:
6103b28642cSBarry Smith +   mat - location to store Jacobian
611ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
612ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
613ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
614ef5ee4d1SLois Curfman McInnes 
615840b8ebdSBarry Smith    Options Database Keys:
616ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
617840b8ebdSBarry Smith 
61815091d37SBarry Smith    Level: intermediate
61915091d37SBarry Smith 
620840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
621840b8ebdSBarry Smith 
622840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
623840b8ebdSBarry Smith @*/
624329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
625840b8ebdSBarry Smith {
626f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
627329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)= (int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
628840b8ebdSBarry Smith   Scalar        dx,mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
629329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
630840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
631840b8ebdSBarry Smith   Vec           w1,w2,w3;
632840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
633f1af5d2fSBarry Smith   PetscTruth    flg;
634840b8ebdSBarry Smith 
6353a40ed3dSBarry Smith   PetscFunctionBegin;
636840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
637840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
638840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
639840b8ebdSBarry Smith 
640840b8ebdSBarry Smith   if (!coloring->w1) {
641840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
642840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
643840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
644840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
645840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
646840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
647840b8ebdSBarry Smith   }
648840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
649840b8ebdSBarry Smith 
650f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
651f1af5d2fSBarry Smith   if (flg) {
652840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
653840b8ebdSBarry Smith   } else {
654840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
655840b8ebdSBarry Smith   }
656840b8ebdSBarry Smith 
657840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
658840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
659840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
660840b8ebdSBarry Smith 
661549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
662840b8ebdSBarry Smith   /*
663840b8ebdSBarry Smith       Loop over each color
664840b8ebdSBarry Smith   */
665840b8ebdSBarry Smith 
6663b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
667840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
668840b8ebdSBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
669840b8ebdSBarry Smith     /*
670840b8ebdSBarry Smith        Loop over each column associated with color adding the
671840b8ebdSBarry Smith        perturbation to the vector w3.
672840b8ebdSBarry Smith     */
673840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
674840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
675840b8ebdSBarry Smith       dx  = xx[col-start];
676840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
677aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
678840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
679840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
680840b8ebdSBarry Smith #else
681329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
682329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
683840b8ebdSBarry Smith #endif
684840b8ebdSBarry Smith       dx          *= epsilon;
685840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
6863b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
687840b8ebdSBarry Smith     }
688840b8ebdSBarry Smith     /*
689840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
690840b8ebdSBarry Smith     */
691840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
692840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
693840b8ebdSBarry Smith     /* Communicate scale to all processors */
6946831982aSBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
695840b8ebdSBarry Smith     /*
696840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
697840b8ebdSBarry Smith     */
6983b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
699840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
700840b8ebdSBarry Smith       row    = coloring->rows[k][l];
701840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
702840b8ebdSBarry Smith       y[row] *= scale[col];
703840b8ebdSBarry Smith       srow   = row + start;
704840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
705840b8ebdSBarry Smith     }
7063b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
707840b8ebdSBarry Smith   }
7083b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
709840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
710840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7113a40ed3dSBarry Smith   PetscFunctionReturn(0);
712840b8ebdSBarry Smith }
7133b28642cSBarry Smith 
7143b28642cSBarry Smith 
7153b28642cSBarry Smith 
716