xref: /petsc/src/mat/matfd/fdmatrix.c (revision d5ba7fb7d3f82433fc93946f26018f4b1c7683c8)
173f4d377SMatthew Knepley /*$Id: fdmatrix.c,v 1.92 2001/08/21 21:03:06 bsmith Exp $*/
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 
118ba1e511SMatthew Knepley /* Logging support */
128ba1e511SMatthew Knepley int MAT_FDCOLORING_COOKIE;
138ba1e511SMatthew Knepley 
144a2ae208SSatish Balay #undef __FUNCT__
153a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
163a7fca6bSBarry Smith int MatFDColoringSetF(MatFDColoring fd,Vec F)
173a7fca6bSBarry Smith {
183a7fca6bSBarry Smith   PetscFunctionBegin;
193a7fca6bSBarry Smith   fd->F = F;
203a7fca6bSBarry Smith   PetscFunctionReturn(0);
213a7fca6bSBarry Smith }
223a7fca6bSBarry Smith 
233a7fca6bSBarry Smith #undef __FUNCT__
244a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
25b0a32e0cSBarry Smith static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
269194eea9SBarry Smith {
279194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
289194eea9SBarry Smith   int           ierr,i,j;
299194eea9SBarry Smith   PetscReal     x,y;
309194eea9SBarry Smith 
319194eea9SBarry Smith   PetscFunctionBegin;
329194eea9SBarry Smith 
339194eea9SBarry Smith   /* loop over colors  */
349194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
359194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
369194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
379194eea9SBarry Smith       x = fd->columnsforrow[i][j];
38b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
399194eea9SBarry Smith     }
409194eea9SBarry Smith   }
419194eea9SBarry Smith   PetscFunctionReturn(0);
429194eea9SBarry Smith }
439194eea9SBarry Smith 
444a2ae208SSatish Balay #undef __FUNCT__
454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
46b0a32e0cSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
47005c665bSBarry Smith {
489194eea9SBarry Smith   int         ierr;
49005c665bSBarry Smith   PetscTruth  isnull;
50b0a32e0cSBarry Smith   PetscDraw   draw;
5154d96fa3SBarry Smith   PetscReal   xr,yr,xl,yl,h,w;
52005c665bSBarry Smith 
533a40ed3dSBarry Smith   PetscFunctionBegin;
54b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
55b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
569194eea9SBarry Smith 
579194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
58005c665bSBarry Smith 
59005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
60005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
61b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
62b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
639194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
643a40ed3dSBarry Smith   PetscFunctionReturn(0);
65005c665bSBarry Smith }
66005c665bSBarry Smith 
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
69bbf0e169SBarry Smith /*@C
70639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
71bbf0e169SBarry Smith 
724c49b128SBarry Smith    Collective on MatFDColoring
73fee21e36SBarry Smith 
74ef5ee4d1SLois Curfman McInnes    Input  Parameters:
75ef5ee4d1SLois Curfman McInnes +  c - the coloring context
76ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
77ef5ee4d1SLois Curfman McInnes 
7815091d37SBarry Smith    Level: intermediate
7915091d37SBarry Smith 
80b4fc646aSLois Curfman McInnes    Notes:
81b4fc646aSLois Curfman McInnes    The available visualization contexts include
82b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
83b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
84ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
85ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
86ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
87b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
88bbf0e169SBarry Smith 
899194eea9SBarry Smith    Notes:
909194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
91b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
929194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
939194eea9SBarry Smith 
94639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
95005c665bSBarry Smith 
96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
97bbf0e169SBarry Smith @*/
98b0a32e0cSBarry Smith int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
99bbf0e169SBarry Smith {
100fb9695e5SSatish Balay   int               i,j,ierr;
1016831982aSBarry Smith   PetscTruth        isdraw,isascii;
102f3ef73ceSBarry Smith   PetscViewerFormat format;
103bbf0e169SBarry Smith 
1043a40ed3dSBarry Smith   PetscFunctionBegin;
105b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
106b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
107b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
1086831982aSBarry Smith   PetscCheckSameComm(c,viewer);
109bbf0e169SBarry Smith 
110fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
111b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
1120f5bd95cSBarry Smith   if (isdraw) {
113b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1140f5bd95cSBarry Smith   } else if (isascii) {
115b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
116b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
117b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
118b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
119ae09f205SBarry Smith 
120b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
121fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
122b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
123b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
124b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
125b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
126b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
127639f9d9dSBarry Smith         }
128b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
129b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
130b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
131b4fc646aSLois Curfman McInnes         }
132bbf0e169SBarry Smith       }
133bbf0e169SBarry Smith     }
134b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1355cd90555SBarry Smith   } else {
13629bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
137bbf0e169SBarry Smith   }
1383a40ed3dSBarry Smith   PetscFunctionReturn(0);
139639f9d9dSBarry Smith }
140639f9d9dSBarry Smith 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
143639f9d9dSBarry Smith /*@
144b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
145b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
146639f9d9dSBarry Smith 
147ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
148ef5ee4d1SLois Curfman McInnes 
149ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
150ef5ee4d1SLois Curfman McInnes .vb
15165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
152f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
153f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
154ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
155ef5ee4d1SLois Curfman McInnes .ve
156639f9d9dSBarry Smith 
157639f9d9dSBarry Smith    Input Parameters:
158ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
159639f9d9dSBarry Smith .  error_rel - relative error
160f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
161fee21e36SBarry Smith 
16215091d37SBarry Smith    Level: advanced
16315091d37SBarry Smith 
164b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
165b4fc646aSLois Curfman McInnes 
166b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
167639f9d9dSBarry Smith @*/
168329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169639f9d9dSBarry Smith {
1703a40ed3dSBarry Smith   PetscFunctionBegin;
171639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
172639f9d9dSBarry Smith 
173639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
174639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
176639f9d9dSBarry Smith }
177639f9d9dSBarry Smith 
1784a2ae208SSatish Balay #undef __FUNCT__
1794a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency"
180005c665bSBarry Smith /*@
181e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
182e0907662SLois Curfman McInnes    matrices.
183005c665bSBarry Smith 
184fee21e36SBarry Smith    Collective on MatFDColoring
185fee21e36SBarry Smith 
186ef5ee4d1SLois Curfman McInnes    Input Parameters:
187ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
188ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
189ef5ee4d1SLois Curfman McInnes 
19015091d37SBarry Smith    Options Database Keys:
19115091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
19215091d37SBarry Smith 
19315091d37SBarry Smith    Level: advanced
19415091d37SBarry Smith 
195e0907662SLois Curfman McInnes    Notes:
196e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
197e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
198e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
199e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
200e0907662SLois Curfman McInnes 
201b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
202ef5ee4d1SLois Curfman McInnes 
20340964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
204005c665bSBarry Smith @*/
205005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
206005c665bSBarry Smith {
2073a40ed3dSBarry Smith   PetscFunctionBegin;
208005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
209005c665bSBarry Smith 
210005c665bSBarry Smith   matfd->freq = freq;
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
212005c665bSBarry Smith }
213005c665bSBarry Smith 
2144a2ae208SSatish Balay #undef __FUNCT__
2154a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency"
216ff0cfa39SBarry Smith /*@
217ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
218ff0cfa39SBarry Smith    matrices.
219ff0cfa39SBarry Smith 
220ef5ee4d1SLois Curfman McInnes    Not Collective
221ef5ee4d1SLois Curfman McInnes 
222ff0cfa39SBarry Smith    Input Parameters:
223ff0cfa39SBarry Smith .  coloring - the coloring context
224ff0cfa39SBarry Smith 
225ff0cfa39SBarry Smith    Output Parameters:
226ff0cfa39SBarry Smith .  freq - frequency (default is 1)
227ff0cfa39SBarry Smith 
22815091d37SBarry Smith    Options Database Keys:
22915091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
23015091d37SBarry Smith 
23115091d37SBarry Smith    Level: advanced
23215091d37SBarry Smith 
233ff0cfa39SBarry Smith    Notes:
234ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
235ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
236ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
237ff0cfa39SBarry Smith    <freq> nonlinear iterations.
238ff0cfa39SBarry Smith 
239ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
240ef5ee4d1SLois Curfman McInnes 
241ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
242ff0cfa39SBarry Smith @*/
243ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
244ff0cfa39SBarry Smith {
2453a40ed3dSBarry Smith   PetscFunctionBegin;
246ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
247ff0cfa39SBarry Smith 
248ff0cfa39SBarry Smith   *freq = matfd->freq;
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
250ff0cfa39SBarry Smith }
251ff0cfa39SBarry Smith 
2524a2ae208SSatish Balay #undef __FUNCT__
2534a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
254d64ed03dSBarry Smith /*@C
255005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
256005c665bSBarry Smith 
257fee21e36SBarry Smith    Collective on MatFDColoring
258fee21e36SBarry Smith 
259ef5ee4d1SLois Curfman McInnes    Input Parameters:
260ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
261ef5ee4d1SLois Curfman McInnes .  f - the function
262ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
263ef5ee4d1SLois Curfman McInnes 
26415091d37SBarry Smith    Level: intermediate
26515091d37SBarry Smith 
266f881d145SBarry Smith    Notes:
267f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
268f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
269f881d145SBarry Smith   with the TS solvers.
270f881d145SBarry Smith 
271b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
272005c665bSBarry Smith @*/
273840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
274005c665bSBarry Smith {
2753a40ed3dSBarry Smith   PetscFunctionBegin;
276005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
277005c665bSBarry Smith 
278005c665bSBarry Smith   matfd->f    = f;
279005c665bSBarry Smith   matfd->fctx = fctx;
280005c665bSBarry Smith 
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
282005c665bSBarry Smith }
283005c665bSBarry Smith 
2844a2ae208SSatish Balay #undef __FUNCT__
2854a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
286639f9d9dSBarry Smith /*@
287b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
288639f9d9dSBarry Smith    the options database.
289639f9d9dSBarry Smith 
290fee21e36SBarry Smith    Collective on MatFDColoring
291fee21e36SBarry Smith 
29265f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
293ef5ee4d1SLois Curfman McInnes .vb
29465f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
295f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
296f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
297ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
298ef5ee4d1SLois Curfman McInnes .ve
299ef5ee4d1SLois Curfman McInnes 
300ef5ee4d1SLois Curfman McInnes    Input Parameter:
301ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
302ef5ee4d1SLois Curfman McInnes 
303b4fc646aSLois Curfman McInnes    Options Database Keys:
304d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
305ef5ee4d1SLois Curfman McInnes            of relative error in the function)
306f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
307ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
308ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
309ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
310ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
311639f9d9dSBarry Smith 
31215091d37SBarry Smith     Level: intermediate
31315091d37SBarry Smith 
314b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
315d1c39f55SBarry Smith 
316d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
317d1c39f55SBarry Smith 
318639f9d9dSBarry Smith @*/
319639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
320639f9d9dSBarry Smith {
321186905e3SBarry Smith   int        ierr;
3223a40ed3dSBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
324639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
325639f9d9dSBarry Smith 
326b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
32787828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
32887828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
329b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
330186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
331b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
332b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
333b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
334b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3353a40ed3dSBarry Smith   PetscFunctionReturn(0);
336005c665bSBarry Smith }
337005c665bSBarry Smith 
3384a2ae208SSatish Balay #undef __FUNCT__
3394a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
340005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
341005c665bSBarry Smith {
342f1af5d2fSBarry Smith   int        ierr;
343f1af5d2fSBarry Smith   PetscTruth flg;
344005c665bSBarry Smith 
3453a40ed3dSBarry Smith   PetscFunctionBegin;
346b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
347005c665bSBarry Smith   if (flg) {
348b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
349005c665bSBarry Smith   }
350b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
351ae09f205SBarry Smith   if (flg) {
352fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
353b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
354b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
355ae09f205SBarry Smith   }
356b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
357005c665bSBarry Smith   if (flg) {
358b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
359b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
360005c665bSBarry Smith   }
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
362bbf0e169SBarry Smith }
363bbf0e169SBarry Smith 
3644a2ae208SSatish Balay #undef __FUNCT__
3654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
366bbf0e169SBarry Smith /*@C
367639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
368639f9d9dSBarry Smith    computation of Jacobians.
369bbf0e169SBarry Smith 
370ef5ee4d1SLois Curfman McInnes    Collective on Mat
371ef5ee4d1SLois Curfman McInnes 
372639f9d9dSBarry Smith    Input Parameters:
373ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
374ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
375bbf0e169SBarry Smith 
376bbf0e169SBarry Smith     Output Parameter:
377639f9d9dSBarry Smith .   color - the new coloring context
378bbf0e169SBarry Smith 
379b4fc646aSLois Curfman McInnes     Options Database Keys:
380ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
381ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
382ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
383639f9d9dSBarry Smith 
38415091d37SBarry Smith     Level: intermediate
38515091d37SBarry Smith 
3866831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
387d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
388d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
389d1c39f55SBarry Smith           MatFDColoringSetParameters()
390bbf0e169SBarry Smith @*/
391639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
392bbf0e169SBarry Smith {
393639f9d9dSBarry Smith   MatFDColoring c;
394639f9d9dSBarry Smith   MPI_Comm      comm;
395639f9d9dSBarry Smith   int           ierr,M,N;
396639f9d9dSBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionBegin;
398*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
399639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
40029bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
401639f9d9dSBarry Smith 
402f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4039194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
404b0a32e0cSBarry Smith   PetscLogObjectCreate(c);
405639f9d9dSBarry Smith 
406f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
407f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
408639f9d9dSBarry Smith   } else {
40929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
410639f9d9dSBarry Smith   }
411639f9d9dSBarry Smith 
41277d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
41377d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
414005c665bSBarry Smith   c->freq              = 1;
41540964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
41640964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
417005c665bSBarry Smith 
418005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
419639f9d9dSBarry Smith 
420639f9d9dSBarry Smith   *color = c;
421*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4223a40ed3dSBarry Smith   PetscFunctionReturn(0);
423bbf0e169SBarry Smith }
424bbf0e169SBarry Smith 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
427bbf0e169SBarry Smith /*@C
428639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
429639f9d9dSBarry Smith     via MatFDColoringCreate().
430bbf0e169SBarry Smith 
431ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
432ef5ee4d1SLois Curfman McInnes 
433b4fc646aSLois Curfman McInnes     Input Parameter:
434639f9d9dSBarry Smith .   c - coloring context
435bbf0e169SBarry Smith 
43615091d37SBarry Smith     Level: intermediate
43715091d37SBarry Smith 
438639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
439bbf0e169SBarry Smith @*/
440639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
441bbf0e169SBarry Smith {
442263760aaSBarry Smith   int i,ierr;
443bbf0e169SBarry Smith 
4443a40ed3dSBarry Smith   PetscFunctionBegin;
4453a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
446d4bb536fSBarry Smith 
447639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
448606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
449606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
450606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
45130b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
452bbf0e169SBarry Smith   }
453606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
454606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
455606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
456606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
457606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
45830b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4594f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
460005c665bSBarry Smith   if (c->w1) {
461005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
462005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
463005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
464005c665bSBarry Smith   }
465b0a32e0cSBarry Smith   PetscLogObjectDestroy(c);
466639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4673a40ed3dSBarry Smith   PetscFunctionReturn(0);
468bbf0e169SBarry Smith }
46943a90d84SBarry Smith 
4704a2ae208SSatish Balay #undef __FUNCT__
4714a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
47243a90d84SBarry Smith /*@
473e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
474e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
47543a90d84SBarry Smith 
476fee21e36SBarry Smith     Collective on MatFDColoring
477fee21e36SBarry Smith 
478ef5ee4d1SLois Curfman McInnes     Input Parameters:
479ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
480ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
481ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
482ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
483ef5ee4d1SLois Curfman McInnes 
4848bba8e72SBarry Smith    Options Database Keys:
485ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
4868bba8e72SBarry Smith 
48715091d37SBarry Smith    Level: intermediate
48815091d37SBarry Smith 
48943a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
49043a90d84SBarry Smith 
49143a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
49243a90d84SBarry Smith @*/
493005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
49443a90d84SBarry Smith {
4955904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
4963a7fca6bSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
497ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
498ea709b57SSatish Balay   PetscScalar   *vscale_array;
499329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
500005c665bSBarry Smith   Vec           w1,w2,w3;
501005c665bSBarry Smith   void          *fctx = coloring->fctx;
502f1af5d2fSBarry Smith   PetscTruth    flg;
503005c665bSBarry Smith 
504a2e34c3dSBarry Smith 
5053a40ed3dSBarry Smith   PetscFunctionBegin;
506e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
507e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
508e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
509e0907662SLois Curfman McInnes 
51040964ad5SBarry Smith   if (coloring->usersetsrecompute) {
51140964ad5SBarry Smith     if (!coloring->recompute) {
51240964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
51376d8deaeSBarry Smith       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
51440964ad5SBarry Smith       PetscFunctionReturn(0);
51540964ad5SBarry Smith     } else {
51640964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
51740964ad5SBarry Smith     }
51840964ad5SBarry Smith   }
51940964ad5SBarry Smith 
520*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
52176d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
52276d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
52376d8deaeSBarry Smith   } else {
52476d8deaeSBarry Smith 
525005c665bSBarry Smith     if (!coloring->w1) {
526005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
527b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w1);
528005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
529b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w2);
530005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
531b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w3);
532005c665bSBarry Smith     }
533005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
53443a90d84SBarry Smith 
5354c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
536b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
537f1af5d2fSBarry Smith     if (flg) {
538b0a32e0cSBarry Smith       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
539e0907662SLois Curfman McInnes     } else {
54043a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
541e0907662SLois Curfman McInnes     }
54243a90d84SBarry Smith 
54343a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
54443a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
545be2fbe1fSBarry Smith 
5463a7fca6bSBarry Smith     /*
5473a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5483a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5493a7fca6bSBarry Smith     */
5503a7fca6bSBarry Smith     if (coloring->F) {
5513a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5523a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5533a7fca6bSBarry Smith       if (m1 != m2) {
5543a7fca6bSBarry Smith 	coloring->F = 0;
5553a7fca6bSBarry Smith       }
5563a7fca6bSBarry Smith     }
5573a7fca6bSBarry Smith 
558be2fbe1fSBarry Smith     if (coloring->F) {
559be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
560be2fbe1fSBarry Smith       coloring->F = 0;
561be2fbe1fSBarry Smith     } else {
56243a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
563be2fbe1fSBarry Smith     }
56443a90d84SBarry Smith 
56543a90d84SBarry Smith     /*
5669782cf98SBarry Smith        Compute all the scale factors and share with other processors
5679782cf98SBarry Smith     */
5689782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
5694f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
5709782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
5719782cf98SBarry Smith       /*
5729782cf98SBarry Smith 	Loop over each column associated with color adding the
5739782cf98SBarry Smith 	perturbation to the vector w3.
5749782cf98SBarry Smith       */
5759782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
5769782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
5779782cf98SBarry Smith 	dx  = xx[col];
5789782cf98SBarry Smith 	if (dx == 0.0) dx = 1.0;
5799782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
5809782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
5819782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
5829782cf98SBarry Smith #else
5839782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
5849782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
5859782cf98SBarry Smith #endif
5869782cf98SBarry Smith 	dx                *= epsilon;
58730b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
5889782cf98SBarry Smith       }
5899782cf98SBarry Smith     }
590a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
59130b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
59230b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5939782cf98SBarry Smith 
594ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
595ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
596b0a32e0cSBarry Smith 
59730b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
59830b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
59930b34957SBarry Smith 
60030b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6019782cf98SBarry Smith     /*
60243a90d84SBarry Smith       Loop over each color
60343a90d84SBarry Smith     */
60443a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
60543a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
60642460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
60743a90d84SBarry Smith       /*
60843a90d84SBarry Smith 	Loop over each column associated with color adding the
60943a90d84SBarry Smith 	perturbation to the vector w3.
61043a90d84SBarry Smith       */
61143a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
61243a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
61342460c72SBarry Smith 	dx  = xx[col];
614ae09f205SBarry Smith 	if (dx == 0.0) dx = 1.0;
615aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
616ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
617ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
61843a90d84SBarry Smith #else
619329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
620329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
62143a90d84SBarry Smith #endif
62243a90d84SBarry Smith 	dx            *= epsilon;
62329bbc08cSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
62442460c72SBarry Smith 	w3_array[col] += dx;
62543a90d84SBarry Smith       }
62642460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6273b28642cSBarry Smith 
62843a90d84SBarry Smith       /*
629e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
63043a90d84SBarry Smith       */
631a2e34c3dSBarry Smith 
63243a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
63343a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6349782cf98SBarry Smith 
63543a90d84SBarry Smith       /*
636e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
63743a90d84SBarry Smith       */
6383b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
63943a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
64043a90d84SBarry Smith 	row    = coloring->rows[k][l];
64143a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6425904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
64343a90d84SBarry Smith 	srow   = row + start;
64443a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
64543a90d84SBarry Smith       }
6463b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
64743a90d84SBarry Smith     }
64830b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
64942460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
65043a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65143a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65276d8deaeSBarry Smith   }
653*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
654a2e34c3dSBarry Smith 
655b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
656a2e34c3dSBarry Smith   if (flg) {
657a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
658a2e34c3dSBarry Smith   }
6593a40ed3dSBarry Smith   PetscFunctionReturn(0);
66043a90d84SBarry Smith }
661840b8ebdSBarry Smith 
6624a2ae208SSatish Balay #undef __FUNCT__
6634a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
664840b8ebdSBarry Smith /*@
665840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
666840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
667840b8ebdSBarry Smith 
668fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
669fee21e36SBarry Smith 
670ef5ee4d1SLois Curfman McInnes     Input Parameters:
6713b28642cSBarry Smith +   mat - location to store Jacobian
672ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
673ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
674ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
675ef5ee4d1SLois Curfman McInnes 
676840b8ebdSBarry Smith    Options Database Keys:
677ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
678840b8ebdSBarry Smith 
67915091d37SBarry Smith    Level: intermediate
68015091d37SBarry Smith 
681840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
682840b8ebdSBarry Smith 
683840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
684840b8ebdSBarry Smith @*/
685329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
686840b8ebdSBarry Smith {
687329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
6885904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
689ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
690ea709b57SSatish Balay   PetscScalar   *vscale_array;
691329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
692840b8ebdSBarry Smith   Vec           w1,w2,w3;
693840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
694f1af5d2fSBarry Smith   PetscTruth    flg;
695840b8ebdSBarry Smith 
6963a40ed3dSBarry Smith   PetscFunctionBegin;
697840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
698840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
699840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
700840b8ebdSBarry Smith 
701*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
702840b8ebdSBarry Smith   if (!coloring->w1) {
703840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
704b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w1);
705840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
706b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w2);
707840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
708b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w3);
709840b8ebdSBarry Smith   }
710840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
711840b8ebdSBarry Smith 
7125904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
713b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
714f1af5d2fSBarry Smith   if (flg) {
715b0a32e0cSBarry Smith     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
716840b8ebdSBarry Smith   } else {
717840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
718840b8ebdSBarry Smith   }
719840b8ebdSBarry Smith 
720840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
721840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
722840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
723840b8ebdSBarry Smith 
724840b8ebdSBarry Smith   /*
7255904e1b1SBarry Smith       Compute all the scale factors and share with other processors
726840b8ebdSBarry Smith   */
7275904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7285904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
729840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
730840b8ebdSBarry Smith     /*
731840b8ebdSBarry Smith        Loop over each column associated with color adding the
732840b8ebdSBarry Smith        perturbation to the vector w3.
733840b8ebdSBarry Smith     */
734840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
735840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7365904e1b1SBarry Smith       dx  = xx[col];
737840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
738aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
739840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
740840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
741840b8ebdSBarry Smith #else
742329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
743329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
744840b8ebdSBarry Smith #endif
745840b8ebdSBarry Smith       dx                *= epsilon;
7465904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
747840b8ebdSBarry Smith     }
7485904e1b1SBarry Smith   }
7495904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7505904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7515904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7525904e1b1SBarry Smith 
7535904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
7545904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
7555904e1b1SBarry Smith 
7565904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7575904e1b1SBarry Smith   /*
7585904e1b1SBarry Smith       Loop over each color
7595904e1b1SBarry Smith   */
7605904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
7615904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
7625904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
7635904e1b1SBarry Smith     /*
7645904e1b1SBarry Smith        Loop over each column associated with color adding the
7655904e1b1SBarry Smith        perturbation to the vector w3.
7665904e1b1SBarry Smith     */
7675904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
7685904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7695904e1b1SBarry Smith       dx  = xx[col];
7705904e1b1SBarry Smith       if (dx == 0.0) dx = 1.0;
7715904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
7725904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
7735904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
7745904e1b1SBarry Smith #else
7755904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
7765904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
7775904e1b1SBarry Smith #endif
7785904e1b1SBarry Smith       dx            *= epsilon;
7795904e1b1SBarry Smith       w3_array[col] += dx;
7805904e1b1SBarry Smith     }
7815904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7825904e1b1SBarry Smith 
783840b8ebdSBarry Smith     /*
784840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
785840b8ebdSBarry Smith     */
786840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
787840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
7885904e1b1SBarry Smith 
789840b8ebdSBarry Smith     /*
790840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
791840b8ebdSBarry Smith     */
7923b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
793840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
794840b8ebdSBarry Smith       row    = coloring->rows[k][l];
795840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
7965904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
797840b8ebdSBarry Smith       srow   = row + start;
798840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
799840b8ebdSBarry Smith     }
8003b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
801840b8ebdSBarry Smith   }
8025904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8035904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
804840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
805840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
806*d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8073a40ed3dSBarry Smith   PetscFunctionReturn(0);
808840b8ebdSBarry Smith }
8093b28642cSBarry Smith 
8103b28642cSBarry Smith 
8114a2ae208SSatish Balay #undef __FUNCT__
8124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
81340964ad5SBarry Smith /*@C
81440964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
81540964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
81640964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
81740964ad5SBarry Smith      called again.
81840964ad5SBarry Smith 
81940964ad5SBarry Smith    Collective on MatFDColoring
82040964ad5SBarry Smith 
82140964ad5SBarry Smith    Input  Parameters:
82240964ad5SBarry Smith .  c - the coloring context
82340964ad5SBarry Smith 
82440964ad5SBarry Smith    Level: intermediate
82540964ad5SBarry Smith 
82640964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
82740964ad5SBarry Smith 
82840964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
82940964ad5SBarry Smith 
83040964ad5SBarry Smith .keywords: Mat, finite differences, coloring
83140964ad5SBarry Smith @*/
83240964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c)
83340964ad5SBarry Smith {
83440964ad5SBarry Smith   PetscFunctionBegin;
83540964ad5SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
83640964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
83740964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
83840964ad5SBarry Smith   PetscFunctionReturn(0);
83940964ad5SBarry Smith }
84040964ad5SBarry Smith 
8413b28642cSBarry Smith 
842