xref: /petsc/src/mat/matfd/fdmatrix.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris 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 */
11*be1d678aSKris Buschelman PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE = 0;
128ba1e511SMatthew Knepley 
134a2ae208SSatish Balay #undef __FUNCT__
143a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
15*be1d678aSKris 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 @*/
98*be1d678aSKris 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 @*/
169*be1d678aSKris 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 @*/
206*be1d678aSKris 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 @*/
244*be1d678aSKris 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 @*/
273*be1d678aSKris 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
306ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
307ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
308ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
309639f9d9dSBarry Smith 
31015091d37SBarry Smith     Level: intermediate
31115091d37SBarry Smith 
312b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
313d1c39f55SBarry Smith 
314d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
315d1c39f55SBarry Smith 
316639f9d9dSBarry Smith @*/
317*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
318639f9d9dSBarry Smith {
319dfbe8321SBarry Smith   PetscErrorCode ierr;
3203a40ed3dSBarry Smith 
3213a40ed3dSBarry Smith   PetscFunctionBegin;
3224482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
323639f9d9dSBarry Smith 
324b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
32587828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
32687828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
327b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
328186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
329b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
330b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
331b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
332b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
334005c665bSBarry Smith }
335005c665bSBarry Smith 
3364a2ae208SSatish Balay #undef __FUNCT__
3374a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
338dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
339005c665bSBarry Smith {
340dfbe8321SBarry Smith   PetscErrorCode ierr;
341f1af5d2fSBarry Smith   PetscTruth     flg;
342005c665bSBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
344b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
345005c665bSBarry Smith   if (flg) {
346b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
347005c665bSBarry Smith   }
348b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
349ae09f205SBarry Smith   if (flg) {
350fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
351b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
352b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
353ae09f205SBarry Smith   }
354b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
355005c665bSBarry Smith   if (flg) {
356b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
357b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
358005c665bSBarry Smith   }
3593a40ed3dSBarry Smith   PetscFunctionReturn(0);
360bbf0e169SBarry Smith }
361bbf0e169SBarry Smith 
3624a2ae208SSatish Balay #undef __FUNCT__
3634a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
364bbf0e169SBarry Smith /*@C
365639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
366639f9d9dSBarry Smith    computation of Jacobians.
367bbf0e169SBarry Smith 
368ef5ee4d1SLois Curfman McInnes    Collective on Mat
369ef5ee4d1SLois Curfman McInnes 
370639f9d9dSBarry Smith    Input Parameters:
371ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
372ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
373bbf0e169SBarry Smith 
374bbf0e169SBarry Smith     Output Parameter:
375639f9d9dSBarry Smith .   color - the new coloring context
376bbf0e169SBarry Smith 
37715091d37SBarry Smith     Level: intermediate
37815091d37SBarry Smith 
3796831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
380d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
381d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
382d1c39f55SBarry Smith           MatFDColoringSetParameters()
383bbf0e169SBarry Smith @*/
384*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
385bbf0e169SBarry Smith {
386639f9d9dSBarry Smith   MatFDColoring  c;
387639f9d9dSBarry Smith   MPI_Comm       comm;
388dfbe8321SBarry Smith   PetscErrorCode ierr;
38938baddfdSBarry Smith   PetscInt       M,N;
390639f9d9dSBarry Smith 
3913a40ed3dSBarry Smith   PetscFunctionBegin;
392d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
393639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
39429bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
395639f9d9dSBarry Smith 
396f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
39752e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
398639f9d9dSBarry Smith 
399f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
400f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
401639f9d9dSBarry Smith   } else {
40229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
403639f9d9dSBarry Smith   }
404639f9d9dSBarry Smith 
40577d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
40677d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
407005c665bSBarry Smith   c->freq              = 1;
40840964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
40940964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
41049b058dcSBarry Smith   c->currentcolor      = -1;
411639f9d9dSBarry Smith 
412639f9d9dSBarry Smith   *color = c;
413d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
415bbf0e169SBarry Smith }
416bbf0e169SBarry Smith 
4174a2ae208SSatish Balay #undef __FUNCT__
4184a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
419bbf0e169SBarry Smith /*@C
420639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
421639f9d9dSBarry Smith     via MatFDColoringCreate().
422bbf0e169SBarry Smith 
423ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
424ef5ee4d1SLois Curfman McInnes 
425b4fc646aSLois Curfman McInnes     Input Parameter:
426639f9d9dSBarry Smith .   c - coloring context
427bbf0e169SBarry Smith 
42815091d37SBarry Smith     Level: intermediate
42915091d37SBarry Smith 
430639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
431bbf0e169SBarry Smith @*/
432*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
433bbf0e169SBarry Smith {
4346849ba73SBarry Smith   PetscErrorCode ierr;
43538baddfdSBarry Smith   PetscInt       i;
436bbf0e169SBarry Smith 
4373a40ed3dSBarry Smith   PetscFunctionBegin;
4383a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
439d4bb536fSBarry Smith 
440639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
441606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
442606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
443606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
44430b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
445bbf0e169SBarry Smith   }
446606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
447606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
448606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
449606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
450606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
45130b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4524f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
453005c665bSBarry Smith   if (c->w1) {
454005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
455005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
456005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
457005c665bSBarry Smith   }
458d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4593a40ed3dSBarry Smith   PetscFunctionReturn(0);
460bbf0e169SBarry Smith }
46143a90d84SBarry Smith 
4624a2ae208SSatish Balay #undef __FUNCT__
46349b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
46449b058dcSBarry Smith /*@C
46549b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
46649b058dcSBarry Smith       that are currently being perturbed.
46749b058dcSBarry Smith 
46849b058dcSBarry Smith     Not Collective
46949b058dcSBarry Smith 
47049b058dcSBarry Smith     Input Parameters:
47149b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
47249b058dcSBarry Smith 
47349b058dcSBarry Smith     Output Parameters:
47449b058dcSBarry Smith +   n - the number of local columns being perturbed
47549b058dcSBarry Smith -   cols - the column indices, in global numbering
47649b058dcSBarry Smith 
47749b058dcSBarry Smith    Level: intermediate
47849b058dcSBarry Smith 
47949b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
48049b058dcSBarry Smith 
48149b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
48249b058dcSBarry Smith @*/
483*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
48449b058dcSBarry Smith {
48549b058dcSBarry Smith   PetscFunctionBegin;
48649b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
48749b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
48849b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
48949b058dcSBarry Smith   } else {
49049b058dcSBarry Smith     *n = 0;
49149b058dcSBarry Smith   }
49249b058dcSBarry Smith   PetscFunctionReturn(0);
49349b058dcSBarry Smith }
49449b058dcSBarry Smith 
49549b058dcSBarry Smith #undef __FUNCT__
4964a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
49743a90d84SBarry Smith /*@
498e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
499e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
50043a90d84SBarry Smith 
501fee21e36SBarry Smith     Collective on MatFDColoring
502fee21e36SBarry Smith 
503ef5ee4d1SLois Curfman McInnes     Input Parameters:
504ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
505ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
506ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
507ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
508ef5ee4d1SLois Curfman McInnes 
5098bba8e72SBarry Smith     Options Database Keys:
510b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
511b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
512b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
513b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5148bba8e72SBarry Smith 
51515091d37SBarry Smith     Level: intermediate
51615091d37SBarry Smith 
51743a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
51843a90d84SBarry Smith 
51943a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
52043a90d84SBarry Smith @*/
521*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
52243a90d84SBarry Smith {
5236849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5246849ba73SBarry Smith   PetscErrorCode ierr;
52538baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
526ea709b57SSatish Balay   PetscScalar    dx,mone = -1.0,*y,*xx,*w3_array;
527ea709b57SSatish Balay   PetscScalar    *vscale_array;
528329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
529005c665bSBarry Smith   Vec            w1,w2,w3;
530005c665bSBarry Smith   void           *fctx = coloring->fctx;
531f1af5d2fSBarry Smith   PetscTruth     flg;
532005c665bSBarry Smith 
533a2e34c3dSBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
5354482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
5364482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
5374482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
538e0907662SLois Curfman McInnes 
53940964ad5SBarry Smith   if (coloring->usersetsrecompute) {
54040964ad5SBarry Smith     if (!coloring->recompute) {
54140964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
54263ba0a88SBarry Smith       ierr = PetscLogInfo((J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"));CHKERRQ(ierr);
54340964ad5SBarry Smith       PetscFunctionReturn(0);
54440964ad5SBarry Smith     } else {
54540964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
54640964ad5SBarry Smith     }
54740964ad5SBarry Smith   }
54840964ad5SBarry Smith 
549d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
55076d8deaeSBarry Smith   if (J->ops->fdcoloringapply) {
55176d8deaeSBarry Smith     ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
55276d8deaeSBarry Smith   } else {
55376d8deaeSBarry Smith 
554005c665bSBarry Smith     if (!coloring->w1) {
555005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
55652e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
557005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
55852e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
559005c665bSBarry Smith       ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
56052e6d16bSBarry Smith       ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
561005c665bSBarry Smith     }
562005c665bSBarry Smith     w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
56343a90d84SBarry Smith 
5644c49b128SBarry Smith     ierr = MatSetUnfactored(J);CHKERRQ(ierr);
565b0a32e0cSBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
566f1af5d2fSBarry Smith     if (flg) {
56763ba0a88SBarry Smith       ierr = PetscLogInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
568e0907662SLois Curfman McInnes     } else {
5690b9b6f31SBarry Smith       PetscTruth assembled;
5700b9b6f31SBarry Smith       ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5710b9b6f31SBarry Smith       if (assembled) {
57243a90d84SBarry Smith 	ierr = MatZeroEntries(J);CHKERRQ(ierr);
573e0907662SLois Curfman McInnes       }
5740b9b6f31SBarry Smith     }
57543a90d84SBarry Smith 
57643a90d84SBarry Smith     ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
57743a90d84SBarry Smith     ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
578be2fbe1fSBarry Smith 
5793a7fca6bSBarry Smith     /*
5803a7fca6bSBarry Smith       This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5813a7fca6bSBarry Smith       coloring->F for the coarser grids from the finest
5823a7fca6bSBarry Smith     */
5833a7fca6bSBarry Smith     if (coloring->F) {
5843a7fca6bSBarry Smith       ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5853a7fca6bSBarry Smith       ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5863a7fca6bSBarry Smith       if (m1 != m2) {
5873a7fca6bSBarry Smith 	coloring->F = 0;
5883a7fca6bSBarry Smith       }
5893a7fca6bSBarry Smith     }
5903a7fca6bSBarry Smith 
591be2fbe1fSBarry Smith     if (coloring->F) {
592be2fbe1fSBarry Smith       w1          = coloring->F; /* use already computed value of function */
593be2fbe1fSBarry Smith       coloring->F = 0;
594be2fbe1fSBarry Smith     } else {
59566f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
59643a90d84SBarry Smith       ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
59766f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
598be2fbe1fSBarry Smith     }
59943a90d84SBarry Smith 
60043a90d84SBarry Smith     /*
6019782cf98SBarry Smith        Compute all the scale factors and share with other processors
6029782cf98SBarry Smith     */
6039782cf98SBarry Smith     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
6044f113abfSBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
6059782cf98SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
6069782cf98SBarry Smith       /*
6079782cf98SBarry Smith 	Loop over each column associated with color adding the
6089782cf98SBarry Smith 	perturbation to the vector w3.
6099782cf98SBarry Smith       */
6109782cf98SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
6119782cf98SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
6129782cf98SBarry Smith 	dx  = xx[col];
6135b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
6149782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6159782cf98SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
6169782cf98SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
6179782cf98SBarry Smith #else
6189782cf98SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6199782cf98SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6209782cf98SBarry Smith #endif
6219782cf98SBarry Smith 	dx                *= epsilon;
62230b34957SBarry Smith 	vscale_array[col] = 1.0/dx;
6239782cf98SBarry Smith       }
6249782cf98SBarry Smith     }
625a2e34c3dSBarry Smith     vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
62630b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62730b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6289782cf98SBarry Smith 
629ce1411ecSBarry Smith     /*  ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
630ce1411ecSBarry Smith 	ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
631b0a32e0cSBarry Smith 
63230b34957SBarry Smith     if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
63330b34957SBarry Smith     else                        vscaleforrow = coloring->columnsforrow;
63430b34957SBarry Smith 
63530b34957SBarry Smith     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6369782cf98SBarry Smith     /*
63743a90d84SBarry Smith       Loop over each color
63843a90d84SBarry Smith     */
63943a90d84SBarry Smith     for (k=0; k<coloring->ncolors; k++) {
64049b058dcSBarry Smith       coloring->currentcolor = k;
64143a90d84SBarry Smith       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
64242460c72SBarry Smith       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
64343a90d84SBarry Smith       /*
64443a90d84SBarry Smith 	Loop over each column associated with color adding the
64543a90d84SBarry Smith 	perturbation to the vector w3.
64643a90d84SBarry Smith       */
64743a90d84SBarry Smith       for (l=0; l<coloring->ncolumns[k]; l++) {
64843a90d84SBarry Smith 	col = coloring->columns[k][l];    /* column of the matrix we are probing for */
64942460c72SBarry Smith 	dx  = xx[col];
6505b8514ebSBarry Smith 	if (dx == 0.0) dx = 1.0;
651aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
652ae09f205SBarry Smith 	if (dx < umin && dx >= 0.0)      dx = umin;
653ae09f205SBarry Smith 	else if (dx < 0.0 && dx > -umin) dx = -umin;
65443a90d84SBarry Smith #else
655329f5518SBarry Smith 	if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
656329f5518SBarry Smith 	else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
65743a90d84SBarry Smith #endif
65843a90d84SBarry Smith 	dx            *= epsilon;
6591302d50aSBarry Smith 	if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
66042460c72SBarry Smith 	w3_array[col] += dx;
66143a90d84SBarry Smith       }
66242460c72SBarry Smith       w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6633b28642cSBarry Smith 
66443a90d84SBarry Smith       /*
665e0907662SLois Curfman McInnes 	Evaluate function at x1 + dx (here dx is a vector of perturbations)
66643a90d84SBarry Smith       */
667a2e34c3dSBarry Smith 
66866f9b7ceSBarry Smith       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
66943a90d84SBarry Smith       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
67066f9b7ceSBarry Smith       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
67143a90d84SBarry Smith       ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6729782cf98SBarry Smith 
67343a90d84SBarry Smith       /*
674e0907662SLois Curfman McInnes 	Loop over rows of vector, putting results into Jacobian matrix
67543a90d84SBarry Smith       */
6763b28642cSBarry Smith       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
67743a90d84SBarry Smith       for (l=0; l<coloring->nrows[k]; l++) {
67843a90d84SBarry Smith 	row    = coloring->rows[k][l];
67943a90d84SBarry Smith 	col    = coloring->columnsforrow[k][l];
6805904e1b1SBarry Smith 	y[row] *= vscale_array[vscaleforrow[k][l]];
68143a90d84SBarry Smith 	srow   = row + start;
68243a90d84SBarry Smith 	ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
68343a90d84SBarry Smith       }
6843b28642cSBarry Smith       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
68543a90d84SBarry Smith     }
68649b058dcSBarry Smith     coloring->currentcolor = -1;
68730b34957SBarry Smith     ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
68842460c72SBarry Smith     xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
68943a90d84SBarry Smith     ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69043a90d84SBarry Smith     ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69176d8deaeSBarry Smith   }
692d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
693a2e34c3dSBarry Smith 
694b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
695a2e34c3dSBarry Smith   if (flg) {
696a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
697a2e34c3dSBarry Smith   }
698b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
699b9722508SLois Curfman McInnes 
7003a40ed3dSBarry Smith   PetscFunctionReturn(0);
70143a90d84SBarry Smith }
702840b8ebdSBarry Smith 
7034a2ae208SSatish Balay #undef __FUNCT__
7044a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
705840b8ebdSBarry Smith /*@
706840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
707840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
708840b8ebdSBarry Smith 
709fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
710fee21e36SBarry Smith 
711ef5ee4d1SLois Curfman McInnes     Input Parameters:
7123b28642cSBarry Smith +   mat - location to store Jacobian
713ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
714ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
715ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
716ef5ee4d1SLois Curfman McInnes 
717840b8ebdSBarry Smith    Options Database Keys:
718ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
719840b8ebdSBarry Smith 
72015091d37SBarry Smith    Level: intermediate
72115091d37SBarry Smith 
722840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
723840b8ebdSBarry Smith 
724840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
725840b8ebdSBarry Smith @*/
726*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
727840b8ebdSBarry Smith {
7286849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7296849ba73SBarry Smith   PetscErrorCode ierr;
73038baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
731ea709b57SSatish Balay   PetscScalar    dx,mone = -1.0,*y,*xx,*w3_array;
732ea709b57SSatish Balay   PetscScalar    *vscale_array;
733329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
734840b8ebdSBarry Smith   Vec            w1,w2,w3;
735840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
736f1af5d2fSBarry Smith   PetscTruth     flg;
737840b8ebdSBarry Smith 
7383a40ed3dSBarry Smith   PetscFunctionBegin;
7394482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
7404482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
7414482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
742840b8ebdSBarry Smith 
743d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
744840b8ebdSBarry Smith   if (!coloring->w1) {
745840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
74652e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
747840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
74852e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
749840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
75052e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
751840b8ebdSBarry Smith   }
752840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
753840b8ebdSBarry Smith 
7545904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
755b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
756f1af5d2fSBarry Smith   if (flg) {
75763ba0a88SBarry Smith     ierr = PetscLogInfo((coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"));CHKERRQ(ierr);
758840b8ebdSBarry Smith   } else {
7590b9b6f31SBarry Smith     PetscTruth assembled;
7600b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
7610b9b6f31SBarry Smith     if (assembled) {
762840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
763840b8ebdSBarry Smith     }
7640b9b6f31SBarry Smith   }
765840b8ebdSBarry Smith 
766840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
767840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
76866f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
769840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
77066f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
771840b8ebdSBarry Smith 
772840b8ebdSBarry Smith   /*
7735904e1b1SBarry Smith       Compute all the scale factors and share with other processors
774840b8ebdSBarry Smith   */
7755904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
7765904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
777840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
778840b8ebdSBarry Smith     /*
779840b8ebdSBarry Smith        Loop over each column associated with color adding the
780840b8ebdSBarry Smith        perturbation to the vector w3.
781840b8ebdSBarry Smith     */
782840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
783840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7845904e1b1SBarry Smith       dx  = xx[col];
7855b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
786aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
787840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
788840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
789840b8ebdSBarry Smith #else
790329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
791329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
792840b8ebdSBarry Smith #endif
793840b8ebdSBarry Smith       dx                *= epsilon;
7945904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
795840b8ebdSBarry Smith     }
7965904e1b1SBarry Smith   }
7975904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7985904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7995904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8005904e1b1SBarry Smith 
8015904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8025904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8035904e1b1SBarry Smith 
8045904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8055904e1b1SBarry Smith   /*
8065904e1b1SBarry Smith       Loop over each color
8075904e1b1SBarry Smith   */
8085904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8095904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8105904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8115904e1b1SBarry Smith     /*
8125904e1b1SBarry Smith        Loop over each column associated with color adding the
8135904e1b1SBarry Smith        perturbation to the vector w3.
8145904e1b1SBarry Smith     */
8155904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8165904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8175904e1b1SBarry Smith       dx  = xx[col];
8185b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8195904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8205904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8215904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8225904e1b1SBarry Smith #else
8235904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8245904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8255904e1b1SBarry Smith #endif
8265904e1b1SBarry Smith       dx            *= epsilon;
8275904e1b1SBarry Smith       w3_array[col] += dx;
8285904e1b1SBarry Smith     }
8295904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8305904e1b1SBarry Smith 
831840b8ebdSBarry Smith     /*
832840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
833840b8ebdSBarry Smith     */
83466f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
835840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
83666f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
837840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
8385904e1b1SBarry Smith 
839840b8ebdSBarry Smith     /*
840840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
841840b8ebdSBarry Smith     */
8423b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
843840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
844840b8ebdSBarry Smith       row    = coloring->rows[k][l];
845840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
8465904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
847840b8ebdSBarry Smith       srow   = row + start;
848840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
849840b8ebdSBarry Smith     }
8503b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
851840b8ebdSBarry Smith   }
8525904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8535904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
854840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
855840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
858840b8ebdSBarry Smith }
8593b28642cSBarry Smith 
8603b28642cSBarry Smith 
8614a2ae208SSatish Balay #undef __FUNCT__
8624a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
86340964ad5SBarry Smith /*@C
86440964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
86540964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
86640964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
86740964ad5SBarry Smith      called again.
86840964ad5SBarry Smith 
86940964ad5SBarry Smith    Collective on MatFDColoring
87040964ad5SBarry Smith 
87140964ad5SBarry Smith    Input  Parameters:
87240964ad5SBarry Smith .  c - the coloring context
87340964ad5SBarry Smith 
87440964ad5SBarry Smith    Level: intermediate
87540964ad5SBarry Smith 
87640964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
87740964ad5SBarry Smith 
87840964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
87940964ad5SBarry Smith 
88040964ad5SBarry Smith .keywords: Mat, finite differences, coloring
88140964ad5SBarry Smith @*/
882*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
88340964ad5SBarry Smith {
88440964ad5SBarry Smith   PetscFunctionBegin;
8854482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
88640964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
88740964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
88840964ad5SBarry Smith   PetscFunctionReturn(0);
88940964ad5SBarry Smith }
89040964ad5SBarry Smith 
8913b28642cSBarry Smith 
892