xref: /petsc/src/mat/matfd/fdmatrix.c (revision 3a7fca6b988f1faea172e46abdf950477dbfaf76)
1*3a7fca6bSBarry Smith /*$Id: fdmatrix.c,v 1.87 2001/05/29 20:13:27 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 
114a2ae208SSatish Balay #undef __FUNCT__
12*3a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
13*3a7fca6bSBarry Smith int MatFDColoringSetF(MatFDColoring fd,Vec F)
14*3a7fca6bSBarry Smith {
15*3a7fca6bSBarry Smith   PetscFunctionBegin;
16*3a7fca6bSBarry Smith   fd->F = F;
17*3a7fca6bSBarry Smith   PetscFunctionReturn(0);
18*3a7fca6bSBarry Smith }
19*3a7fca6bSBarry Smith 
20*3a7fca6bSBarry Smith #undef __FUNCT__
214a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
22b0a32e0cSBarry Smith static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
239194eea9SBarry Smith {
249194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
259194eea9SBarry Smith   int           ierr,i,j;
269194eea9SBarry Smith   PetscReal     x,y;
279194eea9SBarry Smith 
289194eea9SBarry Smith   PetscFunctionBegin;
299194eea9SBarry Smith 
309194eea9SBarry Smith   /* loop over colors  */
319194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
329194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
339194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
349194eea9SBarry Smith       x = fd->columnsforrow[i][j];
35b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
369194eea9SBarry Smith     }
379194eea9SBarry Smith   }
389194eea9SBarry Smith   PetscFunctionReturn(0);
399194eea9SBarry Smith }
409194eea9SBarry Smith 
414a2ae208SSatish Balay #undef __FUNCT__
424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
43b0a32e0cSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
44005c665bSBarry Smith {
459194eea9SBarry Smith   int         ierr;
46005c665bSBarry Smith   PetscTruth  isnull;
47b0a32e0cSBarry Smith   PetscDraw   draw;
4854d96fa3SBarry Smith   PetscReal   xr,yr,xl,yl,h,w;
49005c665bSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
51b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
52b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
539194eea9SBarry Smith 
549194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
55005c665bSBarry Smith 
56005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
57005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
58b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
59b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
609194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
613a40ed3dSBarry Smith   PetscFunctionReturn(0);
62005c665bSBarry Smith }
63005c665bSBarry Smith 
644a2ae208SSatish Balay #undef __FUNCT__
654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
66bbf0e169SBarry Smith /*@C
67639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
68bbf0e169SBarry Smith 
694c49b128SBarry Smith    Collective on MatFDColoring
70fee21e36SBarry Smith 
71ef5ee4d1SLois Curfman McInnes    Input  Parameters:
72ef5ee4d1SLois Curfman McInnes +  c - the coloring context
73ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
74ef5ee4d1SLois Curfman McInnes 
7515091d37SBarry Smith    Level: intermediate
7615091d37SBarry Smith 
77b4fc646aSLois Curfman McInnes    Notes:
78b4fc646aSLois Curfman McInnes    The available visualization contexts include
79b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
80b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
81ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
82ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
83ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
84b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
85bbf0e169SBarry Smith 
869194eea9SBarry Smith    Notes:
879194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
88b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
899194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
909194eea9SBarry Smith 
91639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
92005c665bSBarry Smith 
93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
94bbf0e169SBarry Smith @*/
95b0a32e0cSBarry Smith int MatFDColoringView(MatFDColoring c,PetscViewer viewer)
96bbf0e169SBarry Smith {
97fb9695e5SSatish Balay   int               i,j,ierr;
986831982aSBarry Smith   PetscTruth        isdraw,isascii;
99f3ef73ceSBarry Smith   PetscViewerFormat format;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
102b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
103b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
104b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
1056831982aSBarry Smith   PetscCheckSameComm(c,viewer);
106bbf0e169SBarry Smith 
107fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
108b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
1090f5bd95cSBarry Smith   if (isdraw) {
110b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1110f5bd95cSBarry Smith   } else if (isascii) {
112b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
113b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
114b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
115b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
116ae09f205SBarry Smith 
117b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
118fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
119b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
120b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
121b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
122b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
123b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
124639f9d9dSBarry Smith         }
125b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
126b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
127b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
128b4fc646aSLois Curfman McInnes         }
129bbf0e169SBarry Smith       }
130bbf0e169SBarry Smith     }
131b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1325cd90555SBarry Smith   } else {
13329bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
134bbf0e169SBarry Smith   }
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
136639f9d9dSBarry Smith }
137639f9d9dSBarry Smith 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
140639f9d9dSBarry Smith /*@
141b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
142b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
143639f9d9dSBarry Smith 
144ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
145ef5ee4d1SLois Curfman McInnes 
146ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
147ef5ee4d1SLois Curfman McInnes .vb
14865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
149f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
150f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
151ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
152ef5ee4d1SLois Curfman McInnes .ve
153639f9d9dSBarry Smith 
154639f9d9dSBarry Smith    Input Parameters:
155ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
156639f9d9dSBarry Smith .  error_rel - relative error
157f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
158fee21e36SBarry Smith 
15915091d37SBarry Smith    Level: advanced
16015091d37SBarry Smith 
161b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
162b4fc646aSLois Curfman McInnes 
163b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
164639f9d9dSBarry Smith @*/
165329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
166639f9d9dSBarry Smith {
1673a40ed3dSBarry Smith   PetscFunctionBegin;
168639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
169639f9d9dSBarry Smith 
170639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
171639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
173639f9d9dSBarry Smith }
174639f9d9dSBarry Smith 
1754a2ae208SSatish Balay #undef __FUNCT__
1764a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency"
177005c665bSBarry Smith /*@
178e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
179e0907662SLois Curfman McInnes    matrices.
180005c665bSBarry Smith 
181fee21e36SBarry Smith    Collective on MatFDColoring
182fee21e36SBarry Smith 
183ef5ee4d1SLois Curfman McInnes    Input Parameters:
184ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
185ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
186ef5ee4d1SLois Curfman McInnes 
18715091d37SBarry Smith    Options Database Keys:
18815091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
18915091d37SBarry Smith 
19015091d37SBarry Smith    Level: advanced
19115091d37SBarry Smith 
192e0907662SLois Curfman McInnes    Notes:
193e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
194e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
195e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
196e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
197e0907662SLois Curfman McInnes 
198b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
199ef5ee4d1SLois Curfman McInnes 
20040964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
201005c665bSBarry Smith @*/
202005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
203005c665bSBarry Smith {
2043a40ed3dSBarry Smith   PetscFunctionBegin;
205005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
206005c665bSBarry Smith 
207005c665bSBarry Smith   matfd->freq = freq;
2083a40ed3dSBarry Smith   PetscFunctionReturn(0);
209005c665bSBarry Smith }
210005c665bSBarry Smith 
2114a2ae208SSatish Balay #undef __FUNCT__
2124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency"
213ff0cfa39SBarry Smith /*@
214ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
215ff0cfa39SBarry Smith    matrices.
216ff0cfa39SBarry Smith 
217ef5ee4d1SLois Curfman McInnes    Not Collective
218ef5ee4d1SLois Curfman McInnes 
219ff0cfa39SBarry Smith    Input Parameters:
220ff0cfa39SBarry Smith .  coloring - the coloring context
221ff0cfa39SBarry Smith 
222ff0cfa39SBarry Smith    Output Parameters:
223ff0cfa39SBarry Smith .  freq - frequency (default is 1)
224ff0cfa39SBarry Smith 
22515091d37SBarry Smith    Options Database Keys:
22615091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
22715091d37SBarry Smith 
22815091d37SBarry Smith    Level: advanced
22915091d37SBarry Smith 
230ff0cfa39SBarry Smith    Notes:
231ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
232ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
233ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
234ff0cfa39SBarry Smith    <freq> nonlinear iterations.
235ff0cfa39SBarry Smith 
236ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
237ef5ee4d1SLois Curfman McInnes 
238ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
239ff0cfa39SBarry Smith @*/
240ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
241ff0cfa39SBarry Smith {
2423a40ed3dSBarry Smith   PetscFunctionBegin;
243ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
244ff0cfa39SBarry Smith 
245ff0cfa39SBarry Smith   *freq = matfd->freq;
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
247ff0cfa39SBarry Smith }
248ff0cfa39SBarry Smith 
2494a2ae208SSatish Balay #undef __FUNCT__
2504a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
251d64ed03dSBarry Smith /*@C
252005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
253005c665bSBarry Smith 
254fee21e36SBarry Smith    Collective on MatFDColoring
255fee21e36SBarry Smith 
256ef5ee4d1SLois Curfman McInnes    Input Parameters:
257ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
258ef5ee4d1SLois Curfman McInnes .  f - the function
259ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
260ef5ee4d1SLois Curfman McInnes 
26115091d37SBarry Smith    Level: intermediate
26215091d37SBarry Smith 
263f881d145SBarry Smith    Notes:
264f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
265f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
266f881d145SBarry Smith   with the TS solvers.
267f881d145SBarry Smith 
268b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
269005c665bSBarry Smith @*/
270840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
271005c665bSBarry Smith {
2723a40ed3dSBarry Smith   PetscFunctionBegin;
273005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
274005c665bSBarry Smith 
275005c665bSBarry Smith   matfd->f    = f;
276005c665bSBarry Smith   matfd->fctx = fctx;
277005c665bSBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionReturn(0);
279005c665bSBarry Smith }
280005c665bSBarry Smith 
2814a2ae208SSatish Balay #undef __FUNCT__
2824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
283639f9d9dSBarry Smith /*@
284b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
285639f9d9dSBarry Smith    the options database.
286639f9d9dSBarry Smith 
287fee21e36SBarry Smith    Collective on MatFDColoring
288fee21e36SBarry Smith 
28965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
290ef5ee4d1SLois Curfman McInnes .vb
29165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
292f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
293f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
294ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
295ef5ee4d1SLois Curfman McInnes .ve
296ef5ee4d1SLois Curfman McInnes 
297ef5ee4d1SLois Curfman McInnes    Input Parameter:
298ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
299ef5ee4d1SLois Curfman McInnes 
300b4fc646aSLois Curfman McInnes    Options Database Keys:
301d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
302ef5ee4d1SLois Curfman McInnes            of relative error in the function)
303f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
304ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
305ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
306ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
307ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
308639f9d9dSBarry Smith 
30915091d37SBarry Smith     Level: intermediate
31015091d37SBarry Smith 
311b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
312639f9d9dSBarry Smith @*/
313639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
314639f9d9dSBarry Smith {
315186905e3SBarry Smith   int        ierr;
3163a40ed3dSBarry Smith 
3173a40ed3dSBarry Smith   PetscFunctionBegin;
318639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
319639f9d9dSBarry Smith 
320b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
321b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
322b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
323b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
324186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
325b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
326b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
327b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
328b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3293a40ed3dSBarry Smith   PetscFunctionReturn(0);
330005c665bSBarry Smith }
331005c665bSBarry Smith 
3324a2ae208SSatish Balay #undef __FUNCT__
3334a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
334005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
335005c665bSBarry Smith {
336f1af5d2fSBarry Smith   int        ierr;
337f1af5d2fSBarry Smith   PetscTruth flg;
338005c665bSBarry Smith 
3393a40ed3dSBarry Smith   PetscFunctionBegin;
340b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
341005c665bSBarry Smith   if (flg) {
342b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
343005c665bSBarry Smith   }
344b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
345ae09f205SBarry Smith   if (flg) {
346fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
347b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
348b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
349ae09f205SBarry Smith   }
350b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
351005c665bSBarry Smith   if (flg) {
352b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
353b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
354005c665bSBarry Smith   }
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
356bbf0e169SBarry Smith }
357bbf0e169SBarry Smith 
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
360bbf0e169SBarry Smith /*@C
361639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
362639f9d9dSBarry Smith    computation of Jacobians.
363bbf0e169SBarry Smith 
364ef5ee4d1SLois Curfman McInnes    Collective on Mat
365ef5ee4d1SLois Curfman McInnes 
366639f9d9dSBarry Smith    Input Parameters:
367ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
368ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
369bbf0e169SBarry Smith 
370bbf0e169SBarry Smith     Output Parameter:
371639f9d9dSBarry Smith .   color - the new coloring context
372bbf0e169SBarry Smith 
373b4fc646aSLois Curfman McInnes     Options Database Keys:
374ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
375ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
376ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
377639f9d9dSBarry Smith 
37815091d37SBarry Smith     Level: intermediate
37915091d37SBarry Smith 
3806831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
3816831982aSBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
382bbf0e169SBarry Smith @*/
383639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
384bbf0e169SBarry Smith {
385639f9d9dSBarry Smith   MatFDColoring c;
386639f9d9dSBarry Smith   MPI_Comm      comm;
387639f9d9dSBarry Smith   int           ierr,M,N;
388639f9d9dSBarry Smith 
3893a40ed3dSBarry Smith   PetscFunctionBegin;
390b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
391639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
39229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
393639f9d9dSBarry Smith 
394f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3959194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
396b0a32e0cSBarry Smith   PetscLogObjectCreate(c);
397639f9d9dSBarry Smith 
398f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
399f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
400639f9d9dSBarry Smith   } else {
40129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
402639f9d9dSBarry Smith   }
403639f9d9dSBarry Smith 
404639f9d9dSBarry Smith   c->error_rel         = 1.e-8;
405ae09f205SBarry Smith   c->umin              = 1.e-6;
406005c665bSBarry Smith   c->freq              = 1;
40740964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
40840964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
409005c665bSBarry Smith 
410005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
411639f9d9dSBarry Smith 
412639f9d9dSBarry Smith   *color = c;
413b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
415bbf0e169SBarry Smith }
416bbf0e169SBarry Smith 
4174a2ae208SSatish Balay #undef __FUNCT__
4184a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
419bbf0e169SBarry Smith /*@C
420639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
421639f9d9dSBarry Smith     via MatFDColoringCreate().
422bbf0e169SBarry Smith 
423ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
424ef5ee4d1SLois Curfman McInnes 
425b4fc646aSLois Curfman McInnes     Input Parameter:
426639f9d9dSBarry Smith .   c - coloring context
427bbf0e169SBarry Smith 
42815091d37SBarry Smith     Level: intermediate
42915091d37SBarry Smith 
430639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
431bbf0e169SBarry Smith @*/
432639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
433bbf0e169SBarry Smith {
434263760aaSBarry Smith   int i,ierr;
435bbf0e169SBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
4373a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
438d4bb536fSBarry Smith 
439639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
440606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
441606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
442606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
44330b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
444bbf0e169SBarry Smith   }
445606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
446606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
447606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
448606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
449606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
45030b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4514f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
452005c665bSBarry Smith   if (c->w1) {
453005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
454005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
455005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
456005c665bSBarry Smith   }
457b0a32e0cSBarry Smith   PetscLogObjectDestroy(c);
458639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4593a40ed3dSBarry Smith   PetscFunctionReturn(0);
460bbf0e169SBarry Smith }
46143a90d84SBarry Smith 
4624a2ae208SSatish Balay #undef __FUNCT__
4634a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
46443a90d84SBarry Smith /*@
465e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
466e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
46743a90d84SBarry Smith 
468fee21e36SBarry Smith     Collective on MatFDColoring
469fee21e36SBarry Smith 
470ef5ee4d1SLois Curfman McInnes     Input Parameters:
471ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
472ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
473ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
474ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
475ef5ee4d1SLois Curfman McInnes 
4768bba8e72SBarry Smith    Options Database Keys:
477ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
4788bba8e72SBarry Smith 
47915091d37SBarry Smith    Level: intermediate
48015091d37SBarry Smith 
48143a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
48243a90d84SBarry Smith 
48343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
48443a90d84SBarry Smith @*/
485005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
48643a90d84SBarry Smith {
4875904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
488*3a7fca6bSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
4895904e1b1SBarry Smith   Scalar        dx,mone = -1.0,*y,*xx,*w3_array;
4904f113abfSBarry Smith   Scalar        *vscale_array;
491329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
492005c665bSBarry Smith   Vec           w1,w2,w3;
493005c665bSBarry Smith   void          *fctx = coloring->fctx;
494f1af5d2fSBarry Smith   PetscTruth    flg;
495005c665bSBarry Smith 
496a2e34c3dSBarry Smith 
4973a40ed3dSBarry Smith   PetscFunctionBegin;
498e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
499e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
500e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
501e0907662SLois Curfman McInnes 
50240964ad5SBarry Smith   if (coloring->usersetsrecompute) {
50340964ad5SBarry Smith     if (!coloring->recompute) {
50440964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
50576d8deaeSBarry Smith       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
50640964ad5SBarry Smith       PetscFunctionReturn(0);
50740964ad5SBarry Smith     } else {
50840964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
50940964ad5SBarry Smith     }
51040964ad5SBarry Smith   }
51140964ad5SBarry Smith 
512b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
51376d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
51476d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
51576d8deaeSBarry Smith   } else {
51676d8deaeSBarry Smith 
517005c665bSBarry Smith     if (!coloring->w1) {
518005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
519b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w1);
520005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
521b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w2);
522005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
523b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w3);
524005c665bSBarry Smith     }
525005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
52643a90d84SBarry Smith 
5274c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
528b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
529f1af5d2fSBarry Smith     if (flg) {
530b0a32e0cSBarry Smith       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
531e0907662SLois Curfman McInnes     } else {
53243a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
533e0907662SLois Curfman McInnes     }
53443a90d84SBarry Smith 
53543a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
53643a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
537be2fbe1fSBarry Smith 
538*3a7fca6bSBarry Smith     /*
539*3a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
540*3a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
541*3a7fca6bSBarry Smith     */
542*3a7fca6bSBarry Smith     if (coloring->F) {
543*3a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
544*3a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
545*3a7fca6bSBarry Smith       if (m1 != m2) {
546*3a7fca6bSBarry Smith 	coloring->F = 0;
547*3a7fca6bSBarry Smith       }
548*3a7fca6bSBarry Smith     }
549*3a7fca6bSBarry Smith 
550be2fbe1fSBarry Smith     if (coloring->F) {
551be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
552be2fbe1fSBarry Smith       coloring->F = 0;
553be2fbe1fSBarry Smith     } else {
55443a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
555be2fbe1fSBarry Smith     }
55643a90d84SBarry Smith 
55743a90d84SBarry Smith     /*
5589782cf98SBarry Smith        Compute all the scale factors and share with other processors
5599782cf98SBarry Smith     */
5609782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
5614f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
5629782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
5639782cf98SBarry Smith       /*
5649782cf98SBarry Smith 	Loop over each column associated with color adding the
5659782cf98SBarry Smith 	perturbation to the vector w3.
5669782cf98SBarry Smith       */
5679782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
5689782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
5699782cf98SBarry Smith 	dx  = xx[col];
5709782cf98SBarry Smith 	if (dx == 0.0) dx = 1.0;
5719782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
5729782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
5739782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
5749782cf98SBarry Smith #else
5759782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
5769782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
5779782cf98SBarry Smith #endif
5789782cf98SBarry Smith 	dx                *= epsilon;
57930b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
5809782cf98SBarry Smith       }
5819782cf98SBarry Smith     }
582a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
58330b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
58430b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5859782cf98SBarry Smith 
586ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
587ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
588b0a32e0cSBarry Smith 
58930b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
59030b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
59130b34957SBarry Smith 
59230b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
5939782cf98SBarry Smith     /*
59443a90d84SBarry Smith       Loop over each color
59543a90d84SBarry Smith     */
59643a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
59743a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
59842460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
59943a90d84SBarry Smith       /*
60043a90d84SBarry Smith 	Loop over each column associated with color adding the
60143a90d84SBarry Smith 	perturbation to the vector w3.
60243a90d84SBarry Smith       */
60343a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
60443a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
60542460c72SBarry Smith 	dx  = xx[col];
606ae09f205SBarry Smith 	if (dx == 0.0) dx = 1.0;
607aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
608ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
609ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
61043a90d84SBarry Smith #else
611329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
612329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
61343a90d84SBarry Smith #endif
61443a90d84SBarry Smith 	dx            *= epsilon;
61529bbc08cSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
61642460c72SBarry Smith 	w3_array[col] += dx;
61743a90d84SBarry Smith       }
61842460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6193b28642cSBarry Smith 
62043a90d84SBarry Smith       /*
621e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
62243a90d84SBarry Smith       */
623a2e34c3dSBarry Smith 
62443a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
62543a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6269782cf98SBarry Smith 
62743a90d84SBarry Smith       /*
628e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
62943a90d84SBarry Smith       */
6303b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
63143a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
63243a90d84SBarry Smith 	row    = coloring->rows[k][l];
63343a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6345904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
63543a90d84SBarry Smith 	srow   = row + start;
63643a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
63743a90d84SBarry Smith       }
6383b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
63943a90d84SBarry Smith     }
64030b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
64142460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
64243a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64343a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64476d8deaeSBarry Smith   }
645b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
646a2e34c3dSBarry Smith 
647b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
648a2e34c3dSBarry Smith   if (flg) {
649a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
650a2e34c3dSBarry Smith   }
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
65243a90d84SBarry Smith }
653840b8ebdSBarry Smith 
6544a2ae208SSatish Balay #undef __FUNCT__
6554a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
656840b8ebdSBarry Smith /*@
657840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
658840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
659840b8ebdSBarry Smith 
660fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
661fee21e36SBarry Smith 
662ef5ee4d1SLois Curfman McInnes     Input Parameters:
6633b28642cSBarry Smith +   mat - location to store Jacobian
664ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
665ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
666ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
667ef5ee4d1SLois Curfman McInnes 
668840b8ebdSBarry Smith    Options Database Keys:
669ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
670840b8ebdSBarry Smith 
67115091d37SBarry Smith    Level: intermediate
67215091d37SBarry Smith 
673840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
674840b8ebdSBarry Smith 
675840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
676840b8ebdSBarry Smith @*/
677329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
678840b8ebdSBarry Smith {
679329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
6805904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
6815904e1b1SBarry Smith   Scalar        dx,mone = -1.0,*y,*xx,*w3_array;
6825904e1b1SBarry Smith   Scalar        *vscale_array;
683329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
684840b8ebdSBarry Smith   Vec           w1,w2,w3;
685840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
686f1af5d2fSBarry Smith   PetscTruth    flg;
687840b8ebdSBarry Smith 
6883a40ed3dSBarry Smith   PetscFunctionBegin;
689840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
690840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
691840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
692840b8ebdSBarry Smith 
693b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
694840b8ebdSBarry Smith   if (!coloring->w1) {
695840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
696b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w1);
697840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
698b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w2);
699840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
700b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w3);
701840b8ebdSBarry Smith   }
702840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
703840b8ebdSBarry Smith 
7045904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
705b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
706f1af5d2fSBarry Smith   if (flg) {
707b0a32e0cSBarry Smith     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
708840b8ebdSBarry Smith   } else {
709840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
710840b8ebdSBarry Smith   }
711840b8ebdSBarry Smith 
712840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
713840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
714840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
715840b8ebdSBarry Smith 
716840b8ebdSBarry Smith   /*
7175904e1b1SBarry Smith       Compute all the scale factors and share with other processors
718840b8ebdSBarry Smith   */
7195904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7205904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
721840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
722840b8ebdSBarry Smith     /*
723840b8ebdSBarry Smith        Loop over each column associated with color adding the
724840b8ebdSBarry Smith        perturbation to the vector w3.
725840b8ebdSBarry Smith     */
726840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
727840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7285904e1b1SBarry Smith       dx  = xx[col];
729840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
730aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
731840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
732840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
733840b8ebdSBarry Smith #else
734329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
735329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
736840b8ebdSBarry Smith #endif
737840b8ebdSBarry Smith       dx                *= epsilon;
7385904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
739840b8ebdSBarry Smith     }
7405904e1b1SBarry Smith   }
7415904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7425904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7435904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7445904e1b1SBarry Smith 
7455904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
7465904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
7475904e1b1SBarry Smith 
7485904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7495904e1b1SBarry Smith   /*
7505904e1b1SBarry Smith       Loop over each color
7515904e1b1SBarry Smith   */
7525904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
7535904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
7545904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
7555904e1b1SBarry Smith     /*
7565904e1b1SBarry Smith        Loop over each column associated with color adding the
7575904e1b1SBarry Smith        perturbation to the vector w3.
7585904e1b1SBarry Smith     */
7595904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
7605904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7615904e1b1SBarry Smith       dx  = xx[col];
7625904e1b1SBarry Smith       if (dx == 0.0) dx = 1.0;
7635904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
7645904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
7655904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
7665904e1b1SBarry Smith #else
7675904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
7685904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
7695904e1b1SBarry Smith #endif
7705904e1b1SBarry Smith       dx            *= epsilon;
7715904e1b1SBarry Smith       w3_array[col] += dx;
7725904e1b1SBarry Smith     }
7735904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7745904e1b1SBarry Smith 
775840b8ebdSBarry Smith     /*
776840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
777840b8ebdSBarry Smith     */
778840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
779840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
7805904e1b1SBarry Smith 
781840b8ebdSBarry Smith     /*
782840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
783840b8ebdSBarry Smith     */
7843b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
785840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
786840b8ebdSBarry Smith       row    = coloring->rows[k][l];
787840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
7885904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
789840b8ebdSBarry Smith       srow   = row + start;
790840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
791840b8ebdSBarry Smith     }
7923b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
793840b8ebdSBarry Smith   }
7945904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7955904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
796840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
797840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
798b0a32e0cSBarry Smith   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
7993a40ed3dSBarry Smith   PetscFunctionReturn(0);
800840b8ebdSBarry Smith }
8013b28642cSBarry Smith 
8023b28642cSBarry Smith 
8034a2ae208SSatish Balay #undef __FUNCT__
8044a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
80540964ad5SBarry Smith /*@C
80640964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
80740964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
80840964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
80940964ad5SBarry Smith      called again.
81040964ad5SBarry Smith 
81140964ad5SBarry Smith    Collective on MatFDColoring
81240964ad5SBarry Smith 
81340964ad5SBarry Smith    Input  Parameters:
81440964ad5SBarry Smith .  c - the coloring context
81540964ad5SBarry Smith 
81640964ad5SBarry Smith    Level: intermediate
81740964ad5SBarry Smith 
81840964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
81940964ad5SBarry Smith 
82040964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
82140964ad5SBarry Smith 
82240964ad5SBarry Smith .keywords: Mat, finite differences, coloring
82340964ad5SBarry Smith @*/
82440964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c)
82540964ad5SBarry Smith {
82640964ad5SBarry Smith   PetscFunctionBegin;
82740964ad5SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
82840964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
82940964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
83040964ad5SBarry Smith   PetscFunctionReturn(0);
83140964ad5SBarry Smith }
83240964ad5SBarry Smith 
8333b28642cSBarry Smith 
834