xref: /petsc/src/mat/matfd/fdmatrix.c (revision 5b8514eb5c47e8e3689c8fb9dc9fe3f0f65ee1ea)
1bbf0e169SBarry Smith /*
2639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
3639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
4bbf0e169SBarry Smith */
5bbf0e169SBarry Smith 
6e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
7bbf0e169SBarry Smith 
88ba1e511SMatthew Knepley /* Logging support */
943c77886SBarry Smith int MAT_FDCOLORING_COOKIE = 0;
108ba1e511SMatthew Knepley 
114a2ae208SSatish Balay #undef __FUNCT__
123a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
133a7fca6bSBarry Smith int MatFDColoringSetF(MatFDColoring fd,Vec F)
143a7fca6bSBarry Smith {
153a7fca6bSBarry Smith   PetscFunctionBegin;
163a7fca6bSBarry Smith   fd->F = F;
173a7fca6bSBarry Smith   PetscFunctionReturn(0);
183a7fca6bSBarry Smith }
193a7fca6bSBarry Smith 
203a7fca6bSBarry 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;
9832077d6dSBarry Smith   PetscTruth        isdraw,iascii;
99f3ef73ceSBarry Smith   PetscViewerFormat format;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
1024482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
103b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
1044482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
105c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
106bbf0e169SBarry Smith 
107fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1090f5bd95cSBarry Smith   if (isdraw) {
110b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11132077d6dSBarry Smith   } else if (iascii) {
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;
1684482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
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;
2054482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
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;
2434482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
244ff0cfa39SBarry Smith   *freq = matfd->freq;
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
246ff0cfa39SBarry Smith }
247ff0cfa39SBarry Smith 
2484a2ae208SSatish Balay #undef __FUNCT__
2494a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
250d64ed03dSBarry Smith /*@C
251005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
252005c665bSBarry Smith 
253fee21e36SBarry Smith    Collective on MatFDColoring
254fee21e36SBarry Smith 
255ef5ee4d1SLois Curfman McInnes    Input Parameters:
256ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
257ef5ee4d1SLois Curfman McInnes .  f - the function
258ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
259ef5ee4d1SLois Curfman McInnes 
26015091d37SBarry Smith    Level: intermediate
26115091d37SBarry Smith 
262f881d145SBarry Smith    Notes:
263f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
264f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
265f881d145SBarry Smith   with the TS solvers.
266f881d145SBarry Smith 
267b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
268005c665bSBarry Smith @*/
269840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
270005c665bSBarry Smith {
2713a40ed3dSBarry Smith   PetscFunctionBegin;
2724482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
273005c665bSBarry Smith   matfd->f    = f;
274005c665bSBarry Smith   matfd->fctx = fctx;
2753a40ed3dSBarry Smith   PetscFunctionReturn(0);
276005c665bSBarry Smith }
277005c665bSBarry Smith 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
280639f9d9dSBarry Smith /*@
281b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
282639f9d9dSBarry Smith    the options database.
283639f9d9dSBarry Smith 
284fee21e36SBarry Smith    Collective on MatFDColoring
285fee21e36SBarry Smith 
28665f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
287ef5ee4d1SLois Curfman McInnes .vb
28865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
289f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
290f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
291ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
292ef5ee4d1SLois Curfman McInnes .ve
293ef5ee4d1SLois Curfman McInnes 
294ef5ee4d1SLois Curfman McInnes    Input Parameter:
295ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
296ef5ee4d1SLois Curfman McInnes 
297b4fc646aSLois Curfman McInnes    Options Database Keys:
298d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
299ef5ee4d1SLois Curfman McInnes            of relative error in the function)
300f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
301ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
302ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
303ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
304ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
305639f9d9dSBarry Smith 
30615091d37SBarry Smith     Level: intermediate
30715091d37SBarry Smith 
308b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
309d1c39f55SBarry Smith 
310d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
311d1c39f55SBarry Smith 
312639f9d9dSBarry Smith @*/
313639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
314639f9d9dSBarry Smith {
315186905e3SBarry Smith   int        ierr;
3163a40ed3dSBarry Smith 
3173a40ed3dSBarry Smith   PetscFunctionBegin;
3184482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
319639f9d9dSBarry Smith 
320b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
32187828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
32287828ca2SBarry Smith     ierr = PetscOptionsReal("-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 
37315091d37SBarry Smith     Level: intermediate
37415091d37SBarry Smith 
3756831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
376d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
377d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
378d1c39f55SBarry Smith           MatFDColoringSetParameters()
379bbf0e169SBarry Smith @*/
380639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
381bbf0e169SBarry Smith {
382639f9d9dSBarry Smith   MatFDColoring c;
383639f9d9dSBarry Smith   MPI_Comm      comm;
384639f9d9dSBarry Smith   int           ierr,M,N;
385639f9d9dSBarry Smith 
3863a40ed3dSBarry Smith   PetscFunctionBegin;
387d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
388639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
38929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
390639f9d9dSBarry Smith 
391f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3929194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
393b0a32e0cSBarry Smith   PetscLogObjectCreate(c);
394639f9d9dSBarry Smith 
395f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
396f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
397639f9d9dSBarry Smith   } else {
39829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
399639f9d9dSBarry Smith   }
400639f9d9dSBarry Smith 
40177d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
40277d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
403005c665bSBarry Smith   c->freq              = 1;
40440964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
40540964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
40649b058dcSBarry Smith   c->currentcolor      = -1;
407639f9d9dSBarry Smith 
408639f9d9dSBarry Smith   *color = c;
409d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4103a40ed3dSBarry Smith   PetscFunctionReturn(0);
411bbf0e169SBarry Smith }
412bbf0e169SBarry Smith 
4134a2ae208SSatish Balay #undef __FUNCT__
4144a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
415bbf0e169SBarry Smith /*@C
416639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
417639f9d9dSBarry Smith     via MatFDColoringCreate().
418bbf0e169SBarry Smith 
419ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
420ef5ee4d1SLois Curfman McInnes 
421b4fc646aSLois Curfman McInnes     Input Parameter:
422639f9d9dSBarry Smith .   c - coloring context
423bbf0e169SBarry Smith 
42415091d37SBarry Smith     Level: intermediate
42515091d37SBarry Smith 
426639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
427bbf0e169SBarry Smith @*/
428639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
429bbf0e169SBarry Smith {
430263760aaSBarry Smith   int i,ierr;
431bbf0e169SBarry Smith 
4323a40ed3dSBarry Smith   PetscFunctionBegin;
4333a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
434d4bb536fSBarry Smith 
435639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
436606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
437606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
438606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
43930b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
440bbf0e169SBarry Smith   }
441606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
442606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
443606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
444606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
445606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
44630b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4474f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
448005c665bSBarry Smith   if (c->w1) {
449005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
450005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
451005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
452005c665bSBarry Smith   }
453b0a32e0cSBarry Smith   PetscLogObjectDestroy(c);
454639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4553a40ed3dSBarry Smith   PetscFunctionReturn(0);
456bbf0e169SBarry Smith }
45743a90d84SBarry Smith 
4584a2ae208SSatish Balay #undef __FUNCT__
45949b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
46049b058dcSBarry Smith /*@C
46149b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
46249b058dcSBarry Smith       that are currently being perturbed.
46349b058dcSBarry Smith 
46449b058dcSBarry Smith     Not Collective
46549b058dcSBarry Smith 
46649b058dcSBarry Smith     Input Parameters:
46749b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
46849b058dcSBarry Smith 
46949b058dcSBarry Smith     Output Parameters:
47049b058dcSBarry Smith +   n - the number of local columns being perturbed
47149b058dcSBarry Smith -   cols - the column indices, in global numbering
47249b058dcSBarry Smith 
47349b058dcSBarry Smith    Level: intermediate
47449b058dcSBarry Smith 
47549b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
47649b058dcSBarry Smith 
47749b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
47849b058dcSBarry Smith @*/
4794e7234bfSBarry Smith EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int *cols[])
48049b058dcSBarry Smith {
48149b058dcSBarry Smith   PetscFunctionBegin;
48249b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
48349b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
48449b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
48549b058dcSBarry Smith   } else {
48649b058dcSBarry Smith     *n = 0;
48749b058dcSBarry Smith   }
48849b058dcSBarry Smith   PetscFunctionReturn(0);
48949b058dcSBarry Smith }
49049b058dcSBarry Smith 
49149b058dcSBarry Smith #undef __FUNCT__
4924a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
49343a90d84SBarry Smith /*@
494e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
495e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49643a90d84SBarry Smith 
497fee21e36SBarry Smith     Collective on MatFDColoring
498fee21e36SBarry Smith 
499ef5ee4d1SLois Curfman McInnes     Input Parameters:
500ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
501ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
502ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
503ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
504ef5ee4d1SLois Curfman McInnes 
5058bba8e72SBarry Smith     Options Database Keys:
506b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
507b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
508b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
509b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5108bba8e72SBarry Smith 
51115091d37SBarry Smith     Level: intermediate
51215091d37SBarry Smith 
51343a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
51443a90d84SBarry Smith 
51543a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
51643a90d84SBarry Smith @*/
517005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
51843a90d84SBarry Smith {
5195904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
5203a7fca6bSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
521ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
522ea709b57SSatish Balay   PetscScalar   *vscale_array;
523329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
524005c665bSBarry Smith   Vec           w1,w2,w3;
525005c665bSBarry Smith   void          *fctx = coloring->fctx;
526f1af5d2fSBarry Smith   PetscTruth    flg;
527005c665bSBarry Smith 
528a2e34c3dSBarry Smith 
5293a40ed3dSBarry Smith   PetscFunctionBegin;
5304482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
5314482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
5324482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
533e0907662SLois Curfman McInnes 
53440964ad5SBarry Smith   if (coloring->usersetsrecompute) {
53540964ad5SBarry Smith     if (!coloring->recompute) {
53640964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
53776d8deaeSBarry Smith       PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n");
53840964ad5SBarry Smith       PetscFunctionReturn(0);
53940964ad5SBarry Smith     } else {
54040964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
54140964ad5SBarry Smith     }
54240964ad5SBarry Smith   }
54340964ad5SBarry Smith 
544d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
54576d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
54676d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
54776d8deaeSBarry Smith   } else {
54876d8deaeSBarry Smith 
549005c665bSBarry Smith     if (!coloring->w1) {
550005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
551b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w1);
552005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
553b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w2);
554005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
555b0a32e0cSBarry Smith       PetscLogObjectParent(coloring,coloring->w3);
556005c665bSBarry Smith     }
557005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
55843a90d84SBarry Smith 
5594c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
560b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
561f1af5d2fSBarry Smith     if (flg) {
562b0a32e0cSBarry Smith       PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
563e0907662SLois Curfman McInnes     } else {
56443a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
565e0907662SLois Curfman McInnes     }
56643a90d84SBarry Smith 
56743a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
56843a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
569be2fbe1fSBarry Smith 
5703a7fca6bSBarry Smith     /*
5713a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5723a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5733a7fca6bSBarry Smith     */
5743a7fca6bSBarry Smith     if (coloring->F) {
5753a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5763a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5773a7fca6bSBarry Smith       if (m1 != m2) {
5783a7fca6bSBarry Smith 	coloring->F = 0;
5793a7fca6bSBarry Smith       }
5803a7fca6bSBarry Smith     }
5813a7fca6bSBarry Smith 
582be2fbe1fSBarry Smith     if (coloring->F) {
583be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
584be2fbe1fSBarry Smith       coloring->F = 0;
585be2fbe1fSBarry Smith     } else {
58666f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
58743a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
58866f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
589be2fbe1fSBarry Smith     }
59043a90d84SBarry Smith 
59143a90d84SBarry Smith     /*
5929782cf98SBarry Smith        Compute all the scale factors and share with other processors
5939782cf98SBarry Smith     */
5949782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
5954f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
5969782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
5979782cf98SBarry Smith       /*
5989782cf98SBarry Smith 	Loop over each column associated with color adding the
5999782cf98SBarry Smith 	perturbation to the vector w3.
6009782cf98SBarry Smith       */
6019782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
6029782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
6039782cf98SBarry Smith 	dx  = xx[col];
604*5b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
6059782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6069782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
6079782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
6089782cf98SBarry Smith #else
6099782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6109782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6119782cf98SBarry Smith #endif
6129782cf98SBarry Smith 	dx                *= epsilon;
61330b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
6149782cf98SBarry Smith       }
6159782cf98SBarry Smith     }
616a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
61730b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
61830b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6199782cf98SBarry Smith 
620ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
621ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
622b0a32e0cSBarry Smith 
62330b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
62430b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
62530b34957SBarry Smith 
62630b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6279782cf98SBarry Smith     /*
62843a90d84SBarry Smith       Loop over each color
62943a90d84SBarry Smith     */
63043a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
63149b058dcSBarry Smith       coloring->currentcolor = k;
63243a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
63342460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
63443a90d84SBarry Smith       /*
63543a90d84SBarry Smith 	Loop over each column associated with color adding the
63643a90d84SBarry Smith 	perturbation to the vector w3.
63743a90d84SBarry Smith       */
63843a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
63943a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
64042460c72SBarry Smith 	dx  = xx[col];
641*5b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
642aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
643ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
644ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
64543a90d84SBarry Smith #else
646329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
647329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
64843a90d84SBarry Smith #endif
64943a90d84SBarry Smith 	dx            *= epsilon;
65029bbc08cSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
65142460c72SBarry Smith 	w3_array[col] += dx;
65243a90d84SBarry Smith       }
65342460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6543b28642cSBarry Smith 
65543a90d84SBarry Smith       /*
656e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
65743a90d84SBarry Smith       */
658a2e34c3dSBarry Smith 
65966f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66043a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
66166f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66243a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6639782cf98SBarry Smith 
66443a90d84SBarry Smith       /*
665e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
66643a90d84SBarry Smith       */
6673b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
66843a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
66943a90d84SBarry Smith 	row    = coloring->rows[k][l];
67043a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6715904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
67243a90d84SBarry Smith 	srow   = row + start;
67343a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
67443a90d84SBarry Smith       }
6753b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
67643a90d84SBarry Smith     }
67749b058dcSBarry Smith     coloring->currentcolor = -1;
67830b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
67942460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
68043a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68143a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68276d8deaeSBarry Smith   }
683d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
684a2e34c3dSBarry Smith 
685b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
686a2e34c3dSBarry Smith   if (flg) {
687a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
688a2e34c3dSBarry Smith   }
689b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
690b9722508SLois Curfman McInnes 
6913a40ed3dSBarry Smith   PetscFunctionReturn(0);
69243a90d84SBarry Smith }
693840b8ebdSBarry Smith 
6944a2ae208SSatish Balay #undef __FUNCT__
6954a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
696840b8ebdSBarry Smith /*@
697840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
698840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
699840b8ebdSBarry Smith 
700fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
701fee21e36SBarry Smith 
702ef5ee4d1SLois Curfman McInnes     Input Parameters:
7033b28642cSBarry Smith +   mat - location to store Jacobian
704ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
705ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
706ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
707ef5ee4d1SLois Curfman McInnes 
708840b8ebdSBarry Smith    Options Database Keys:
709ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
710840b8ebdSBarry Smith 
71115091d37SBarry Smith    Level: intermediate
71215091d37SBarry Smith 
713840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
714840b8ebdSBarry Smith 
715840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
716840b8ebdSBarry Smith @*/
717329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
718840b8ebdSBarry Smith {
719329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
7205904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
721ea709b57SSatish Balay   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
722ea709b57SSatish Balay   PetscScalar   *vscale_array;
723329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
724840b8ebdSBarry Smith   Vec           w1,w2,w3;
725840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
726f1af5d2fSBarry Smith   PetscTruth    flg;
727840b8ebdSBarry Smith 
7283a40ed3dSBarry Smith   PetscFunctionBegin;
7294482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
7304482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
7314482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
732840b8ebdSBarry Smith 
733d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
734840b8ebdSBarry Smith   if (!coloring->w1) {
735840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
736b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w1);
737840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
738b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w2);
739840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
740b0a32e0cSBarry Smith     PetscLogObjectParent(coloring,coloring->w3);
741840b8ebdSBarry Smith   }
742840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
743840b8ebdSBarry Smith 
7445904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
745b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
746f1af5d2fSBarry Smith   if (flg) {
747b0a32e0cSBarry Smith     PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
748840b8ebdSBarry Smith   } else {
749840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
750840b8ebdSBarry Smith   }
751840b8ebdSBarry Smith 
752840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
753840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
75466f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
755840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
75666f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
757840b8ebdSBarry Smith 
758840b8ebdSBarry Smith   /*
7595904e1b1SBarry Smith       Compute all the scale factors and share with other processors
760840b8ebdSBarry Smith   */
7615904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7625904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
763840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
764840b8ebdSBarry Smith     /*
765840b8ebdSBarry Smith        Loop over each column associated with color adding the
766840b8ebdSBarry Smith        perturbation to the vector w3.
767840b8ebdSBarry Smith     */
768840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
769840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7705904e1b1SBarry Smith       dx  = xx[col];
771*5b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
772aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
773840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
774840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
775840b8ebdSBarry Smith #else
776329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
777329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
778840b8ebdSBarry Smith #endif
779840b8ebdSBarry Smith       dx                *= epsilon;
7805904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
781840b8ebdSBarry Smith     }
7825904e1b1SBarry Smith   }
7835904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7845904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7855904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7865904e1b1SBarry Smith 
7875904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
7885904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
7895904e1b1SBarry Smith 
7905904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7915904e1b1SBarry Smith   /*
7925904e1b1SBarry Smith       Loop over each color
7935904e1b1SBarry Smith   */
7945904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
7955904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
7965904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
7975904e1b1SBarry Smith     /*
7985904e1b1SBarry Smith        Loop over each column associated with color adding the
7995904e1b1SBarry Smith        perturbation to the vector w3.
8005904e1b1SBarry Smith     */
8015904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8025904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8035904e1b1SBarry Smith       dx  = xx[col];
804*5b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8055904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8065904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8075904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8085904e1b1SBarry Smith #else
8095904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8105904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8115904e1b1SBarry Smith #endif
8125904e1b1SBarry Smith       dx            *= epsilon;
8135904e1b1SBarry Smith       w3_array[col] += dx;
8145904e1b1SBarry Smith     }
8155904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8165904e1b1SBarry Smith 
817840b8ebdSBarry Smith     /*
818840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
819840b8ebdSBarry Smith     */
82066f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
821840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
82266f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
823840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
8245904e1b1SBarry Smith 
825840b8ebdSBarry Smith     /*
826840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
827840b8ebdSBarry Smith     */
8283b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
829840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
830840b8ebdSBarry Smith       row    = coloring->rows[k][l];
831840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8325904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
833840b8ebdSBarry Smith       srow   = row + start;
834840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
835840b8ebdSBarry Smith     }
8363b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
837840b8ebdSBarry Smith   }
8385904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8395904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
840840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
841840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
842d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8433a40ed3dSBarry Smith   PetscFunctionReturn(0);
844840b8ebdSBarry Smith }
8453b28642cSBarry Smith 
8463b28642cSBarry Smith 
8474a2ae208SSatish Balay #undef __FUNCT__
8484a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
84940964ad5SBarry Smith /*@C
85040964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
85140964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
85240964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
85340964ad5SBarry Smith      called again.
85440964ad5SBarry Smith 
85540964ad5SBarry Smith    Collective on MatFDColoring
85640964ad5SBarry Smith 
85740964ad5SBarry Smith    Input  Parameters:
85840964ad5SBarry Smith .  c - the coloring context
85940964ad5SBarry Smith 
86040964ad5SBarry Smith    Level: intermediate
86140964ad5SBarry Smith 
86240964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
86340964ad5SBarry Smith 
86440964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
86540964ad5SBarry Smith 
86640964ad5SBarry Smith .keywords: Mat, finite differences, coloring
86740964ad5SBarry Smith @*/
86840964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c)
86940964ad5SBarry Smith {
87040964ad5SBarry Smith   PetscFunctionBegin;
8714482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
87240964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
87340964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
87440964ad5SBarry Smith   PetscFunctionReturn(0);
87540964ad5SBarry Smith }
87640964ad5SBarry Smith 
8773b28642cSBarry Smith 
878