xref: /petsc/src/mat/matfd/fdmatrix.c (revision efb308899f58eff35e31f4136d1431f0c81bf4bc)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
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 
8e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
9bbf0e169SBarry Smith 
108ba1e511SMatthew Knepley /* Logging support */
11be1d678aSKris Buschelman PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE = 0;
128ba1e511SMatthew Knepley 
134a2ae208SSatish Balay #undef __FUNCT__
143a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
15be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F)
163a7fca6bSBarry Smith {
173a7fca6bSBarry Smith   PetscFunctionBegin;
183a7fca6bSBarry Smith   fd->F = F;
193a7fca6bSBarry Smith   PetscFunctionReturn(0);
203a7fca6bSBarry Smith }
213a7fca6bSBarry Smith 
223a7fca6bSBarry Smith #undef __FUNCT__
234a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
246849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
259194eea9SBarry Smith {
269194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
27dfbe8321SBarry Smith   PetscErrorCode ierr;
2838baddfdSBarry Smith   PetscInt       i,j;
299194eea9SBarry Smith   PetscReal      x,y;
309194eea9SBarry Smith 
319194eea9SBarry Smith   PetscFunctionBegin;
329194eea9SBarry Smith 
339194eea9SBarry Smith   /* loop over colors  */
349194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
359194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
369194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
379194eea9SBarry Smith       x = fd->columnsforrow[i][j];
38b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
399194eea9SBarry Smith     }
409194eea9SBarry Smith   }
419194eea9SBarry Smith   PetscFunctionReturn(0);
429194eea9SBarry Smith }
439194eea9SBarry Smith 
444a2ae208SSatish Balay #undef __FUNCT__
454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
466849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
47005c665bSBarry Smith {
48dfbe8321SBarry Smith   PetscErrorCode ierr;
49005c665bSBarry Smith   PetscTruth     isnull;
50b0a32e0cSBarry Smith   PetscDraw      draw;
5154d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
52005c665bSBarry Smith 
533a40ed3dSBarry Smith   PetscFunctionBegin;
54b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
55b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
569194eea9SBarry Smith 
579194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
58005c665bSBarry Smith 
59005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
60005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
61b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
62b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
639194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
643a40ed3dSBarry Smith   PetscFunctionReturn(0);
65005c665bSBarry Smith }
66005c665bSBarry Smith 
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
69bbf0e169SBarry Smith /*@C
70639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
71bbf0e169SBarry Smith 
724c49b128SBarry Smith    Collective on MatFDColoring
73fee21e36SBarry Smith 
74ef5ee4d1SLois Curfman McInnes    Input  Parameters:
75ef5ee4d1SLois Curfman McInnes +  c - the coloring context
76ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
77ef5ee4d1SLois Curfman McInnes 
7815091d37SBarry Smith    Level: intermediate
7915091d37SBarry Smith 
80b4fc646aSLois Curfman McInnes    Notes:
81b4fc646aSLois Curfman McInnes    The available visualization contexts include
82b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
83b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
84ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
85ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
86ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
87b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
88bbf0e169SBarry Smith 
899194eea9SBarry Smith    Notes:
909194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
91b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
929194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
939194eea9SBarry Smith 
94639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
95005c665bSBarry Smith 
96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
97bbf0e169SBarry Smith @*/
98be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer)
99bbf0e169SBarry Smith {
1006849ba73SBarry Smith   PetscErrorCode    ierr;
10138baddfdSBarry Smith   PetscInt          i,j;
10232077d6dSBarry Smith   PetscTruth        isdraw,iascii;
103f3ef73ceSBarry Smith   PetscViewerFormat format;
104bbf0e169SBarry Smith 
1053a40ed3dSBarry Smith   PetscFunctionBegin;
1064482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
107b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm);
1084482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
109c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
110bbf0e169SBarry Smith 
111fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
11232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1130f5bd95cSBarry Smith   if (isdraw) {
114b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11532077d6dSBarry Smith   } else if (iascii) {
116b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
117b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
118b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
11977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
120ae09f205SBarry Smith 
121b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
122fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
123b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12577431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
126b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12777431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
128639f9d9dSBarry Smith         }
12977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13177431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
132b4fc646aSLois Curfman McInnes         }
133bbf0e169SBarry Smith       }
134bbf0e169SBarry Smith     }
135b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1365cd90555SBarry Smith   } else {
13779a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
138bbf0e169SBarry Smith   }
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
140639f9d9dSBarry Smith }
141639f9d9dSBarry Smith 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
144639f9d9dSBarry Smith /*@
145b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
146b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
147639f9d9dSBarry Smith 
148ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
149ef5ee4d1SLois Curfman McInnes 
150ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
151ef5ee4d1SLois Curfman McInnes .vb
15265f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
153f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
154f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
155ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
156ef5ee4d1SLois Curfman McInnes .ve
157639f9d9dSBarry Smith 
158639f9d9dSBarry Smith    Input Parameters:
159ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
160639f9d9dSBarry Smith .  error_rel - relative error
161f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
162fee21e36SBarry Smith 
16315091d37SBarry Smith    Level: advanced
16415091d37SBarry Smith 
165b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
166b4fc646aSLois Curfman McInnes 
167b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
168639f9d9dSBarry Smith @*/
169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
170639f9d9dSBarry Smith {
1713a40ed3dSBarry Smith   PetscFunctionBegin;
1724482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
173639f9d9dSBarry Smith 
174639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
175639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177639f9d9dSBarry Smith }
178639f9d9dSBarry Smith 
1794a2ae208SSatish Balay #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency"
181005c665bSBarry Smith /*@
182e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
183e0907662SLois Curfman McInnes    matrices.
184005c665bSBarry Smith 
185fee21e36SBarry Smith    Collective on MatFDColoring
186fee21e36SBarry Smith 
187ef5ee4d1SLois Curfman McInnes    Input Parameters:
188ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
189ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
190ef5ee4d1SLois Curfman McInnes 
19115091d37SBarry Smith    Options Database Keys:
19215091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
19315091d37SBarry Smith 
19415091d37SBarry Smith    Level: advanced
19515091d37SBarry Smith 
196e0907662SLois Curfman McInnes    Notes:
197e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
198e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
199e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
200e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
201e0907662SLois Curfman McInnes 
202b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
203ef5ee4d1SLois Curfman McInnes 
20440964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
205005c665bSBarry Smith @*/
206be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq)
207005c665bSBarry Smith {
2083a40ed3dSBarry Smith   PetscFunctionBegin;
2094482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
210005c665bSBarry Smith 
211005c665bSBarry Smith   matfd->freq = freq;
2123a40ed3dSBarry Smith   PetscFunctionReturn(0);
213005c665bSBarry Smith }
214005c665bSBarry Smith 
2154a2ae208SSatish Balay #undef __FUNCT__
2164a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency"
217ff0cfa39SBarry Smith /*@
218ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
219ff0cfa39SBarry Smith    matrices.
220ff0cfa39SBarry Smith 
221ef5ee4d1SLois Curfman McInnes    Not Collective
222ef5ee4d1SLois Curfman McInnes 
223ff0cfa39SBarry Smith    Input Parameters:
224ff0cfa39SBarry Smith .  coloring - the coloring context
225ff0cfa39SBarry Smith 
226ff0cfa39SBarry Smith    Output Parameters:
227ff0cfa39SBarry Smith .  freq - frequency (default is 1)
228ff0cfa39SBarry Smith 
22915091d37SBarry Smith    Options Database Keys:
23015091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
23115091d37SBarry Smith 
23215091d37SBarry Smith    Level: advanced
23315091d37SBarry Smith 
234ff0cfa39SBarry Smith    Notes:
235ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
236ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
237ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
238ff0cfa39SBarry Smith    <freq> nonlinear iterations.
239ff0cfa39SBarry Smith 
240ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
241ef5ee4d1SLois Curfman McInnes 
242ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
243ff0cfa39SBarry Smith @*/
244be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq)
245ff0cfa39SBarry Smith {
2463a40ed3dSBarry Smith   PetscFunctionBegin;
2474482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
248ff0cfa39SBarry Smith   *freq = matfd->freq;
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
250ff0cfa39SBarry Smith }
251ff0cfa39SBarry Smith 
2524a2ae208SSatish Balay #undef __FUNCT__
2534a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
254d64ed03dSBarry Smith /*@C
255005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
256005c665bSBarry Smith 
257fee21e36SBarry Smith    Collective on MatFDColoring
258fee21e36SBarry Smith 
259ef5ee4d1SLois Curfman McInnes    Input Parameters:
260ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
261ef5ee4d1SLois Curfman McInnes .  f - the function
262ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
263ef5ee4d1SLois Curfman McInnes 
26415091d37SBarry Smith    Level: intermediate
26515091d37SBarry Smith 
266f881d145SBarry Smith    Notes:
267f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
268f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
269f881d145SBarry Smith   with the TS solvers.
270f881d145SBarry Smith 
271b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
272005c665bSBarry Smith @*/
273be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
274005c665bSBarry Smith {
2753a40ed3dSBarry Smith   PetscFunctionBegin;
2764482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
277005c665bSBarry Smith   matfd->f    = f;
278005c665bSBarry Smith   matfd->fctx = fctx;
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
280005c665bSBarry Smith }
281005c665bSBarry Smith 
2824a2ae208SSatish Balay #undef __FUNCT__
2834a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
284639f9d9dSBarry Smith /*@
285b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
286639f9d9dSBarry Smith    the options database.
287639f9d9dSBarry Smith 
288fee21e36SBarry Smith    Collective on MatFDColoring
289fee21e36SBarry Smith 
29065f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
291ef5ee4d1SLois Curfman McInnes .vb
29265f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
293f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
294f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
295ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
296ef5ee4d1SLois Curfman McInnes .ve
297ef5ee4d1SLois Curfman McInnes 
298ef5ee4d1SLois Curfman McInnes    Input Parameter:
299ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
300ef5ee4d1SLois Curfman McInnes 
301b4fc646aSLois Curfman McInnes    Options Database Keys:
302d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
303ef5ee4d1SLois Curfman McInnes            of relative error in the function)
304f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
305ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
306*efb30889SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
307ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
308ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
309ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
310639f9d9dSBarry Smith 
31115091d37SBarry Smith     Level: intermediate
31215091d37SBarry Smith 
313b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
314d1c39f55SBarry Smith 
315d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
316d1c39f55SBarry Smith 
317639f9d9dSBarry Smith @*/
318be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
319639f9d9dSBarry Smith {
320dfbe8321SBarry Smith   PetscErrorCode ierr;
321*efb30889SBarry Smith   PetscTruth     flg;
322*efb30889SBarry Smith   char           value[3];
3233a40ed3dSBarry Smith 
3243a40ed3dSBarry Smith   PetscFunctionBegin;
3254482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
326639f9d9dSBarry Smith 
327b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
32887828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
32987828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
330b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
331*efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
332*efb30889SBarry Smith     if (flg) {
333*efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
334*efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
335*efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
336*efb30889SBarry Smith     }
337186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
338b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
339b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
340b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
341b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
343005c665bSBarry Smith }
344005c665bSBarry Smith 
3454a2ae208SSatish Balay #undef __FUNCT__
3464a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
347dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
348005c665bSBarry Smith {
349dfbe8321SBarry Smith   PetscErrorCode ierr;
350f1af5d2fSBarry Smith   PetscTruth     flg;
351005c665bSBarry Smith 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
353b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
354005c665bSBarry Smith   if (flg) {
355b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
356005c665bSBarry Smith   }
357b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
358ae09f205SBarry Smith   if (flg) {
359fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
360b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
361b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
362ae09f205SBarry Smith   }
363b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
364005c665bSBarry Smith   if (flg) {
365b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
366b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
367005c665bSBarry Smith   }
3683a40ed3dSBarry Smith   PetscFunctionReturn(0);
369bbf0e169SBarry Smith }
370bbf0e169SBarry Smith 
3714a2ae208SSatish Balay #undef __FUNCT__
3724a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
37305869f15SSatish Balay /*@
374639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
375639f9d9dSBarry Smith    computation of Jacobians.
376bbf0e169SBarry Smith 
377ef5ee4d1SLois Curfman McInnes    Collective on Mat
378ef5ee4d1SLois Curfman McInnes 
379639f9d9dSBarry Smith    Input Parameters:
380ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
381ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
382bbf0e169SBarry Smith 
383bbf0e169SBarry Smith     Output Parameter:
384639f9d9dSBarry Smith .   color - the new coloring context
385bbf0e169SBarry Smith 
38615091d37SBarry Smith     Level: intermediate
38715091d37SBarry Smith 
3886831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
389d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
390d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
391d1c39f55SBarry Smith           MatFDColoringSetParameters()
392bbf0e169SBarry Smith @*/
393be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
394bbf0e169SBarry Smith {
395639f9d9dSBarry Smith   MatFDColoring  c;
396639f9d9dSBarry Smith   MPI_Comm       comm;
397dfbe8321SBarry Smith   PetscErrorCode ierr;
39838baddfdSBarry Smith   PetscInt       M,N;
399639f9d9dSBarry Smith 
4003a40ed3dSBarry Smith   PetscFunctionBegin;
401d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
402639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
40329bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
404639f9d9dSBarry Smith 
405f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
40652e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
407639f9d9dSBarry Smith 
408f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
409f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
410639f9d9dSBarry Smith   } else {
41129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
412639f9d9dSBarry Smith   }
413639f9d9dSBarry Smith 
41477d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
41577d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
416005c665bSBarry Smith   c->freq              = 1;
41740964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
41840964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
41949b058dcSBarry Smith   c->currentcolor      = -1;
420*efb30889SBarry Smith   c->htype             = "wp";
421639f9d9dSBarry Smith 
422639f9d9dSBarry Smith   *color = c;
423d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4243a40ed3dSBarry Smith   PetscFunctionReturn(0);
425bbf0e169SBarry Smith }
426bbf0e169SBarry Smith 
4274a2ae208SSatish Balay #undef __FUNCT__
4284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
42905869f15SSatish Balay /*@
430639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
431639f9d9dSBarry Smith     via MatFDColoringCreate().
432bbf0e169SBarry Smith 
433ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
434ef5ee4d1SLois Curfman McInnes 
435b4fc646aSLois Curfman McInnes     Input Parameter:
436639f9d9dSBarry Smith .   c - coloring context
437bbf0e169SBarry Smith 
43815091d37SBarry Smith     Level: intermediate
43915091d37SBarry Smith 
440639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
441bbf0e169SBarry Smith @*/
442be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
443bbf0e169SBarry Smith {
4446849ba73SBarry Smith   PetscErrorCode ierr;
44538baddfdSBarry Smith   PetscInt       i;
446bbf0e169SBarry Smith 
4473a40ed3dSBarry Smith   PetscFunctionBegin;
4483a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
449d4bb536fSBarry Smith 
450639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
451606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
452606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
453606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
45430b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
455bbf0e169SBarry Smith   }
456606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
457606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
458606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
459606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
460606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
46130b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4624f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
463005c665bSBarry Smith   if (c->w1) {
464005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
465005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
466005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
467005c665bSBarry Smith   }
468d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4693a40ed3dSBarry Smith   PetscFunctionReturn(0);
470bbf0e169SBarry Smith }
47143a90d84SBarry Smith 
4724a2ae208SSatish Balay #undef __FUNCT__
47349b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
47449b058dcSBarry Smith /*@C
47549b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
47649b058dcSBarry Smith       that are currently being perturbed.
47749b058dcSBarry Smith 
47849b058dcSBarry Smith     Not Collective
47949b058dcSBarry Smith 
48049b058dcSBarry Smith     Input Parameters:
48149b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
48249b058dcSBarry Smith 
48349b058dcSBarry Smith     Output Parameters:
48449b058dcSBarry Smith +   n - the number of local columns being perturbed
48549b058dcSBarry Smith -   cols - the column indices, in global numbering
48649b058dcSBarry Smith 
48749b058dcSBarry Smith    Level: intermediate
48849b058dcSBarry Smith 
48949b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
49049b058dcSBarry Smith 
49149b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
49249b058dcSBarry Smith @*/
493be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
49449b058dcSBarry Smith {
49549b058dcSBarry Smith   PetscFunctionBegin;
49649b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
49749b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
49849b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
49949b058dcSBarry Smith   } else {
50049b058dcSBarry Smith     *n = 0;
50149b058dcSBarry Smith   }
50249b058dcSBarry Smith   PetscFunctionReturn(0);
50349b058dcSBarry Smith }
50449b058dcSBarry Smith 
50549b058dcSBarry Smith #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
50743a90d84SBarry Smith /*@
508e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
509e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
51043a90d84SBarry Smith 
511fee21e36SBarry Smith     Collective on MatFDColoring
512fee21e36SBarry Smith 
513ef5ee4d1SLois Curfman McInnes     Input Parameters:
514ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
515ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
516ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
517ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
518ef5ee4d1SLois Curfman McInnes 
5198bba8e72SBarry Smith     Options Database Keys:
520b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
521*efb30889SBarry Smith .    -mat_fd_type - "wp" or "ds"  (see MATSNESMF_WP or MATSNESMF_DS)
522b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
523b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
524b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5258bba8e72SBarry Smith 
52615091d37SBarry Smith     Level: intermediate
52715091d37SBarry Smith 
52843a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
52943a90d84SBarry Smith 
53043a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
53143a90d84SBarry Smith @*/
532be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
53343a90d84SBarry Smith {
5346849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5356849ba73SBarry Smith   PetscErrorCode ierr;
53638baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
537*efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
538ea709b57SSatish Balay   PetscScalar    *vscale_array;
539*efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
540005c665bSBarry Smith   Vec            w1,w2,w3;
541005c665bSBarry Smith   void           *fctx = coloring->fctx;
542f1af5d2fSBarry Smith   PetscTruth     flg;
543005c665bSBarry Smith 
544a2e34c3dSBarry Smith 
5453a40ed3dSBarry Smith   PetscFunctionBegin;
5464482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
5474482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
5484482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
549e0907662SLois Curfman McInnes 
55040964ad5SBarry Smith   if (coloring->usersetsrecompute) {
55140964ad5SBarry Smith     if (!coloring->recompute) {
55240964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
55363ba0a88SBarry Smith       ierr = PetscLogInfo((J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"));CHKERRQ(ierr);
55440964ad5SBarry Smith       PetscFunctionReturn(0);
55540964ad5SBarry Smith     } else {
55640964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
55740964ad5SBarry Smith     }
55840964ad5SBarry Smith   }
55940964ad5SBarry Smith 
560d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
56176d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
56276d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
56376d8deaeSBarry Smith   } else {
56476d8deaeSBarry Smith 
565005c665bSBarry Smith     if (!coloring->w1) {
566005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
56752e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
568005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
56952e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
570005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
57152e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
572005c665bSBarry Smith     }
573005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
57443a90d84SBarry Smith 
5754c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
576b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
577f1af5d2fSBarry Smith     if (flg) {
57863ba0a88SBarry Smith       ierr = PetscLogInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
579e0907662SLois Curfman McInnes     } else {
5800b9b6f31SBarry Smith       PetscTruth assembled;
5810b9b6f31SBarry Smith       ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5820b9b6f31SBarry Smith       if (assembled) {
58343a90d84SBarry Smith 	ierr = MatZeroEntries(J);CHKERRQ(ierr);
584e0907662SLois Curfman McInnes       }
5850b9b6f31SBarry Smith     }
58643a90d84SBarry Smith 
58743a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
58843a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
589be2fbe1fSBarry Smith 
5903a7fca6bSBarry Smith     /*
5913a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5923a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5933a7fca6bSBarry Smith     */
5943a7fca6bSBarry Smith     if (coloring->F) {
5953a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5963a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5973a7fca6bSBarry Smith       if (m1 != m2) {
5983a7fca6bSBarry Smith 	coloring->F = 0;
5993a7fca6bSBarry Smith       }
6003a7fca6bSBarry Smith     }
6013a7fca6bSBarry Smith 
602be2fbe1fSBarry Smith     if (coloring->F) {
603be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
604be2fbe1fSBarry Smith       coloring->F = 0;
605be2fbe1fSBarry Smith     } else {
60666f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
60743a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
60866f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
609be2fbe1fSBarry Smith     }
61043a90d84SBarry Smith 
611*efb30889SBarry Smith     if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
612*efb30889SBarry Smith       ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
613*efb30889SBarry Smith     }
614*efb30889SBarry Smith 
61543a90d84SBarry Smith     /*
6169782cf98SBarry Smith        Compute all the scale factors and share with other processors
6179782cf98SBarry Smith     */
6189782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
6194f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
6209782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
6219782cf98SBarry Smith       /*
6229782cf98SBarry Smith 	Loop over each column associated with color adding the
6239782cf98SBarry Smith 	perturbation to the vector w3.
6249782cf98SBarry Smith       */
6259782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
6269782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
627*efb30889SBarry Smith         if (coloring->htype[0] == 'w') {
628*efb30889SBarry Smith           dx = 1.0 + unorm;
629*efb30889SBarry Smith         } else {
6309782cf98SBarry Smith   	  dx  = xx[col];
631*efb30889SBarry Smith         }
6325b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
6339782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6349782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
6359782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
6369782cf98SBarry Smith #else
6379782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6389782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6399782cf98SBarry Smith #endif
6409782cf98SBarry Smith 	dx                *= epsilon;
64130b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
6429782cf98SBarry Smith       }
6439782cf98SBarry Smith     }
644*efb30889SBarry Smith     vscale_array = vscale_array + start;
645*efb30889SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
64630b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64730b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6489782cf98SBarry Smith 
649ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
650ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
651b0a32e0cSBarry Smith 
65230b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
65330b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
65430b34957SBarry Smith 
65530b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6569782cf98SBarry Smith     /*
65743a90d84SBarry Smith       Loop over each color
65843a90d84SBarry Smith     */
65943a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
66049b058dcSBarry Smith       coloring->currentcolor = k;
66143a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
66242460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
66343a90d84SBarry Smith       /*
66443a90d84SBarry Smith 	Loop over each column associated with color adding the
66543a90d84SBarry Smith 	perturbation to the vector w3.
66643a90d84SBarry Smith       */
66743a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
66843a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
669*efb30889SBarry Smith         if (coloring->htype[0] == 'w') {
670*efb30889SBarry Smith           dx = 1.0 + unorm;
671*efb30889SBarry Smith         } else {
67242460c72SBarry Smith   	  dx  = xx[col];
673*efb30889SBarry Smith         }
6745b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
675aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
676ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
677ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
67843a90d84SBarry Smith #else
679329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
680329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
68143a90d84SBarry Smith #endif
68243a90d84SBarry Smith 	dx            *= epsilon;
6831302d50aSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
68442460c72SBarry Smith 	w3_array[col] += dx;
68543a90d84SBarry Smith       }
68642460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6873b28642cSBarry Smith 
68843a90d84SBarry Smith       /*
689e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
69043a90d84SBarry Smith       */
691a2e34c3dSBarry Smith 
69266f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
69343a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
69466f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
695*efb30889SBarry Smith       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6969782cf98SBarry Smith 
69743a90d84SBarry Smith       /*
698e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
69943a90d84SBarry Smith       */
7003b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
70143a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
70243a90d84SBarry Smith 	row    = coloring->rows[k][l];
70343a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
7045904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
70543a90d84SBarry Smith 	srow   = row + start;
70643a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
70743a90d84SBarry Smith       }
7083b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
70943a90d84SBarry Smith     }
71049b058dcSBarry Smith     coloring->currentcolor = -1;
71130b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
71242460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
71343a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71443a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71576d8deaeSBarry Smith   }
716d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
717a2e34c3dSBarry Smith 
718b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
719a2e34c3dSBarry Smith   if (flg) {
720a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
721a2e34c3dSBarry Smith   }
722b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
723b9722508SLois Curfman McInnes 
7243a40ed3dSBarry Smith   PetscFunctionReturn(0);
72543a90d84SBarry Smith }
726840b8ebdSBarry Smith 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
729840b8ebdSBarry Smith /*@
730840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
731840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
732840b8ebdSBarry Smith 
733fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
734fee21e36SBarry Smith 
735ef5ee4d1SLois Curfman McInnes     Input Parameters:
7363b28642cSBarry Smith +   mat - location to store Jacobian
737ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
738ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
739ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
740ef5ee4d1SLois Curfman McInnes 
741840b8ebdSBarry Smith    Options Database Keys:
742ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
743840b8ebdSBarry Smith 
74415091d37SBarry Smith    Level: intermediate
74515091d37SBarry Smith 
746840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
747840b8ebdSBarry Smith 
748840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
749840b8ebdSBarry Smith @*/
750be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
751840b8ebdSBarry Smith {
7526849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7536849ba73SBarry Smith   PetscErrorCode ierr;
75438baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
755*efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
756ea709b57SSatish Balay   PetscScalar    *vscale_array;
757329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
758840b8ebdSBarry Smith   Vec            w1,w2,w3;
759840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
760f1af5d2fSBarry Smith   PetscTruth     flg;
761840b8ebdSBarry Smith 
7623a40ed3dSBarry Smith   PetscFunctionBegin;
7634482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
7644482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
7654482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
766840b8ebdSBarry Smith 
767d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
768840b8ebdSBarry Smith   if (!coloring->w1) {
769840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
77052e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
771840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
77252e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
773840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
77452e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
775840b8ebdSBarry Smith   }
776840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
777840b8ebdSBarry Smith 
7785904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
779b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
780f1af5d2fSBarry Smith   if (flg) {
78163ba0a88SBarry Smith     ierr = PetscLogInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
782840b8ebdSBarry Smith   } else {
7830b9b6f31SBarry Smith     PetscTruth assembled;
7840b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7850b9b6f31SBarry Smith     if (assembled) {
786840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
787840b8ebdSBarry Smith     }
7880b9b6f31SBarry Smith   }
789840b8ebdSBarry Smith 
790840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
791840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
79266f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
793840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
79466f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
795840b8ebdSBarry Smith 
796840b8ebdSBarry Smith   /*
7975904e1b1SBarry Smith       Compute all the scale factors and share with other processors
798840b8ebdSBarry Smith   */
7995904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
8005904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
801840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
802840b8ebdSBarry Smith     /*
803840b8ebdSBarry Smith        Loop over each column associated with color adding the
804840b8ebdSBarry Smith        perturbation to the vector w3.
805840b8ebdSBarry Smith     */
806840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
807840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8085904e1b1SBarry Smith       dx  = xx[col];
8095b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
810aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
811840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
812840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
813840b8ebdSBarry Smith #else
814329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
815329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
816840b8ebdSBarry Smith #endif
817840b8ebdSBarry Smith       dx                *= epsilon;
8185904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
819840b8ebdSBarry Smith     }
8205904e1b1SBarry Smith   }
8215904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8225904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8235904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8245904e1b1SBarry Smith 
8255904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8265904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8275904e1b1SBarry Smith 
8285904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8295904e1b1SBarry Smith   /*
8305904e1b1SBarry Smith       Loop over each color
8315904e1b1SBarry Smith   */
8325904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8335904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8345904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8355904e1b1SBarry Smith     /*
8365904e1b1SBarry Smith        Loop over each column associated with color adding the
8375904e1b1SBarry Smith        perturbation to the vector w3.
8385904e1b1SBarry Smith     */
8395904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8405904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8415904e1b1SBarry Smith       dx  = xx[col];
8425b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8435904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8445904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8455904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8465904e1b1SBarry Smith #else
8475904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8485904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8495904e1b1SBarry Smith #endif
8505904e1b1SBarry Smith       dx            *= epsilon;
8515904e1b1SBarry Smith       w3_array[col] += dx;
8525904e1b1SBarry Smith     }
8535904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8545904e1b1SBarry Smith 
855840b8ebdSBarry Smith     /*
856840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
857840b8ebdSBarry Smith     */
85866f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
859840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
86066f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
861*efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
8625904e1b1SBarry Smith 
863840b8ebdSBarry Smith     /*
864840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
865840b8ebdSBarry Smith     */
8663b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
867840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
868840b8ebdSBarry Smith       row    = coloring->rows[k][l];
869840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8705904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
871840b8ebdSBarry Smith       srow   = row + start;
872840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
873840b8ebdSBarry Smith     }
8743b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
875840b8ebdSBarry Smith   }
8765904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8775904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
878840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
879840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
880d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8813a40ed3dSBarry Smith   PetscFunctionReturn(0);
882840b8ebdSBarry Smith }
8833b28642cSBarry Smith 
8843b28642cSBarry Smith 
8854a2ae208SSatish Balay #undef __FUNCT__
8864a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
88740964ad5SBarry Smith /*@C
88840964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
88940964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
89040964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
89140964ad5SBarry Smith      called again.
89240964ad5SBarry Smith 
89340964ad5SBarry Smith    Collective on MatFDColoring
89440964ad5SBarry Smith 
89540964ad5SBarry Smith    Input  Parameters:
89640964ad5SBarry Smith .  c - the coloring context
89740964ad5SBarry Smith 
89840964ad5SBarry Smith    Level: intermediate
89940964ad5SBarry Smith 
90040964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
90140964ad5SBarry Smith 
90240964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
90340964ad5SBarry Smith 
90440964ad5SBarry Smith .keywords: Mat, finite differences, coloring
90540964ad5SBarry Smith @*/
906be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
90740964ad5SBarry Smith {
90840964ad5SBarry Smith   PetscFunctionBegin;
9094482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
91040964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
91140964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
91240964ad5SBarry Smith   PetscFunctionReturn(0);
91340964ad5SBarry Smith }
91440964ad5SBarry Smith 
9153b28642cSBarry Smith 
916