xref: /petsc/src/mat/matfd/fdmatrix.c (revision aeb7e98aec3d9a965d430920b809464becf809f2)
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);
117a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
118a83599f4SBarry 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__
2534a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
2544a9d489dSBarry Smith /*@C
2554a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2564a9d489dSBarry Smith 
2574a9d489dSBarry Smith    Collective on MatFDColoring
2584a9d489dSBarry Smith 
2594a9d489dSBarry Smith    Input Parameters:
2604a9d489dSBarry Smith .  coloring - the coloring context
2614a9d489dSBarry Smith 
2624a9d489dSBarry Smith    Output Parameters:
2634a9d489dSBarry Smith +  f - the function
2644a9d489dSBarry Smith -  fctx - the optional user-defined function context
2654a9d489dSBarry Smith 
2664a9d489dSBarry Smith    Level: intermediate
2674a9d489dSBarry Smith 
2684a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
2694a9d489dSBarry Smith @*/
2704a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2714a9d489dSBarry Smith {
2724a9d489dSBarry Smith   PetscFunctionBegin;
2734a9d489dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
2744a9d489dSBarry Smith   if (f) *f = matfd->f;
2754a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2764a9d489dSBarry Smith   PetscFunctionReturn(0);
2774a9d489dSBarry Smith }
2784a9d489dSBarry Smith 
2794a9d489dSBarry Smith #undef __FUNCT__
2804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
281d64ed03dSBarry Smith /*@C
282005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
283005c665bSBarry Smith 
284fee21e36SBarry Smith    Collective on MatFDColoring
285fee21e36SBarry Smith 
286ef5ee4d1SLois Curfman McInnes    Input Parameters:
287ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
288ef5ee4d1SLois Curfman McInnes .  f - the function
289ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
290ef5ee4d1SLois Curfman McInnes 
29115091d37SBarry Smith    Level: intermediate
29215091d37SBarry Smith 
293f881d145SBarry Smith    Notes:
294f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
295f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
296f881d145SBarry Smith   with the TS solvers.
297f881d145SBarry Smith 
298b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
299005c665bSBarry Smith @*/
300be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
301005c665bSBarry Smith {
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3034482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
304005c665bSBarry Smith   matfd->f    = f;
305005c665bSBarry Smith   matfd->fctx = fctx;
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
307005c665bSBarry Smith }
308005c665bSBarry Smith 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
311639f9d9dSBarry Smith /*@
312b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
313639f9d9dSBarry Smith    the options database.
314639f9d9dSBarry Smith 
315fee21e36SBarry Smith    Collective on MatFDColoring
316fee21e36SBarry Smith 
31765f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
318ef5ee4d1SLois Curfman McInnes .vb
31965f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
320f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
321f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
322ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
323ef5ee4d1SLois Curfman McInnes .ve
324ef5ee4d1SLois Curfman McInnes 
325ef5ee4d1SLois Curfman McInnes    Input Parameter:
326ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
327ef5ee4d1SLois Curfman McInnes 
328b4fc646aSLois Curfman McInnes    Options Database Keys:
329d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
330ef5ee4d1SLois Curfman McInnes            of relative error in the function)
331f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
332ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
333efb30889SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
334ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
335ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
336ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
337639f9d9dSBarry Smith 
33815091d37SBarry Smith     Level: intermediate
33915091d37SBarry Smith 
340b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
341d1c39f55SBarry Smith 
342d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
343d1c39f55SBarry Smith 
344639f9d9dSBarry Smith @*/
345be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
346639f9d9dSBarry Smith {
347dfbe8321SBarry Smith   PetscErrorCode ierr;
348efb30889SBarry Smith   PetscTruth     flg;
349efb30889SBarry Smith   char           value[3];
350b8f8c88eSHong Zhang   PetscMPIInt    size;
351b8f8c88eSHong Zhang   const char     *isctypes[] = {"IS_COLORING_LOCAL","IS_COLORING_GHOSTED"};
352b8f8c88eSHong Zhang   PetscInt       isctype;
3533a40ed3dSBarry Smith 
3543a40ed3dSBarry Smith   PetscFunctionBegin;
3554482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
356639f9d9dSBarry Smith 
357b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
35887828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
35987828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
360b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
361b8f8c88eSHong Zhang 
362b8f8c88eSHong Zhang     ierr = MPI_Comm_size(matfd->comm,&size);CHKERRQ(ierr);
363b8f8c88eSHong Zhang     /* set default coloring type */
364b8f8c88eSHong Zhang     if (size == 1){
365b8f8c88eSHong Zhang       isctype      = 0; /* IS_COLORING_LOCAL */
366b8f8c88eSHong Zhang     } else {
367b8f8c88eSHong Zhang       isctype      = 0; /* should be 1, IS_COLORING_GHOSTED */
368b8f8c88eSHong Zhang     }
369b8f8c88eSHong Zhang     ierr = PetscOptionsEList("-mat_fd_coloring_type","Type of MatFDColoring","None",isctypes,2,isctypes[isctype],&isctype,PETSC_NULL);CHKERRQ(ierr);
370b8f8c88eSHong Zhang     matfd->ctype = (ISColoringType)isctype;
371b8f8c88eSHong Zhang 
372efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
373efb30889SBarry Smith     if (flg) {
374efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
375efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
376efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
377efb30889SBarry Smith     }
378186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
379b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
380b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
381b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
382b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3833a40ed3dSBarry Smith   PetscFunctionReturn(0);
384005c665bSBarry Smith }
385005c665bSBarry Smith 
3864a2ae208SSatish Balay #undef __FUNCT__
3874a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
388dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
389005c665bSBarry Smith {
390dfbe8321SBarry Smith   PetscErrorCode ierr;
391f1af5d2fSBarry Smith   PetscTruth     flg;
392005c665bSBarry Smith 
3933a40ed3dSBarry Smith   PetscFunctionBegin;
394b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
395005c665bSBarry Smith   if (flg) {
396b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
397005c665bSBarry Smith   }
398b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
399ae09f205SBarry Smith   if (flg) {
400fb9695e5SSatish Balay     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
401b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
402b0a32e0cSBarry Smith     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
403ae09f205SBarry Smith   }
404b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
405005c665bSBarry Smith   if (flg) {
406b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
407b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
408005c665bSBarry Smith   }
4093a40ed3dSBarry Smith   PetscFunctionReturn(0);
410bbf0e169SBarry Smith }
411bbf0e169SBarry Smith 
4124a2ae208SSatish Balay #undef __FUNCT__
4134a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
41405869f15SSatish Balay /*@
415639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
416639f9d9dSBarry Smith    computation of Jacobians.
417bbf0e169SBarry Smith 
418ef5ee4d1SLois Curfman McInnes    Collective on Mat
419ef5ee4d1SLois Curfman McInnes 
420639f9d9dSBarry Smith    Input Parameters:
421ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
422ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
423bbf0e169SBarry Smith 
424bbf0e169SBarry Smith     Output Parameter:
425639f9d9dSBarry Smith .   color - the new coloring context
426bbf0e169SBarry Smith 
42715091d37SBarry Smith     Level: intermediate
42815091d37SBarry Smith 
4296831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
430d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
431d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
432d1c39f55SBarry Smith           MatFDColoringSetParameters()
433bbf0e169SBarry Smith @*/
434be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
435bbf0e169SBarry Smith {
436639f9d9dSBarry Smith   MatFDColoring  c;
437639f9d9dSBarry Smith   MPI_Comm       comm;
438dfbe8321SBarry Smith   PetscErrorCode ierr;
43938baddfdSBarry Smith   PetscInt       M,N;
440b8f8c88eSHong Zhang   PetscMPIInt    size;
441639f9d9dSBarry Smith 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
443d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
444639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
44529bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
446639f9d9dSBarry Smith 
447f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
44852e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
449639f9d9dSBarry Smith 
450b8f8c88eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
451b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
452b8f8c88eSHong Zhang   if (size == 1) c->ctype = iscoloring->ctype = IS_COLORING_LOCAL;
453b8f8c88eSHong Zhang 
454f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
455f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
456639f9d9dSBarry Smith   } else {
45729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
458639f9d9dSBarry Smith   }
459639f9d9dSBarry Smith 
460b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
461b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
462b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
463b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
464b8f8c88eSHong Zhang 
46577d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
46677d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
467005c665bSBarry Smith   c->freq              = 1;
46840964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
46940964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
47049b058dcSBarry Smith   c->currentcolor      = -1;
471efb30889SBarry Smith   c->htype             = "wp";
472639f9d9dSBarry Smith 
473639f9d9dSBarry Smith   *color = c;
474d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4753a40ed3dSBarry Smith   PetscFunctionReturn(0);
476bbf0e169SBarry Smith }
477bbf0e169SBarry Smith 
4784a2ae208SSatish Balay #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
48005869f15SSatish Balay /*@
481639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
482639f9d9dSBarry Smith     via MatFDColoringCreate().
483bbf0e169SBarry Smith 
484ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
485ef5ee4d1SLois Curfman McInnes 
486b4fc646aSLois Curfman McInnes     Input Parameter:
487639f9d9dSBarry Smith .   c - coloring context
488bbf0e169SBarry Smith 
48915091d37SBarry Smith     Level: intermediate
49015091d37SBarry Smith 
491639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
492bbf0e169SBarry Smith @*/
493be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
494bbf0e169SBarry Smith {
4956849ba73SBarry Smith   PetscErrorCode ierr;
49638baddfdSBarry Smith   PetscInt       i;
497bbf0e169SBarry Smith 
4983a40ed3dSBarry Smith   PetscFunctionBegin;
4993a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
500d4bb536fSBarry Smith 
501639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
50205b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
50305b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
50405b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
50505b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
506bbf0e169SBarry Smith   }
507606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
508606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
509606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
510606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
511606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
51205b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
5134f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
514005c665bSBarry Smith   if (c->w1) {
515005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
516005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
517b8f8c88eSHong Zhang   }
518b8f8c88eSHong Zhang   if (c->w3){
519005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
520005c665bSBarry Smith   }
521d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5223a40ed3dSBarry Smith   PetscFunctionReturn(0);
523bbf0e169SBarry Smith }
52443a90d84SBarry Smith 
5254a2ae208SSatish Balay #undef __FUNCT__
52649b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
52749b058dcSBarry Smith /*@C
52849b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52949b058dcSBarry Smith       that are currently being perturbed.
53049b058dcSBarry Smith 
53149b058dcSBarry Smith     Not Collective
53249b058dcSBarry Smith 
53349b058dcSBarry Smith     Input Parameters:
53449b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
53549b058dcSBarry Smith 
53649b058dcSBarry Smith     Output Parameters:
53749b058dcSBarry Smith +   n - the number of local columns being perturbed
53849b058dcSBarry Smith -   cols - the column indices, in global numbering
53949b058dcSBarry Smith 
54049b058dcSBarry Smith    Level: intermediate
54149b058dcSBarry Smith 
54249b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
54349b058dcSBarry Smith 
54449b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
54549b058dcSBarry Smith @*/
546be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
54749b058dcSBarry Smith {
54849b058dcSBarry Smith   PetscFunctionBegin;
54949b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
55049b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
55149b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55249b058dcSBarry Smith   } else {
55349b058dcSBarry Smith     *n = 0;
55449b058dcSBarry Smith   }
55549b058dcSBarry Smith   PetscFunctionReturn(0);
55649b058dcSBarry Smith }
55749b058dcSBarry Smith 
558b8f8c88eSHong Zhang #include "petscda.h"      /*I      "petscda.h"    I*/
559b8f8c88eSHong Zhang 
56049b058dcSBarry Smith #undef __FUNCT__
5614a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
56243a90d84SBarry Smith /*@
563e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
564e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
56543a90d84SBarry Smith 
566fee21e36SBarry Smith     Collective on MatFDColoring
567fee21e36SBarry Smith 
568ef5ee4d1SLois Curfman McInnes     Input Parameters:
569ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
570ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
571ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
572ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
573ef5ee4d1SLois Curfman McInnes 
5748bba8e72SBarry Smith     Options Database Keys:
575b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
576efb30889SBarry Smith .    -mat_fd_type - "wp" or "ds"  (see MATSNESMF_WP or MATSNESMF_DS)
577b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
578b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
579b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5808bba8e72SBarry Smith 
58115091d37SBarry Smith     Level: intermediate
58215091d37SBarry Smith 
58343a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
58443a90d84SBarry Smith 
58543a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
58643a90d84SBarry Smith @*/
587be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
58843a90d84SBarry Smith {
5896849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5906849ba73SBarry Smith   PetscErrorCode ierr;
591b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
592efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
593ea709b57SSatish Balay   PetscScalar    *vscale_array;
594efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
595b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
596005c665bSBarry Smith   void           *fctx = coloring->fctx;
597f1af5d2fSBarry Smith   PetscTruth     flg;
598b88e453bSHong Zhang   PetscInt       ctype=coloring->ctype,N,col_start,col_end;
599b8f8c88eSHong Zhang   Vec            x1_tmp;
600*aeb7e98aSMatthew Knepley   /* remove ! */
601b8f8c88eSHong Zhang   PetscMPIInt rank;
602b8f8c88eSHong Zhang   PetscInt    prid=10;
60314b491afSHong Zhang   /*  ex5
604b8f8c88eSHong Zhang   PetscTruth  fd_jacobian_ghost=PETSC_FALSE;
605b8f8c88eSHong Zhang   DA          da;
60614b491afSHong Zhang   */
607a2e34c3dSBarry Smith 
6083a40ed3dSBarry Smith   PetscFunctionBegin;
609b8f8c88eSHong Zhang 
610b8f8c88eSHong Zhang     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
611b8f8c88eSHong Zhang     ierr = PetscOptionsGetInt(PETSC_NULL,"-prid",&prid,PETSC_NULL);CHKERRQ(ierr);
612b8f8c88eSHong Zhang 
6134482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
6144482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
6154482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
616feba9168SBarry Smith   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
617e0907662SLois Curfman McInnes 
61840964ad5SBarry Smith   if (coloring->usersetsrecompute) {
61940964ad5SBarry Smith     if (!coloring->recompute) {
62040964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
621ae15b995SBarry Smith       ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr);
62240964ad5SBarry Smith       PetscFunctionReturn(0);
62340964ad5SBarry Smith     } else {
62440964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
62540964ad5SBarry Smith     }
62640964ad5SBarry Smith   }
62740964ad5SBarry Smith 
628d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6294c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
630b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
631f1af5d2fSBarry Smith   if (flg) {
632ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
633e0907662SLois Curfman McInnes   } else {
6340b9b6f31SBarry Smith     PetscTruth assembled;
6350b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
6360b9b6f31SBarry Smith     if (assembled) {
63743a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
638e0907662SLois Curfman McInnes     }
6390b9b6f31SBarry Smith   }
64043a90d84SBarry Smith 
641b8f8c88eSHong Zhang   x1_tmp = x1;
642b8f8c88eSHong Zhang 
643b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED && !coloring->vscale){
644b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
645b8f8c88eSHong Zhang   }
646be2fbe1fSBarry Smith 
6473a7fca6bSBarry Smith   /*
6483a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
6493a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
6503a7fca6bSBarry Smith   */
6513a7fca6bSBarry Smith   if (coloring->F) {
6523a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
6533a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
6543a7fca6bSBarry Smith     if (m1 != m2) {
6553a7fca6bSBarry Smith       coloring->F = 0;
6563a7fca6bSBarry Smith       }
6573a7fca6bSBarry Smith     }
6583a7fca6bSBarry Smith 
659b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
660b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
661b8f8c88eSHong Zhang   }
662b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
663b8f8c88eSHong Zhang 
664b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
665be2fbe1fSBarry Smith   if (coloring->F) {
666be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
667be2fbe1fSBarry Smith     coloring->F = 0;
668be2fbe1fSBarry Smith   } else {
66966f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
670b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
67166f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
672be2fbe1fSBarry Smith   }
67343a90d84SBarry Smith 
674b8f8c88eSHong Zhang   if (!coloring->w3) {
675b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
676b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
677efb30889SBarry Smith   }
678b8f8c88eSHong Zhang   w3 = coloring->w3;
679efb30889SBarry Smith 
680b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
681b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
682b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
683b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
684b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
685b8f8c88eSHong Zhang     col_start = 0; col_end = N;
686b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_LOCAL){
687b8f8c88eSHong Zhang     xx = xx - start;
688b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
689b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
690b8f8c88eSHong Zhang   }
691b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
692b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
693efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
694efb30889SBarry Smith       dx = 1.0 + unorm;
695efb30889SBarry Smith     } else {
6969782cf98SBarry Smith       dx  = xx[col];
697efb30889SBarry Smith     }
6985b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6999782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
7009782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
7019782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
7029782cf98SBarry Smith #else
7039782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
7049782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
7059782cf98SBarry Smith #endif
7069782cf98SBarry Smith     dx               *= epsilon;
70730b34957SBarry Smith     vscale_array[col] = 1.0/dx;
7089782cf98SBarry Smith   }
709b8f8c88eSHong Zhang   if (ctype == IS_COLORING_LOCAL)  vscale_array = vscale_array + start;
710efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
711b8f8c88eSHong Zhang   if (ctype == IS_COLORING_LOCAL){
71230b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71330b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
714b8f8c88eSHong Zhang   }
7159782cf98SBarry Smith 
716b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
717b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
718b8f8c88eSHong Zhang   } else {
719b8f8c88eSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
720b8f8c88eSHong Zhang   }
721b0a32e0cSBarry Smith 
7229782cf98SBarry Smith   /*
72343a90d84SBarry Smith     Loop over each color
72443a90d84SBarry Smith   */
725b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
72643a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
72749b058dcSBarry Smith     coloring->currentcolor = k;
728b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
729b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
730b8f8c88eSHong Zhang     if (ctype == IS_COLORING_LOCAL) w3_array = w3_array - start;
731b8f8c88eSHong Zhang 
732b8f8c88eSHong Zhang     if (prid == rank) printf(" [%d] color %d \n -----------\n",rank,k);
73343a90d84SBarry Smith     /*
734b8f8c88eSHong Zhang       Loop over each column associated with color
735b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
73643a90d84SBarry Smith     */
73743a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
738b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
739efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
740efb30889SBarry Smith         dx = 1.0 + unorm;
741efb30889SBarry Smith       } else {
74242460c72SBarry Smith         dx  = xx[col];
743efb30889SBarry Smith       }
7445b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
745aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
746ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
747ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
74843a90d84SBarry Smith #else
749329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
750329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
75143a90d84SBarry Smith #endif
75243a90d84SBarry Smith       dx            *= epsilon;
7531302d50aSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
75442460c72SBarry Smith       w3_array[col] += dx;
75543a90d84SBarry Smith     }
756b8f8c88eSHong Zhang     if (ctype == IS_COLORING_LOCAL) w3_array = w3_array + start;
757b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7583b28642cSBarry Smith 
75943a90d84SBarry Smith     /*
760b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
761b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
76243a90d84SBarry Smith     */
76366f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
76443a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
76566f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
766efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
7679782cf98SBarry Smith 
76843a90d84SBarry Smith     /*
769e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
77043a90d84SBarry Smith     */
7713b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
77243a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
773b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
774b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
7755904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
77643a90d84SBarry Smith       srow   = row + start;
77743a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
77843a90d84SBarry Smith     }
7793b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
780*aeb7e98aSMatthew Knepley   } /* endof for each color */
781b8f8c88eSHong Zhang   if (ctype == IS_COLORING_LOCAL) xx = xx + start;
78230b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
783b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
784b8f8c88eSHong Zhang 
785b8f8c88eSHong Zhang   coloring->currentcolor = -1;
78643a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78743a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
788*aeb7e98aSMatthew Knepley   /*ierr = MatView(J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
789d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
790a2e34c3dSBarry Smith 
791b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
792a2e34c3dSBarry Smith   if (flg) {
793a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
794a2e34c3dSBarry Smith   }
795b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
796b9722508SLois Curfman McInnes 
7973a40ed3dSBarry Smith   PetscFunctionReturn(0);
79843a90d84SBarry Smith }
799840b8ebdSBarry Smith 
8004a2ae208SSatish Balay #undef __FUNCT__
8014a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
802840b8ebdSBarry Smith /*@
803840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
804840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
805840b8ebdSBarry Smith 
806fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
807fee21e36SBarry Smith 
808ef5ee4d1SLois Curfman McInnes     Input Parameters:
8093b28642cSBarry Smith +   mat - location to store Jacobian
810ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
811ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
812ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
813ef5ee4d1SLois Curfman McInnes 
814840b8ebdSBarry Smith    Options Database Keys:
815ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
816840b8ebdSBarry Smith 
81715091d37SBarry Smith    Level: intermediate
81815091d37SBarry Smith 
819840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
820840b8ebdSBarry Smith 
821840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
822840b8ebdSBarry Smith @*/
823be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
824840b8ebdSBarry Smith {
8256849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
8266849ba73SBarry Smith   PetscErrorCode ierr;
82738baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
828efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
829ea709b57SSatish Balay   PetscScalar    *vscale_array;
830329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
831840b8ebdSBarry Smith   Vec            w1,w2,w3;
832840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
833f1af5d2fSBarry Smith   PetscTruth     flg;
834840b8ebdSBarry Smith 
8353a40ed3dSBarry Smith   PetscFunctionBegin;
8364482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
8374482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
8384482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
839840b8ebdSBarry Smith 
840d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
841840b8ebdSBarry Smith   if (!coloring->w1) {
842840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
84352e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr);
844840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
84552e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr);
846840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
84752e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
848840b8ebdSBarry Smith   }
849840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
850840b8ebdSBarry Smith 
8515904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
852b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
853f1af5d2fSBarry Smith   if (flg) {
854ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
855840b8ebdSBarry Smith   } else {
8560b9b6f31SBarry Smith     PetscTruth assembled;
8570b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
8580b9b6f31SBarry Smith     if (assembled) {
859840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
860840b8ebdSBarry Smith     }
8610b9b6f31SBarry Smith   }
862840b8ebdSBarry Smith 
863840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
864840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
86566f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
866840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
86766f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
868840b8ebdSBarry Smith 
869840b8ebdSBarry Smith   /*
8705904e1b1SBarry Smith       Compute all the scale factors and share with other processors
871840b8ebdSBarry Smith   */
8725904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
8735904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
874840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
875840b8ebdSBarry Smith     /*
876840b8ebdSBarry Smith        Loop over each column associated with color adding the
877840b8ebdSBarry Smith        perturbation to the vector w3.
878840b8ebdSBarry Smith     */
879840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
880840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8815904e1b1SBarry Smith       dx  = xx[col];
8825b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
883aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
884840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
885840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
886840b8ebdSBarry Smith #else
887329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
888329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
889840b8ebdSBarry Smith #endif
890840b8ebdSBarry Smith       dx                *= epsilon;
8915904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
892840b8ebdSBarry Smith     }
8935904e1b1SBarry Smith   }
8945904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8955904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8965904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8975904e1b1SBarry Smith 
8985904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8995904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
9005904e1b1SBarry Smith 
9015904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
9025904e1b1SBarry Smith   /*
9035904e1b1SBarry Smith       Loop over each color
9045904e1b1SBarry Smith   */
9055904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
9065904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
9075904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
9085904e1b1SBarry Smith     /*
9095904e1b1SBarry Smith        Loop over each column associated with color adding the
9105904e1b1SBarry Smith        perturbation to the vector w3.
9115904e1b1SBarry Smith     */
9125904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
9135904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
9145904e1b1SBarry Smith       dx  = xx[col];
9155b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
9165904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9175904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
9185904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
9195904e1b1SBarry Smith #else
9205904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
9215904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
9225904e1b1SBarry Smith #endif
9235904e1b1SBarry Smith       dx            *= epsilon;
9245904e1b1SBarry Smith       w3_array[col] += dx;
9255904e1b1SBarry Smith     }
9265904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
9275904e1b1SBarry Smith 
928840b8ebdSBarry Smith     /*
929840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
930840b8ebdSBarry Smith     */
93166f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
932840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
93366f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
934efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
9355904e1b1SBarry Smith 
936840b8ebdSBarry Smith     /*
937840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
938840b8ebdSBarry Smith     */
9393b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
940840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
941840b8ebdSBarry Smith       row    = coloring->rows[k][l];
942840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
9435904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
944840b8ebdSBarry Smith       srow   = row + start;
945840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
946840b8ebdSBarry Smith     }
9473b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
948840b8ebdSBarry Smith   }
9495904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
9505904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
951840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
952840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
953d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
9543a40ed3dSBarry Smith   PetscFunctionReturn(0);
955840b8ebdSBarry Smith }
9563b28642cSBarry Smith 
9573b28642cSBarry Smith 
9584a2ae208SSatish Balay #undef __FUNCT__
9594a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
96040964ad5SBarry Smith /*@C
96140964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
96240964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
96340964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
96440964ad5SBarry Smith      called again.
96540964ad5SBarry Smith 
96640964ad5SBarry Smith    Collective on MatFDColoring
96740964ad5SBarry Smith 
96840964ad5SBarry Smith    Input  Parameters:
96940964ad5SBarry Smith .  c - the coloring context
97040964ad5SBarry Smith 
97140964ad5SBarry Smith    Level: intermediate
97240964ad5SBarry Smith 
97340964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
97440964ad5SBarry Smith 
97540964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
97640964ad5SBarry Smith 
97740964ad5SBarry Smith .keywords: Mat, finite differences, coloring
97840964ad5SBarry Smith @*/
979be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
98040964ad5SBarry Smith {
98140964ad5SBarry Smith   PetscFunctionBegin;
9824482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
98340964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
98440964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
98540964ad5SBarry Smith   PetscFunctionReturn(0);
98640964ad5SBarry Smith }
98740964ad5SBarry Smith 
9883b28642cSBarry Smith 
989