xref: /petsc/src/mat/matfd/fdmatrix.c (revision 3050cee28ff83b093d983d7c5197e87b438ca09b)
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 
8b9147fbbSdalcinl #include "include/private/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);
107*3050cee2SBarry Smith   if (!viewer) {
108*3050cee2SBarry Smith     ierr = PetscViewerASCIIGetStdout(c->comm,&viewer);CHKERRQ(ierr);
109*3050cee2SBarry Smith   }
1104482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
111c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
112bbf0e169SBarry Smith 
113fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
11432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1150f5bd95cSBarry Smith   if (isdraw) {
116b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11732077d6dSBarry Smith   } else if (iascii) {
118b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
119a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
120a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
12177431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
122ae09f205SBarry Smith 
123b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
124fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
125b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12677431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
128b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12977431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
130639f9d9dSBarry Smith         }
13177431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
132b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
13377431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
134b4fc646aSLois Curfman McInnes         }
135bbf0e169SBarry Smith       }
136bbf0e169SBarry Smith     }
137b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1385cd90555SBarry Smith   } else {
13979a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
140bbf0e169SBarry Smith   }
1413a40ed3dSBarry Smith   PetscFunctionReturn(0);
142639f9d9dSBarry Smith }
143639f9d9dSBarry Smith 
1444a2ae208SSatish Balay #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
146639f9d9dSBarry Smith /*@
147b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
148b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
149639f9d9dSBarry Smith 
150ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
151ef5ee4d1SLois Curfman McInnes 
152ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
153ef5ee4d1SLois Curfman McInnes .vb
15465f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
155f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
156f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
157ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
158ef5ee4d1SLois Curfman McInnes .ve
159639f9d9dSBarry Smith 
160639f9d9dSBarry Smith    Input Parameters:
161ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
162639f9d9dSBarry Smith .  error_rel - relative error
163f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
164fee21e36SBarry Smith 
16515091d37SBarry Smith    Level: advanced
16615091d37SBarry Smith 
167b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
168b4fc646aSLois Curfman McInnes 
169b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
170639f9d9dSBarry Smith @*/
171be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
172639f9d9dSBarry Smith {
1733a40ed3dSBarry Smith   PetscFunctionBegin;
1744482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
175639f9d9dSBarry Smith 
176639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
177639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1783a40ed3dSBarry Smith   PetscFunctionReturn(0);
179639f9d9dSBarry Smith }
180639f9d9dSBarry Smith 
1814a2ae208SSatish Balay #undef __FUNCT__
1824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency"
183005c665bSBarry Smith /*@
184e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
185e0907662SLois Curfman McInnes    matrices.
186005c665bSBarry Smith 
187fee21e36SBarry Smith    Collective on MatFDColoring
188fee21e36SBarry Smith 
189ef5ee4d1SLois Curfman McInnes    Input Parameters:
190ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
191ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
192ef5ee4d1SLois Curfman McInnes 
19315091d37SBarry Smith    Options Database Keys:
19415091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
19515091d37SBarry Smith 
19615091d37SBarry Smith    Level: advanced
19715091d37SBarry Smith 
198e0907662SLois Curfman McInnes    Notes:
199e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
200e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
201e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
202e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
203e0907662SLois Curfman McInnes 
204b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
205ef5ee4d1SLois Curfman McInnes 
20640964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute()
207005c665bSBarry Smith @*/
208be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq)
209005c665bSBarry Smith {
2103a40ed3dSBarry Smith   PetscFunctionBegin;
2114482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
212005c665bSBarry Smith 
213005c665bSBarry Smith   matfd->freq = freq;
2143a40ed3dSBarry Smith   PetscFunctionReturn(0);
215005c665bSBarry Smith }
216005c665bSBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency"
219ff0cfa39SBarry Smith /*@
220ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
221ff0cfa39SBarry Smith    matrices.
222ff0cfa39SBarry Smith 
223ef5ee4d1SLois Curfman McInnes    Not Collective
224ef5ee4d1SLois Curfman McInnes 
225ff0cfa39SBarry Smith    Input Parameters:
226ff0cfa39SBarry Smith .  coloring - the coloring context
227ff0cfa39SBarry Smith 
228ff0cfa39SBarry Smith    Output Parameters:
229ff0cfa39SBarry Smith .  freq - frequency (default is 1)
230ff0cfa39SBarry Smith 
23115091d37SBarry Smith    Options Database Keys:
23215091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
23315091d37SBarry Smith 
23415091d37SBarry Smith    Level: advanced
23515091d37SBarry Smith 
236ff0cfa39SBarry Smith    Notes:
237ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
238ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
239ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
240ff0cfa39SBarry Smith    <freq> nonlinear iterations.
241ff0cfa39SBarry Smith 
242ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
243ef5ee4d1SLois Curfman McInnes 
244ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
245ff0cfa39SBarry Smith @*/
246be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq)
247ff0cfa39SBarry Smith {
2483a40ed3dSBarry Smith   PetscFunctionBegin;
2494482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
250ff0cfa39SBarry Smith   *freq = matfd->freq;
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
252ff0cfa39SBarry Smith }
253ff0cfa39SBarry Smith 
2544a2ae208SSatish Balay #undef __FUNCT__
2554a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
2564a9d489dSBarry Smith /*@C
2574a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2584a9d489dSBarry Smith 
2594a9d489dSBarry Smith    Collective on MatFDColoring
2604a9d489dSBarry Smith 
2614a9d489dSBarry Smith    Input Parameters:
2624a9d489dSBarry Smith .  coloring - the coloring context
2634a9d489dSBarry Smith 
2644a9d489dSBarry Smith    Output Parameters:
2654a9d489dSBarry Smith +  f - the function
2664a9d489dSBarry Smith -  fctx - the optional user-defined function context
2674a9d489dSBarry Smith 
2684a9d489dSBarry Smith    Level: intermediate
2694a9d489dSBarry Smith 
2704a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
2714a9d489dSBarry Smith @*/
2724a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2734a9d489dSBarry Smith {
2744a9d489dSBarry Smith   PetscFunctionBegin;
2754a9d489dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
2764a9d489dSBarry Smith   if (f) *f = matfd->f;
2774a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2784a9d489dSBarry Smith   PetscFunctionReturn(0);
2794a9d489dSBarry Smith }
2804a9d489dSBarry Smith 
2814a9d489dSBarry Smith #undef __FUNCT__
2824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
283d64ed03dSBarry Smith /*@C
284005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
285005c665bSBarry Smith 
286fee21e36SBarry Smith    Collective on MatFDColoring
287fee21e36SBarry Smith 
288ef5ee4d1SLois Curfman McInnes    Input Parameters:
289ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
290ef5ee4d1SLois Curfman McInnes .  f - the function
291ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
292ef5ee4d1SLois Curfman McInnes 
29315091d37SBarry Smith    Level: intermediate
29415091d37SBarry Smith 
295f881d145SBarry Smith    Notes:
296f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
297f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
298f881d145SBarry Smith   with the TS solvers.
299f881d145SBarry Smith 
300b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
301005c665bSBarry Smith @*/
302be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
303005c665bSBarry Smith {
3043a40ed3dSBarry Smith   PetscFunctionBegin;
3054482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
306005c665bSBarry Smith   matfd->f    = f;
307005c665bSBarry Smith   matfd->fctx = fctx;
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
309005c665bSBarry Smith }
310005c665bSBarry Smith 
3114a2ae208SSatish Balay #undef __FUNCT__
3124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
313639f9d9dSBarry Smith /*@
314b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
315639f9d9dSBarry Smith    the options database.
316639f9d9dSBarry Smith 
317fee21e36SBarry Smith    Collective on MatFDColoring
318fee21e36SBarry Smith 
31965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
320ef5ee4d1SLois Curfman McInnes .vb
32165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
322f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
323f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
324ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
325ef5ee4d1SLois Curfman McInnes .ve
326ef5ee4d1SLois Curfman McInnes 
327ef5ee4d1SLois Curfman McInnes    Input Parameter:
328ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
329ef5ee4d1SLois Curfman McInnes 
330b4fc646aSLois Curfman McInnes    Options Database Keys:
331d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
332ef5ee4d1SLois Curfman McInnes            of relative error in the function)
333f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
334ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
335efb30889SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS)
336ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
337ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
338ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
339639f9d9dSBarry Smith 
34015091d37SBarry Smith     Level: intermediate
34115091d37SBarry Smith 
342b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
343d1c39f55SBarry Smith 
344d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
345d1c39f55SBarry Smith 
346639f9d9dSBarry Smith @*/
347be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
348639f9d9dSBarry Smith {
349dfbe8321SBarry Smith   PetscErrorCode ierr;
350efb30889SBarry Smith   PetscTruth     flg;
351efb30889SBarry Smith   char           value[3];
3523a40ed3dSBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
3544482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
355639f9d9dSBarry Smith 
356b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
35787828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
35887828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
359b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
360efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
361efb30889SBarry Smith     if (flg) {
362efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
363efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
364efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
365efb30889SBarry Smith     }
366186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
367b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
368b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
369b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
370b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3713a40ed3dSBarry Smith   PetscFunctionReturn(0);
372005c665bSBarry Smith }
373005c665bSBarry Smith 
3744a2ae208SSatish Balay #undef __FUNCT__
3754a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
376dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
377005c665bSBarry Smith {
378dfbe8321SBarry Smith   PetscErrorCode ierr;
379f1af5d2fSBarry Smith   PetscTruth     flg;
380*3050cee2SBarry Smith   PetscViewer    viewer;
381005c665bSBarry Smith 
3823a40ed3dSBarry Smith   PetscFunctionBegin;
383*3050cee2SBarry Smith   ierr = PetscViewerASCIIGetStdout(fd->comm,&viewer);CHKERRQ(ierr);
384b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
385005c665bSBarry Smith   if (flg) {
386*3050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
387005c665bSBarry Smith   }
388b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
389ae09f205SBarry Smith   if (flg) {
390*3050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
391*3050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
392*3050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
393ae09f205SBarry Smith   }
394b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
395005c665bSBarry Smith   if (flg) {
396b0a32e0cSBarry Smith     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
397b0a32e0cSBarry Smith     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
398005c665bSBarry Smith   }
3993a40ed3dSBarry Smith   PetscFunctionReturn(0);
400bbf0e169SBarry Smith }
401bbf0e169SBarry Smith 
4024a2ae208SSatish Balay #undef __FUNCT__
4034a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
40405869f15SSatish Balay /*@
405639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
406639f9d9dSBarry Smith    computation of Jacobians.
407bbf0e169SBarry Smith 
408ef5ee4d1SLois Curfman McInnes    Collective on Mat
409ef5ee4d1SLois Curfman McInnes 
410639f9d9dSBarry Smith    Input Parameters:
411ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
412ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
413bbf0e169SBarry Smith 
414bbf0e169SBarry Smith     Output Parameter:
415639f9d9dSBarry Smith .   color - the new coloring context
416bbf0e169SBarry Smith 
41715091d37SBarry Smith     Level: intermediate
41815091d37SBarry Smith 
4196831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
420d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
421d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
422d1c39f55SBarry Smith           MatFDColoringSetParameters()
423bbf0e169SBarry Smith @*/
424be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
425bbf0e169SBarry Smith {
426639f9d9dSBarry Smith   MatFDColoring  c;
427639f9d9dSBarry Smith   MPI_Comm       comm;
428dfbe8321SBarry Smith   PetscErrorCode ierr;
42938baddfdSBarry Smith   PetscInt       M,N;
430b8f8c88eSHong Zhang   PetscMPIInt    size;
431639f9d9dSBarry Smith 
4323a40ed3dSBarry Smith   PetscFunctionBegin;
433d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
434639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
43529bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
436639f9d9dSBarry Smith 
437f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
43852e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
439639f9d9dSBarry Smith 
440b8f8c88eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
441b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
442b8f8c88eSHong Zhang 
443f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
444f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
445639f9d9dSBarry Smith   } else {
44629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
447639f9d9dSBarry Smith   }
448639f9d9dSBarry Smith 
449b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
450b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
451b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
452b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
453b8f8c88eSHong Zhang 
45477d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
45577d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
456005c665bSBarry Smith   c->freq              = 1;
45740964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
45840964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
45949b058dcSBarry Smith   c->currentcolor      = -1;
460efb30889SBarry Smith   c->htype             = "wp";
461639f9d9dSBarry Smith 
462639f9d9dSBarry Smith   *color = c;
463d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
465bbf0e169SBarry Smith }
466bbf0e169SBarry Smith 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
46905869f15SSatish Balay /*@
470639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
471639f9d9dSBarry Smith     via MatFDColoringCreate().
472bbf0e169SBarry Smith 
473ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
474ef5ee4d1SLois Curfman McInnes 
475b4fc646aSLois Curfman McInnes     Input Parameter:
476639f9d9dSBarry Smith .   c - coloring context
477bbf0e169SBarry Smith 
47815091d37SBarry Smith     Level: intermediate
47915091d37SBarry Smith 
480639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
481bbf0e169SBarry Smith @*/
482be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
483bbf0e169SBarry Smith {
4846849ba73SBarry Smith   PetscErrorCode ierr;
48538baddfdSBarry Smith   PetscInt       i;
486bbf0e169SBarry Smith 
4873a40ed3dSBarry Smith   PetscFunctionBegin;
4883a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
489d4bb536fSBarry Smith 
490639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
49105b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
49205b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
49305b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
49405b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
495bbf0e169SBarry Smith   }
496606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
497606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
498606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
499606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
500606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
50105b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
5024f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
503005c665bSBarry Smith   if (c->w1) {
504005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
505005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
506b8f8c88eSHong Zhang   }
507b8f8c88eSHong Zhang   if (c->w3){
508005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
509005c665bSBarry Smith   }
510d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
512bbf0e169SBarry Smith }
51343a90d84SBarry Smith 
5144a2ae208SSatish Balay #undef __FUNCT__
51549b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
51649b058dcSBarry Smith /*@C
51749b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
51849b058dcSBarry Smith       that are currently being perturbed.
51949b058dcSBarry Smith 
52049b058dcSBarry Smith     Not Collective
52149b058dcSBarry Smith 
52249b058dcSBarry Smith     Input Parameters:
52349b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
52449b058dcSBarry Smith 
52549b058dcSBarry Smith     Output Parameters:
52649b058dcSBarry Smith +   n - the number of local columns being perturbed
52749b058dcSBarry Smith -   cols - the column indices, in global numbering
52849b058dcSBarry Smith 
52949b058dcSBarry Smith    Level: intermediate
53049b058dcSBarry Smith 
53149b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
53249b058dcSBarry Smith 
53349b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
53449b058dcSBarry Smith @*/
535be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
53649b058dcSBarry Smith {
53749b058dcSBarry Smith   PetscFunctionBegin;
53849b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
53949b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
54049b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
54149b058dcSBarry Smith   } else {
54249b058dcSBarry Smith     *n = 0;
54349b058dcSBarry Smith   }
54449b058dcSBarry Smith   PetscFunctionReturn(0);
54549b058dcSBarry Smith }
54649b058dcSBarry Smith 
547b8f8c88eSHong Zhang #include "petscda.h"      /*I      "petscda.h"    I*/
548b8f8c88eSHong Zhang 
54949b058dcSBarry Smith #undef __FUNCT__
5504a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
55143a90d84SBarry Smith /*@
552e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
553e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
55443a90d84SBarry Smith 
555fee21e36SBarry Smith     Collective on MatFDColoring
556fee21e36SBarry Smith 
557ef5ee4d1SLois Curfman McInnes     Input Parameters:
558ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
559ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
560ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
561ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
562ef5ee4d1SLois Curfman McInnes 
5638bba8e72SBarry Smith     Options Database Keys:
564b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
565efb30889SBarry Smith .    -mat_fd_type - "wp" or "ds"  (see MATSNESMF_WP or MATSNESMF_DS)
566b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
567b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
568b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5698bba8e72SBarry Smith 
57015091d37SBarry Smith     Level: intermediate
57115091d37SBarry Smith 
57243a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
57343a90d84SBarry Smith 
57443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
57543a90d84SBarry Smith @*/
576be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
57743a90d84SBarry Smith {
5786849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5796849ba73SBarry Smith   PetscErrorCode ierr;
580b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
581efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
582ea709b57SSatish Balay   PetscScalar    *vscale_array;
583efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
584b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
585005c665bSBarry Smith   void           *fctx = coloring->fctx;
586f1af5d2fSBarry Smith   PetscTruth     flg;
587705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
588b8f8c88eSHong Zhang   Vec            x1_tmp;
589a2e34c3dSBarry Smith 
5903a40ed3dSBarry Smith   PetscFunctionBegin;
5914482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
5924482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
5934482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
594feba9168SBarry Smith   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
595e0907662SLois Curfman McInnes 
59640964ad5SBarry Smith   if (coloring->usersetsrecompute) {
59740964ad5SBarry Smith     if (!coloring->recompute) {
59840964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
599ae15b995SBarry Smith       ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr);
60040964ad5SBarry Smith       PetscFunctionReturn(0);
60140964ad5SBarry Smith     } else {
60240964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
60340964ad5SBarry Smith     }
60440964ad5SBarry Smith   }
60540964ad5SBarry Smith 
606d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6074c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
608b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
609f1af5d2fSBarry Smith   if (flg) {
610ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
611e0907662SLois Curfman McInnes   } else {
6120b9b6f31SBarry Smith     PetscTruth assembled;
6130b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
6140b9b6f31SBarry Smith     if (assembled) {
61543a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
616e0907662SLois Curfman McInnes     }
6170b9b6f31SBarry Smith   }
61843a90d84SBarry Smith 
619b8f8c88eSHong Zhang   x1_tmp = x1;
620dac7f5fdSBarry Smith   if (!coloring->vscale){
621b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
622b8f8c88eSHong Zhang   }
623be2fbe1fSBarry Smith 
6243a7fca6bSBarry Smith   /*
6253a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
6263a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
6273a7fca6bSBarry Smith   */
6283a7fca6bSBarry Smith   if (coloring->F) {
6293a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
6303a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
6313a7fca6bSBarry Smith     if (m1 != m2) {
6323a7fca6bSBarry Smith       coloring->F = 0;
6333a7fca6bSBarry Smith       }
6343a7fca6bSBarry Smith     }
6353a7fca6bSBarry Smith 
636b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
637b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
638b8f8c88eSHong Zhang   }
639b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
640b8f8c88eSHong Zhang 
641b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
642be2fbe1fSBarry Smith   if (coloring->F) {
643be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
644be2fbe1fSBarry Smith     coloring->F = 0;
645be2fbe1fSBarry Smith   } else {
64666f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
647b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
64866f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
649be2fbe1fSBarry Smith   }
65043a90d84SBarry Smith 
651b8f8c88eSHong Zhang   if (!coloring->w3) {
652b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
653b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
654efb30889SBarry Smith   }
655b8f8c88eSHong Zhang   w3 = coloring->w3;
656efb30889SBarry Smith 
657b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
658b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
659b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
660b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
661b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
662b8f8c88eSHong Zhang     col_start = 0; col_end = N;
6638ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
664b8f8c88eSHong Zhang     xx = xx - start;
665b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
666b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
667b8f8c88eSHong Zhang   }
668b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
669b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
670efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
671efb30889SBarry Smith       dx = 1.0 + unorm;
672efb30889SBarry Smith     } else {
6739782cf98SBarry Smith       dx  = xx[col];
674efb30889SBarry Smith     }
6755b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6769782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6779782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6789782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6799782cf98SBarry Smith #else
6809782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6819782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6829782cf98SBarry Smith #endif
6839782cf98SBarry Smith     dx               *= epsilon;
68430b34957SBarry Smith     vscale_array[col] = 1.0/dx;
6859782cf98SBarry Smith   }
6868ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
687efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6888ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
68930b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69030b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
691b8f8c88eSHong Zhang   }
6929782cf98SBarry Smith 
693b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
694b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
695b8f8c88eSHong Zhang   } else {
696b8f8c88eSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
697b8f8c88eSHong Zhang   }
698b0a32e0cSBarry Smith 
6999782cf98SBarry Smith   /*
70043a90d84SBarry Smith     Loop over each color
70143a90d84SBarry Smith   */
702b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
70343a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
70449b058dcSBarry Smith     coloring->currentcolor = k;
705b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
706b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
7078ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
70843a90d84SBarry Smith     /*
709b8f8c88eSHong Zhang       Loop over each column associated with color
710b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
71143a90d84SBarry Smith     */
71243a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
713b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
714efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
715efb30889SBarry Smith         dx = 1.0 + unorm;
716efb30889SBarry Smith       } else {
71742460c72SBarry Smith         dx  = xx[col];
718efb30889SBarry Smith       }
7195b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
720aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
721ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
722ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
72343a90d84SBarry Smith #else
724329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
725329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
72643a90d84SBarry Smith #endif
72743a90d84SBarry Smith       dx            *= epsilon;
7281302d50aSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
72942460c72SBarry Smith       w3_array[col] += dx;
73043a90d84SBarry Smith     }
7318ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
732b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7333b28642cSBarry Smith 
73443a90d84SBarry Smith     /*
735b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
736b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
73743a90d84SBarry Smith     */
73866f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
73943a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
74066f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
741efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
7429782cf98SBarry Smith 
74343a90d84SBarry Smith     /*
744e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
74543a90d84SBarry Smith     */
7463b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
74743a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
748b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
749b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
7505904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
75143a90d84SBarry Smith       srow   = row + start;
75243a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
75343a90d84SBarry Smith     }
7543b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
755aeb7e98aSMatthew Knepley   } /* endof for each color */
7568ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
75730b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
758b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
759b8f8c88eSHong Zhang 
760b8f8c88eSHong Zhang   coloring->currentcolor = -1;
76143a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76243a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
763d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
764a2e34c3dSBarry Smith 
765b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
766a2e34c3dSBarry Smith   if (flg) {
767a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
768a2e34c3dSBarry Smith   }
769b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7703a40ed3dSBarry Smith   PetscFunctionReturn(0);
77143a90d84SBarry Smith }
772840b8ebdSBarry Smith 
7734a2ae208SSatish Balay #undef __FUNCT__
7744a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
775840b8ebdSBarry Smith /*@
776840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
777840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
778840b8ebdSBarry Smith 
779fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
780fee21e36SBarry Smith 
781ef5ee4d1SLois Curfman McInnes     Input Parameters:
7823b28642cSBarry Smith +   mat - location to store Jacobian
783ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
784ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
785ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
786ef5ee4d1SLois Curfman McInnes 
787840b8ebdSBarry Smith    Options Database Keys:
788ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
789840b8ebdSBarry Smith 
79015091d37SBarry Smith    Level: intermediate
79115091d37SBarry Smith 
792840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
793840b8ebdSBarry Smith 
794840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
795840b8ebdSBarry Smith @*/
796be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
797840b8ebdSBarry Smith {
7986849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
7996849ba73SBarry Smith   PetscErrorCode ierr;
80038baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
801efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
802ea709b57SSatish Balay   PetscScalar    *vscale_array;
803329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
804ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
805840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
806f1af5d2fSBarry Smith   PetscTruth     flg;
807840b8ebdSBarry Smith 
8083a40ed3dSBarry Smith   PetscFunctionBegin;
8094482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
8104482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
8114482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
812840b8ebdSBarry Smith 
813d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
814ced766e8SHong Zhang   if (!coloring->w3) {
815840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
81652e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
817840b8ebdSBarry Smith   }
818ced766e8SHong Zhang   w3 = coloring->w3;
819840b8ebdSBarry Smith 
8205904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
821b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
822f1af5d2fSBarry Smith   if (flg) {
823ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
824840b8ebdSBarry Smith   } else {
8250b9b6f31SBarry Smith     PetscTruth assembled;
8260b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
8270b9b6f31SBarry Smith     if (assembled) {
828840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
829840b8ebdSBarry Smith     }
8300b9b6f31SBarry Smith   }
831840b8ebdSBarry Smith 
832840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
833840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
83466f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
835840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
83666f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
837840b8ebdSBarry Smith 
838840b8ebdSBarry Smith   /*
8395904e1b1SBarry Smith       Compute all the scale factors and share with other processors
840840b8ebdSBarry Smith   */
8415904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
8425904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
843840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
844840b8ebdSBarry Smith     /*
845840b8ebdSBarry Smith        Loop over each column associated with color adding the
846840b8ebdSBarry Smith        perturbation to the vector w3.
847840b8ebdSBarry Smith     */
848840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
849840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8505904e1b1SBarry Smith       dx  = xx[col];
8515b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
852aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
853840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
854840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
855840b8ebdSBarry Smith #else
856329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
857329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
858840b8ebdSBarry Smith #endif
859840b8ebdSBarry Smith       dx                *= epsilon;
8605904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
861840b8ebdSBarry Smith     }
8625904e1b1SBarry Smith   }
8635904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8645904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8655904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8665904e1b1SBarry Smith 
8675904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8685904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8695904e1b1SBarry Smith 
8705904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8715904e1b1SBarry Smith   /*
8725904e1b1SBarry Smith       Loop over each color
8735904e1b1SBarry Smith   */
8745904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8755904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8765904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8775904e1b1SBarry Smith     /*
8785904e1b1SBarry Smith        Loop over each column associated with color adding the
8795904e1b1SBarry Smith        perturbation to the vector w3.
8805904e1b1SBarry Smith     */
8815904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8825904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8835904e1b1SBarry Smith       dx  = xx[col];
8845b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8855904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8865904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8875904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8885904e1b1SBarry Smith #else
8895904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8905904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
8915904e1b1SBarry Smith #endif
8925904e1b1SBarry Smith       dx            *= epsilon;
8935904e1b1SBarry Smith       w3_array[col] += dx;
8945904e1b1SBarry Smith     }
8955904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
8965904e1b1SBarry Smith 
897840b8ebdSBarry Smith     /*
898840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
899840b8ebdSBarry Smith     */
90066f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
901840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
90266f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
903efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
9045904e1b1SBarry Smith 
905840b8ebdSBarry Smith     /*
906840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
907840b8ebdSBarry Smith     */
9083b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
909840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
910840b8ebdSBarry Smith       row    = coloring->rows[k][l];
911840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
9125904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
913840b8ebdSBarry Smith       srow   = row + start;
914840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
915840b8ebdSBarry Smith     }
9163b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
917840b8ebdSBarry Smith   }
9185904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
9195904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
920840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
921840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
922d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
924840b8ebdSBarry Smith }
9253b28642cSBarry Smith 
9263b28642cSBarry Smith 
9274a2ae208SSatish Balay #undef __FUNCT__
9284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
92940964ad5SBarry Smith /*@C
93040964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
93140964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
93240964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
93340964ad5SBarry Smith      called again.
93440964ad5SBarry Smith 
93540964ad5SBarry Smith    Collective on MatFDColoring
93640964ad5SBarry Smith 
93740964ad5SBarry Smith    Input  Parameters:
93840964ad5SBarry Smith .  c - the coloring context
93940964ad5SBarry Smith 
94040964ad5SBarry Smith    Level: intermediate
94140964ad5SBarry Smith 
94240964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
94340964ad5SBarry Smith 
94440964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
94540964ad5SBarry Smith 
94640964ad5SBarry Smith .keywords: Mat, finite differences, coloring
94740964ad5SBarry Smith @*/
948be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
94940964ad5SBarry Smith {
95040964ad5SBarry Smith   PetscFunctionBegin;
9514482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
95240964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
95340964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
95440964ad5SBarry Smith   PetscFunctionReturn(0);
95540964ad5SBarry Smith }
95640964ad5SBarry Smith 
9573b28642cSBarry Smith 
958