xref: /petsc/src/mat/matfd/fdmatrix.c (revision 166c7f255c91bb243aa54e4217727fe58a0f23be)
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 */
11*166c7f25SBarry Smith PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
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);
1073050cee2SBarry Smith   if (!viewer) {
1087adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1093050cee2SBarry 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 
2937850c7c0SBarry Smith    Calling sequence of (*f) function:
2947850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
2957850c7c0SBarry Smith     For TS:      PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*)
2967850c7c0SBarry Smith     If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
29715091d37SBarry Smith 
2987850c7c0SBarry Smith    Level: advanced
2997850c7c0SBarry Smith 
3007850c7c0SBarry Smith    Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument
3017850c7c0SBarry Smith      SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used
3027850c7c0SBarry Smith      by someone computing a matrix via coloring directly by calling MatFDColoringApply()
3037850c7c0SBarry Smith 
3047850c7c0SBarry Smith    Fortran Notes:
3057850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
3067850c7c0SBarry Smith   be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
3077850c7c0SBarry Smith   within the TS solvers.
308f881d145SBarry Smith 
309b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
310005c665bSBarry Smith @*/
311be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
312005c665bSBarry Smith {
3133a40ed3dSBarry Smith   PetscFunctionBegin;
3144482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
315005c665bSBarry Smith   matfd->f    = f;
316005c665bSBarry Smith   matfd->fctx = fctx;
3173a40ed3dSBarry Smith   PetscFunctionReturn(0);
318005c665bSBarry Smith }
319005c665bSBarry Smith 
3204a2ae208SSatish Balay #undef __FUNCT__
3214a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
322639f9d9dSBarry Smith /*@
323b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
324639f9d9dSBarry Smith    the options database.
325639f9d9dSBarry Smith 
326fee21e36SBarry Smith    Collective on MatFDColoring
327fee21e36SBarry Smith 
32865f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
329ef5ee4d1SLois Curfman McInnes .vb
33065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
331f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
332f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
333ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
334ef5ee4d1SLois Curfman McInnes .ve
335ef5ee4d1SLois Curfman McInnes 
336ef5ee4d1SLois Curfman McInnes    Input Parameter:
337ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
338ef5ee4d1SLois Curfman McInnes 
339b4fc646aSLois Curfman McInnes    Options Database Keys:
340d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
341ef5ee4d1SLois Curfman McInnes            of relative error in the function)
342f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
343ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
3443ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
345ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
346ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
347ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
348639f9d9dSBarry Smith 
34915091d37SBarry Smith     Level: intermediate
35015091d37SBarry Smith 
351b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
352d1c39f55SBarry Smith 
353d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
354d1c39f55SBarry Smith 
355639f9d9dSBarry Smith @*/
356be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd)
357639f9d9dSBarry Smith {
358dfbe8321SBarry Smith   PetscErrorCode ierr;
359efb30889SBarry Smith   PetscTruth     flg;
360efb30889SBarry Smith   char           value[3];
3613a40ed3dSBarry Smith 
3623a40ed3dSBarry Smith   PetscFunctionBegin;
3634482741eSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1);
364639f9d9dSBarry Smith 
3657adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)matfd)->comm,((PetscObject)matfd)->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
36687828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
36787828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
368b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
369efb30889SBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr);
370efb30889SBarry Smith     if (flg) {
371efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
372efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
373efb30889SBarry Smith       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
374efb30889SBarry Smith     }
375186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
376b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
377b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
378b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
379b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
381005c665bSBarry Smith }
382005c665bSBarry Smith 
3834a2ae208SSatish Balay #undef __FUNCT__
3844a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
385dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
386005c665bSBarry Smith {
387dfbe8321SBarry Smith   PetscErrorCode ierr;
388f1af5d2fSBarry Smith   PetscTruth     flg;
3893050cee2SBarry Smith   PetscViewer    viewer;
390005c665bSBarry Smith 
3913a40ed3dSBarry Smith   PetscFunctionBegin;
3927adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
393b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
394005c665bSBarry Smith   if (flg) {
3953050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
396005c665bSBarry Smith   }
397b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
398ae09f205SBarry Smith   if (flg) {
3993050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
4003050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
4013050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
402ae09f205SBarry Smith   }
403b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
404005c665bSBarry Smith   if (flg) {
4057adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
4067adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
407005c665bSBarry Smith   }
4083a40ed3dSBarry Smith   PetscFunctionReturn(0);
409bbf0e169SBarry Smith }
410bbf0e169SBarry Smith 
4114a2ae208SSatish Balay #undef __FUNCT__
4124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
41305869f15SSatish Balay /*@
414639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
415639f9d9dSBarry Smith    computation of Jacobians.
416bbf0e169SBarry Smith 
417ef5ee4d1SLois Curfman McInnes    Collective on Mat
418ef5ee4d1SLois Curfman McInnes 
419639f9d9dSBarry Smith    Input Parameters:
420ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
421ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
422bbf0e169SBarry Smith 
423bbf0e169SBarry Smith     Output Parameter:
424639f9d9dSBarry Smith .   color - the new coloring context
425bbf0e169SBarry Smith 
42615091d37SBarry Smith     Level: intermediate
42715091d37SBarry Smith 
4286831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
429d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
430d1c39f55SBarry Smith           MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(),
431d1c39f55SBarry Smith           MatFDColoringSetParameters()
432bbf0e169SBarry Smith @*/
433be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
434bbf0e169SBarry Smith {
435639f9d9dSBarry Smith   MatFDColoring  c;
436639f9d9dSBarry Smith   MPI_Comm       comm;
437dfbe8321SBarry Smith   PetscErrorCode ierr;
43838baddfdSBarry Smith   PetscInt       M,N;
439b8f8c88eSHong Zhang   PetscMPIInt    size;
440639f9d9dSBarry Smith 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
442d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
443639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
44429bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
445639f9d9dSBarry Smith 
446f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
44752e6d16bSBarry Smith   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
448639f9d9dSBarry Smith 
449b8f8c88eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
450b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
451b8f8c88eSHong Zhang 
452f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
453f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
454639f9d9dSBarry Smith   } else {
45529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
456639f9d9dSBarry Smith   }
457639f9d9dSBarry Smith 
458b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
459b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
460b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
461b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
462b8f8c88eSHong Zhang 
46377d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
46477d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
465005c665bSBarry Smith   c->freq              = 1;
46640964ad5SBarry Smith   c->usersetsrecompute = PETSC_FALSE;
46740964ad5SBarry Smith   c->recompute         = PETSC_FALSE;
46849b058dcSBarry Smith   c->currentcolor      = -1;
469efb30889SBarry Smith   c->htype             = "wp";
470639f9d9dSBarry Smith 
471639f9d9dSBarry Smith   *color = c;
472d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
474bbf0e169SBarry Smith }
475bbf0e169SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
47805869f15SSatish Balay /*@
479639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
480639f9d9dSBarry Smith     via MatFDColoringCreate().
481bbf0e169SBarry Smith 
482ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
483ef5ee4d1SLois Curfman McInnes 
484b4fc646aSLois Curfman McInnes     Input Parameter:
485639f9d9dSBarry Smith .   c - coloring context
486bbf0e169SBarry Smith 
48715091d37SBarry Smith     Level: intermediate
48815091d37SBarry Smith 
489639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
490bbf0e169SBarry Smith @*/
491be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c)
492bbf0e169SBarry Smith {
4936849ba73SBarry Smith   PetscErrorCode ierr;
49438baddfdSBarry Smith   PetscInt       i;
495bbf0e169SBarry Smith 
4963a40ed3dSBarry Smith   PetscFunctionBegin;
4977adad957SLisandro Dalcin   if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0);
498d4bb536fSBarry Smith 
499639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
50005b42c5fSBarry Smith     ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);
50105b42c5fSBarry Smith     ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);
50205b42c5fSBarry Smith     ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);
50305b42c5fSBarry Smith     if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
504bbf0e169SBarry Smith   }
505606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
506606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
507606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
508606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
509606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
51005b42c5fSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
5114f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
512005c665bSBarry Smith   if (c->w1) {
513005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
514005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
515b8f8c88eSHong Zhang   }
516b8f8c88eSHong Zhang   if (c->w3){
517005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
518005c665bSBarry Smith   }
519d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5203a40ed3dSBarry Smith   PetscFunctionReturn(0);
521bbf0e169SBarry Smith }
52243a90d84SBarry Smith 
5234a2ae208SSatish Balay #undef __FUNCT__
52449b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
52549b058dcSBarry Smith /*@C
52649b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52749b058dcSBarry Smith       that are currently being perturbed.
52849b058dcSBarry Smith 
52949b058dcSBarry Smith     Not Collective
53049b058dcSBarry Smith 
53149b058dcSBarry Smith     Input Parameters:
53249b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
53349b058dcSBarry Smith 
53449b058dcSBarry Smith     Output Parameters:
53549b058dcSBarry Smith +   n - the number of local columns being perturbed
53649b058dcSBarry Smith -   cols - the column indices, in global numbering
53749b058dcSBarry Smith 
53849b058dcSBarry Smith    Level: intermediate
53949b058dcSBarry Smith 
54049b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
54149b058dcSBarry Smith 
54249b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
54349b058dcSBarry Smith @*/
544be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
54549b058dcSBarry Smith {
54649b058dcSBarry Smith   PetscFunctionBegin;
54749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
54849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
54949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55049b058dcSBarry Smith   } else {
55149b058dcSBarry Smith     *n = 0;
55249b058dcSBarry Smith   }
55349b058dcSBarry Smith   PetscFunctionReturn(0);
55449b058dcSBarry Smith }
55549b058dcSBarry Smith 
556b8f8c88eSHong Zhang #include "petscda.h"      /*I      "petscda.h"    I*/
557b8f8c88eSHong Zhang 
55849b058dcSBarry Smith #undef __FUNCT__
5594a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
56043a90d84SBarry Smith /*@
561e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
562e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
56343a90d84SBarry Smith 
564fee21e36SBarry Smith     Collective on MatFDColoring
565fee21e36SBarry Smith 
566ef5ee4d1SLois Curfman McInnes     Input Parameters:
567ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
568ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
569ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5707850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
571ef5ee4d1SLois Curfman McInnes 
5728bba8e72SBarry Smith     Options Database Keys:
573b9722508SLois Curfman McInnes +    -mat_fd_coloring_freq <freq> - Sets coloring frequency
5743ec795f1SBarry Smith .    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
575b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
576b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
577b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
5788bba8e72SBarry Smith 
57915091d37SBarry Smith     Level: intermediate
58015091d37SBarry Smith 
5817850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
58243a90d84SBarry Smith 
58343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
58443a90d84SBarry Smith @*/
585be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
58643a90d84SBarry Smith {
5876849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5886849ba73SBarry Smith   PetscErrorCode ierr;
589b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
590efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
591ea709b57SSatish Balay   PetscScalar    *vscale_array;
592efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
593b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
594005c665bSBarry Smith   void           *fctx = coloring->fctx;
595f1af5d2fSBarry Smith   PetscTruth     flg;
596705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
597b8f8c88eSHong Zhang   Vec            x1_tmp;
598a2e34c3dSBarry Smith 
5993a40ed3dSBarry Smith   PetscFunctionBegin;
6004482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
6014482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
6024482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,3);
603feba9168SBarry Smith   if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
604e0907662SLois Curfman McInnes 
60540964ad5SBarry Smith   if (coloring->usersetsrecompute) {
60640964ad5SBarry Smith     if (!coloring->recompute) {
60740964ad5SBarry Smith       *flag = SAME_PRECONDITIONER;
608ae15b995SBarry Smith       ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr);
60940964ad5SBarry Smith       PetscFunctionReturn(0);
61040964ad5SBarry Smith     } else {
61140964ad5SBarry Smith       coloring->recompute = PETSC_FALSE;
61240964ad5SBarry Smith     }
61340964ad5SBarry Smith   }
61440964ad5SBarry Smith 
615d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6164c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
617b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
618f1af5d2fSBarry Smith   if (flg) {
619ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
620e0907662SLois Curfman McInnes   } else {
6210b9b6f31SBarry Smith     PetscTruth assembled;
6220b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
6230b9b6f31SBarry Smith     if (assembled) {
62443a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
625e0907662SLois Curfman McInnes     }
6260b9b6f31SBarry Smith   }
62743a90d84SBarry Smith 
628b8f8c88eSHong Zhang   x1_tmp = x1;
629dac7f5fdSBarry Smith   if (!coloring->vscale){
630b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
631b8f8c88eSHong Zhang   }
632be2fbe1fSBarry Smith 
6333a7fca6bSBarry Smith   /*
6343a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
6353a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
6363a7fca6bSBarry Smith   */
6373a7fca6bSBarry Smith   if (coloring->F) {
6383a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
6393a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
6403a7fca6bSBarry Smith     if (m1 != m2) {
6413a7fca6bSBarry Smith       coloring->F = 0;
6423a7fca6bSBarry Smith       }
6433a7fca6bSBarry Smith     }
6443a7fca6bSBarry Smith 
645b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
646b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
647b8f8c88eSHong Zhang   }
648b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
649b8f8c88eSHong Zhang 
650b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
651be2fbe1fSBarry Smith   if (coloring->F) {
652be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
653be2fbe1fSBarry Smith     coloring->F = 0;
654be2fbe1fSBarry Smith   } else {
65566f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
656b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
65766f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
658be2fbe1fSBarry Smith   }
65943a90d84SBarry Smith 
660b8f8c88eSHong Zhang   if (!coloring->w3) {
661b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
662b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
663efb30889SBarry Smith   }
664b8f8c88eSHong Zhang   w3 = coloring->w3;
665efb30889SBarry Smith 
666b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
667b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
668b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
669b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
670b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
671b8f8c88eSHong Zhang     col_start = 0; col_end = N;
6728ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
673b8f8c88eSHong Zhang     xx = xx - start;
674b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
675b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
676b8f8c88eSHong Zhang   }
677b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
678b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
679efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
680efb30889SBarry Smith       dx = 1.0 + unorm;
681efb30889SBarry Smith     } else {
6829782cf98SBarry Smith       dx  = xx[col];
683efb30889SBarry Smith     }
6845b8514ebSBarry Smith     if (dx == 0.0) dx = 1.0;
6859782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6869782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6879782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6889782cf98SBarry Smith #else
6899782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6909782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6919782cf98SBarry Smith #endif
6929782cf98SBarry Smith     dx               *= epsilon;
69330b34957SBarry Smith     vscale_array[col] = 1.0/dx;
6949782cf98SBarry Smith   }
6958ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
696efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6978ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
69830b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69930b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
700b8f8c88eSHong Zhang   }
7019782cf98SBarry Smith 
702b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
703b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
704b8f8c88eSHong Zhang   } else {
705b8f8c88eSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
706b8f8c88eSHong Zhang   }
707b0a32e0cSBarry Smith 
7089782cf98SBarry Smith   /*
70943a90d84SBarry Smith     Loop over each color
71043a90d84SBarry Smith   */
711b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
71243a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
71349b058dcSBarry Smith     coloring->currentcolor = k;
714b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
715b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
7168ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
71743a90d84SBarry Smith     /*
718b8f8c88eSHong Zhang       Loop over each column associated with color
719b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
72043a90d84SBarry Smith     */
72143a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
722b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
723efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
724efb30889SBarry Smith         dx = 1.0 + unorm;
725efb30889SBarry Smith       } else {
72642460c72SBarry Smith         dx  = xx[col];
727efb30889SBarry Smith       }
7285b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
729aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
730ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
731ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
73243a90d84SBarry Smith #else
733329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
734329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
73543a90d84SBarry Smith #endif
73643a90d84SBarry Smith       dx            *= epsilon;
7371302d50aSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
73842460c72SBarry Smith       w3_array[col] += dx;
73943a90d84SBarry Smith     }
7408ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
741b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7423b28642cSBarry Smith 
74343a90d84SBarry Smith     /*
744b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
745b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
74643a90d84SBarry Smith     */
74766f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
74843a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
74966f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
750efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
7519782cf98SBarry Smith 
75243a90d84SBarry Smith     /*
753e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
75443a90d84SBarry Smith     */
7553b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
75643a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
757b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
758b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
7595904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
76043a90d84SBarry Smith       srow   = row + start;
76143a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
76243a90d84SBarry Smith     }
7633b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
764aeb7e98aSMatthew Knepley   } /* endof for each color */
7658ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
76630b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
767b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
768b8f8c88eSHong Zhang 
769b8f8c88eSHong Zhang   coloring->currentcolor = -1;
77043a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77143a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
772d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
773a2e34c3dSBarry Smith 
774b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
775a2e34c3dSBarry Smith   if (flg) {
776a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
777a2e34c3dSBarry Smith   }
778b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7793a40ed3dSBarry Smith   PetscFunctionReturn(0);
78043a90d84SBarry Smith }
781840b8ebdSBarry Smith 
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS"
784840b8ebdSBarry Smith /*@
785840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
786840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
787840b8ebdSBarry Smith 
788fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
789fee21e36SBarry Smith 
790ef5ee4d1SLois Curfman McInnes     Input Parameters:
7913b28642cSBarry Smith +   mat - location to store Jacobian
792ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
793ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
7947850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null
795ef5ee4d1SLois Curfman McInnes 
796840b8ebdSBarry Smith    Options Database Keys:
797ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
798840b8ebdSBarry Smith 
79915091d37SBarry Smith    Level: intermediate
80015091d37SBarry Smith 
8017850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
802840b8ebdSBarry Smith 
803840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
804840b8ebdSBarry Smith @*/
805be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
806840b8ebdSBarry Smith {
8076849ba73SBarry Smith   PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f;
8086849ba73SBarry Smith   PetscErrorCode ierr;
80938baddfdSBarry Smith   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow;
810efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
811ea709b57SSatish Balay   PetscScalar    *vscale_array;
812329f5518SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
813ced766e8SHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
814840b8ebdSBarry Smith   void           *fctx = coloring->fctx;
815f1af5d2fSBarry Smith   PetscTruth     flg;
816840b8ebdSBarry Smith 
8173a40ed3dSBarry Smith   PetscFunctionBegin;
8184482741eSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
8194482741eSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2);
8204482741eSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE,4);
821840b8ebdSBarry Smith 
822d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
823ced766e8SHong Zhang   if (!coloring->w3) {
824840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
82552e6d16bSBarry Smith     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
826840b8ebdSBarry Smith   }
827ced766e8SHong Zhang   w3 = coloring->w3;
828840b8ebdSBarry Smith 
8295904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
830b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
831f1af5d2fSBarry Smith   if (flg) {
832ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
833840b8ebdSBarry Smith   } else {
8340b9b6f31SBarry Smith     PetscTruth assembled;
8350b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
8360b9b6f31SBarry Smith     if (assembled) {
837840b8ebdSBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
838840b8ebdSBarry Smith     }
8390b9b6f31SBarry Smith   }
840840b8ebdSBarry Smith 
841840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
842840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
84366f9b7ceSBarry Smith   ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
844840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
84566f9b7ceSBarry Smith   ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
846840b8ebdSBarry Smith 
847840b8ebdSBarry Smith   /*
8485904e1b1SBarry Smith       Compute all the scale factors and share with other processors
849840b8ebdSBarry Smith   */
8505904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
8515904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
852840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
853840b8ebdSBarry Smith     /*
854840b8ebdSBarry Smith        Loop over each column associated with color adding the
855840b8ebdSBarry Smith        perturbation to the vector w3.
856840b8ebdSBarry Smith     */
857840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
858840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8595904e1b1SBarry Smith       dx  = xx[col];
8605b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
861aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
862840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
863840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
864840b8ebdSBarry Smith #else
865329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
866329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
867840b8ebdSBarry Smith #endif
868840b8ebdSBarry Smith       dx                *= epsilon;
8695904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
870840b8ebdSBarry Smith     }
8715904e1b1SBarry Smith   }
8725904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8735904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8745904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8755904e1b1SBarry Smith 
8765904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
8775904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
8785904e1b1SBarry Smith 
8795904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
8805904e1b1SBarry Smith   /*
8815904e1b1SBarry Smith       Loop over each color
8825904e1b1SBarry Smith   */
8835904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
8845904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
8855904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
8865904e1b1SBarry Smith     /*
8875904e1b1SBarry Smith        Loop over each column associated with color adding the
8885904e1b1SBarry Smith        perturbation to the vector w3.
8895904e1b1SBarry Smith     */
8905904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
8915904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
8925904e1b1SBarry Smith       dx  = xx[col];
8935b8514ebSBarry Smith       if (dx == 0.0) dx = 1.0;
8945904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8955904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
8965904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
8975904e1b1SBarry Smith #else
8985904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
8995904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
9005904e1b1SBarry Smith #endif
9015904e1b1SBarry Smith       dx            *= epsilon;
9025904e1b1SBarry Smith       w3_array[col] += dx;
9035904e1b1SBarry Smith     }
9045904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
9055904e1b1SBarry Smith 
906840b8ebdSBarry Smith     /*
907840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
908840b8ebdSBarry Smith     */
90966f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
910840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
91166f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
912efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
9135904e1b1SBarry Smith 
914840b8ebdSBarry Smith     /*
915840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
916840b8ebdSBarry Smith     */
9173b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
918840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
919840b8ebdSBarry Smith       row    = coloring->rows[k][l];
920840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
9215904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
922840b8ebdSBarry Smith       srow   = row + start;
923840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
924840b8ebdSBarry Smith     }
9253b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
926840b8ebdSBarry Smith   }
9275904e1b1SBarry Smith   ierr  = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
9285904e1b1SBarry Smith   xx    = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
929840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
930840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
931d5ba7fb7SMatthew Knepley   ierr  = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
9323a40ed3dSBarry Smith   PetscFunctionReturn(0);
933840b8ebdSBarry Smith }
9343b28642cSBarry Smith 
9353b28642cSBarry Smith 
9364a2ae208SSatish Balay #undef __FUNCT__
9374a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()"
93840964ad5SBarry Smith /*@C
93940964ad5SBarry Smith    MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner
94040964ad5SBarry Smith      is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed
94140964ad5SBarry Smith      no additional Jacobian's will be computed (the same one will be used) until this is
94240964ad5SBarry Smith      called again.
94340964ad5SBarry Smith 
94440964ad5SBarry Smith    Collective on MatFDColoring
94540964ad5SBarry Smith 
94640964ad5SBarry Smith    Input  Parameters:
94740964ad5SBarry Smith .  c - the coloring context
94840964ad5SBarry Smith 
94940964ad5SBarry Smith    Level: intermediate
95040964ad5SBarry Smith 
95140964ad5SBarry Smith    Notes: The MatFDColoringSetFrequency() is ignored once this is called
95240964ad5SBarry Smith 
95340964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency()
95440964ad5SBarry Smith 
95540964ad5SBarry Smith .keywords: Mat, finite differences, coloring
95640964ad5SBarry Smith @*/
957be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c)
95840964ad5SBarry Smith {
95940964ad5SBarry Smith   PetscFunctionBegin;
9604482741eSBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1);
96140964ad5SBarry Smith   c->usersetsrecompute = PETSC_TRUE;
96240964ad5SBarry Smith   c->recompute         = PETSC_TRUE;
96340964ad5SBarry Smith   PetscFunctionReturn(0);
96440964ad5SBarry Smith }
96540964ad5SBarry Smith 
9663b28642cSBarry Smith 
967