xref: /petsc/src/mat/matfd/fdmatrix.c (revision e7e920445bafd8671dd3ea579405cd582ea97ef2)
1be1d678aSKris Buschelman 
2bbf0e169SBarry Smith /*
3639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
4639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
5bbf0e169SBarry Smith */
6bbf0e169SBarry Smith 
7af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
8af0996ceSBarry Smith #include <petsc/private/isimpl.h>
9bbf0e169SBarry Smith 
107087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
113a7fca6bSBarry Smith {
124e269d77SPeter Brune   PetscErrorCode ierr;
136e111a19SKarl Rupp 
143a7fca6bSBarry Smith   PetscFunctionBegin;
154e269d77SPeter Brune   if (F) {
164e269d77SPeter Brune     ierr     = VecCopy(F,fd->w1);CHKERRQ(ierr);
174e269d77SPeter Brune     fd->fset = PETSC_TRUE;
184e269d77SPeter Brune   } else {
194e269d77SPeter Brune     fd->fset = PETSC_FALSE;
204e269d77SPeter Brune   }
213a7fca6bSBarry Smith   PetscFunctionReturn(0);
223a7fca6bSBarry Smith }
233a7fca6bSBarry Smith 
249804daf3SBarry Smith #include <petscdraw.h>
256849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
269194eea9SBarry Smith {
279194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
28dfbe8321SBarry Smith   PetscErrorCode ierr;
29b312708bSHong Zhang   PetscInt       i,j,nz,row;
309194eea9SBarry Smith   PetscReal      x,y;
31b312708bSHong Zhang   MatEntry       *Jentry=fd->matentry;
329194eea9SBarry Smith 
339194eea9SBarry Smith   PetscFunctionBegin;
349194eea9SBarry Smith   /* loop over colors  */
35b312708bSHong Zhang   nz = 0;
369194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
379194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
38b312708bSHong Zhang       row = Jentry[nz].row;
39b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
40b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
41b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
429194eea9SBarry Smith     }
439194eea9SBarry Smith   }
449194eea9SBarry Smith   PetscFunctionReturn(0);
459194eea9SBarry Smith }
469194eea9SBarry Smith 
476849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
48005c665bSBarry Smith {
49dfbe8321SBarry Smith   PetscErrorCode ierr;
50ace3abfcSBarry Smith   PetscBool      isnull;
51b0a32e0cSBarry Smith   PetscDraw      draw;
5254d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
53005c665bSBarry Smith 
543a40ed3dSBarry Smith   PetscFunctionBegin;
55b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
56832b7cebSLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
57832b7cebSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
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);
62832b7cebSLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
63b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
64f77a5242SKarl Rupp   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr);
65832b7cebSLisandro Dalcin   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
663a40ed3dSBarry Smith   PetscFunctionReturn(0);
67005c665bSBarry Smith }
68005c665bSBarry Smith 
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 @*/
987087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
99bbf0e169SBarry Smith {
1006849ba73SBarry Smith   PetscErrorCode    ierr;
10138baddfdSBarry Smith   PetscInt          i,j;
102ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
103f3ef73ceSBarry Smith   PetscViewerFormat format;
104bbf0e169SBarry Smith 
1053a40ed3dSBarry Smith   PetscFunctionBegin;
1060700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1073050cee2SBarry Smith   if (!viewer) {
108ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
1093050cee2SBarry Smith   }
1100700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
111c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
112bbf0e169SBarry Smith 
113251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
114251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1150f5bd95cSBarry Smith   if (isdraw) {
116b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11732077d6dSBarry Smith   } else if (iascii) {
118dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
11957622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr);
12057622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",(double)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) {
125b312708bSHong Zhang       PetscInt row,col,nz;
126b312708bSHong Zhang       nz = 0;
127b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13177431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
132639f9d9dSBarry Smith         }
13377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
1345bdb020cSBarry Smith         if (c->matentry) {
135b4fc646aSLois Curfman McInnes           for (j=0; j<c->nrows[i]; j++) {
136b312708bSHong Zhang             row  = c->matentry[nz].row;
137b312708bSHong Zhang             col  = c->matentry[nz++].col;
138b312708bSHong Zhang             ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);CHKERRQ(ierr);
139b4fc646aSLois Curfman McInnes           }
140bbf0e169SBarry Smith         }
141bbf0e169SBarry Smith       }
1425bdb020cSBarry Smith     }
143b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
144bbf0e169SBarry Smith   }
1453a40ed3dSBarry Smith   PetscFunctionReturn(0);
146639f9d9dSBarry Smith }
147639f9d9dSBarry Smith 
148639f9d9dSBarry Smith /*@
149b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
150b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
151639f9d9dSBarry Smith 
1523f9fe445SBarry Smith    Logically Collective on MatFDColoring
153ef5ee4d1SLois Curfman McInnes 
154ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
155ef5ee4d1SLois Curfman McInnes .vb
15665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
157d461c19dSHong Zhang        htype = 'ds':
158f23b5b22SLois Curfman McInnes          h = error_rel*u[i]                 if  abs(u[i]) > umin
159f23b5b22SLois Curfman McInnes            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
160ef5ee4d1SLois Curfman McInnes          dx_{i} = (0, ... 1, .... 0)
161d461c19dSHong Zhang 
162d461c19dSHong Zhang        htype = 'wp':
163d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
164ef5ee4d1SLois Curfman McInnes .ve
165639f9d9dSBarry Smith 
166639f9d9dSBarry Smith    Input Parameters:
167ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
168639f9d9dSBarry Smith .  error_rel - relative error
169f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
170fee21e36SBarry Smith 
17115091d37SBarry Smith    Level: advanced
17215091d37SBarry Smith 
173b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
174b4fc646aSLois Curfman McInnes 
17595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17695b89fc3SBarry Smith 
177639f9d9dSBarry Smith @*/
1787087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
179639f9d9dSBarry Smith {
1803a40ed3dSBarry Smith   PetscFunctionBegin;
1810700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
182c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
183c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
184639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
185639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1863a40ed3dSBarry Smith   PetscFunctionReturn(0);
187639f9d9dSBarry Smith }
188639f9d9dSBarry Smith 
189f86b9fbaSHong Zhang /*@
190c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
191005c665bSBarry Smith 
192f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
193f86b9fbaSHong Zhang 
194f86b9fbaSHong Zhang    Input Parameters:
195f86b9fbaSHong Zhang +  coloring - the coloring context
196f86b9fbaSHong Zhang .  brows - number of rows in the block
197f86b9fbaSHong Zhang -  bcols - number of columns in the block
198f86b9fbaSHong Zhang 
199f86b9fbaSHong Zhang    Level: intermediate
200f86b9fbaSHong Zhang 
201f86b9fbaSHong Zhang .keywords: Mat, coloring
202f86b9fbaSHong Zhang 
203f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
204f86b9fbaSHong Zhang 
205f86b9fbaSHong Zhang @*/
206f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
207f86b9fbaSHong Zhang {
208f86b9fbaSHong Zhang   PetscFunctionBegin;
209f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
210f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,brows,2);
211f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,bcols,3);
212f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
213f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
214f86b9fbaSHong Zhang   PetscFunctionReturn(0);
215f86b9fbaSHong Zhang }
216f86b9fbaSHong Zhang 
217f86b9fbaSHong Zhang /*@
218f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
219f86b9fbaSHong Zhang 
220f86b9fbaSHong Zhang    Collective on Mat
221f86b9fbaSHong Zhang 
222f86b9fbaSHong Zhang    Input Parameters:
223f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
224f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
225f86b9fbaSHong Zhang -  color - the matrix coloring context
226f86b9fbaSHong Zhang 
227f86b9fbaSHong Zhang    Level: beginner
228f86b9fbaSHong Zhang 
229f86b9fbaSHong Zhang .keywords: MatFDColoring, setup
230f86b9fbaSHong Zhang 
231f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
232f86b9fbaSHong Zhang @*/
233f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
234f86b9fbaSHong Zhang {
235f86b9fbaSHong Zhang   PetscErrorCode ierr;
236f86b9fbaSHong Zhang 
237f86b9fbaSHong Zhang   PetscFunctionBegin;
238c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
239c8a9c622SHong Zhang   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
240c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
241c8a9c622SHong Zhang 
2420df34763SHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
243f86b9fbaSHong Zhang   if (mat->ops->fdcoloringsetup) {
244f86b9fbaSHong Zhang     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
245f86b9fbaSHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
246c8a9c622SHong Zhang 
247c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2480df34763SHong Zhang    ierr   = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
249f86b9fbaSHong Zhang   PetscFunctionReturn(0);
250f86b9fbaSHong Zhang }
251ff0cfa39SBarry Smith 
2524a9d489dSBarry Smith /*@C
2534a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2544a9d489dSBarry Smith 
2553f9fe445SBarry Smith    Not Collective
2564a9d489dSBarry Smith 
2574a9d489dSBarry Smith    Input Parameters:
2584a9d489dSBarry Smith .  coloring - the coloring context
2594a9d489dSBarry Smith 
2604a9d489dSBarry Smith    Output Parameters:
2614a9d489dSBarry Smith +  f - the function
2624a9d489dSBarry Smith -  fctx - the optional user-defined function context
2634a9d489dSBarry Smith 
2644a9d489dSBarry Smith    Level: intermediate
2654a9d489dSBarry Smith 
2664a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
26795b89fc3SBarry Smith 
26895b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
26995b89fc3SBarry Smith 
2704a9d489dSBarry Smith @*/
2717087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2724a9d489dSBarry Smith {
2734a9d489dSBarry Smith   PetscFunctionBegin;
2740700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2754a9d489dSBarry Smith   if (f) *f = matfd->f;
2764a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2774a9d489dSBarry Smith   PetscFunctionReturn(0);
2784a9d489dSBarry Smith }
2794a9d489dSBarry Smith 
280d64ed03dSBarry Smith /*@C
281005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
282005c665bSBarry Smith 
2833f9fe445SBarry Smith    Logically Collective on MatFDColoring
284fee21e36SBarry Smith 
285ef5ee4d1SLois Curfman McInnes    Input Parameters:
286ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
287ef5ee4d1SLois Curfman McInnes .  f - the function
288ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
289ef5ee4d1SLois Curfman McInnes 
2907850c7c0SBarry Smith    Calling sequence of (*f) function:
2917850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
292ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
29315091d37SBarry Smith 
2947850c7c0SBarry Smith    Level: advanced
2957850c7c0SBarry Smith 
29695452b02SPatrick Sanan    Notes:
29795452b02SPatrick Sanan     This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2988d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
299ab637aeaSJed Brown      calling MatFDColoringApply()
3007850c7c0SBarry Smith 
3017850c7c0SBarry Smith    Fortran Notes:
3027850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
303ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
304f881d145SBarry Smith 
305b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
30695b89fc3SBarry Smith 
30795b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
30895b89fc3SBarry Smith 
309005c665bSBarry Smith @*/
3107087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
311005c665bSBarry Smith {
3123a40ed3dSBarry Smith   PetscFunctionBegin;
3130700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
314005c665bSBarry Smith   matfd->f    = f;
315005c665bSBarry Smith   matfd->fctx = fctx;
3163a40ed3dSBarry Smith   PetscFunctionReturn(0);
317005c665bSBarry Smith }
318005c665bSBarry Smith 
319639f9d9dSBarry Smith /*@
320b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
321639f9d9dSBarry Smith    the options database.
322639f9d9dSBarry Smith 
323fee21e36SBarry Smith    Collective on MatFDColoring
324fee21e36SBarry Smith 
32565f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
326ef5ee4d1SLois Curfman McInnes .vb
32765f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
328f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
329f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
330ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
331ef5ee4d1SLois Curfman McInnes .ve
332ef5ee4d1SLois Curfman McInnes 
333ef5ee4d1SLois Curfman McInnes    Input Parameter:
334ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
335ef5ee4d1SLois Curfman McInnes 
336b4fc646aSLois Curfman McInnes    Options Database Keys:
3375620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
338f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3393ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
340ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
341e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
342e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
343639f9d9dSBarry Smith 
34415091d37SBarry Smith     Level: intermediate
34515091d37SBarry Smith 
346b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
347d1c39f55SBarry Smith 
348d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
349d1c39f55SBarry Smith 
350639f9d9dSBarry Smith @*/
3517087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
352639f9d9dSBarry Smith {
353dfbe8321SBarry Smith   PetscErrorCode ierr;
354ace3abfcSBarry Smith   PetscBool      flg;
355efb30889SBarry Smith   char           value[3];
3563a40ed3dSBarry Smith 
3573a40ed3dSBarry Smith   PetscFunctionBegin;
3580700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
359639f9d9dSBarry Smith 
3603194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
36187828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
36287828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3631bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
364efb30889SBarry Smith   if (flg) {
365efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
366efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
367e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
368efb30889SBarry Smith   }
369f86b9fbaSHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
37093dfae19SHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
37193dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
37293dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
37393dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
37493dfae19SHong Zhang   }
375f86b9fbaSHong Zhang 
3765d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3770633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr);
378b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
380005c665bSBarry Smith }
381005c665bSBarry Smith 
38261ab5df0SBarry Smith /*@C
38361ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
38461ab5df0SBarry Smith 
38561ab5df0SBarry Smith    Collective on MatFDColoring
38661ab5df0SBarry Smith 
38761ab5df0SBarry Smith    Input Parameters:
38861ab5df0SBarry Smith +  coloring - the coloring context
38961ab5df0SBarry Smith -  type - either MATMFFD_WP or MATMFFD_DS
39061ab5df0SBarry Smith 
39161ab5df0SBarry Smith    Options Database Keys:
39261ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
39361ab5df0SBarry Smith 
39461ab5df0SBarry Smith    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
39561ab5df0SBarry Smith          but the process of computing the entries is the same as as with the MatMFFD operation so we should reuse the names instead of
39661ab5df0SBarry Smith          introducing another one.
39761ab5df0SBarry Smith 
39861ab5df0SBarry Smith    Level: intermediate
39961ab5df0SBarry Smith 
40061ab5df0SBarry Smith .keywords: Mat, finite differences, parameters
40161ab5df0SBarry Smith 
40261ab5df0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
40361ab5df0SBarry Smith 
40461ab5df0SBarry Smith @*/
40561ab5df0SBarry Smith PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
40661ab5df0SBarry Smith {
40761ab5df0SBarry Smith   PetscFunctionBegin;
40861ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
40961ab5df0SBarry Smith   /*
41061ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
41161ab5df0SBarry Smith      and this function is being provided as patch for a release so it shouldn't change the implementaton
41261ab5df0SBarry Smith   */
41361ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
41461ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
41561ab5df0SBarry Smith   else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
41661ab5df0SBarry Smith   PetscFunctionReturn(0);
41761ab5df0SBarry Smith }
41861ab5df0SBarry Smith 
419146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
420005c665bSBarry Smith {
421dfbe8321SBarry Smith   PetscErrorCode    ierr;
422e350db43SBarry Smith   PetscBool         flg;
4233050cee2SBarry Smith   PetscViewer       viewer;
424cffb1e40SBarry Smith   PetscViewerFormat format;
425005c665bSBarry Smith 
4263a40ed3dSBarry Smith   PetscFunctionBegin;
427146574abSBarry Smith   if (prefix) {
42816413a6aSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
429146574abSBarry Smith   } else {
43016413a6aSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
431146574abSBarry Smith   }
432005c665bSBarry Smith   if (flg) {
433cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
4343050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
435cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
436cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
437005c665bSBarry Smith   }
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
439bbf0e169SBarry Smith }
440bbf0e169SBarry Smith 
44105869f15SSatish Balay /*@
442639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
443639f9d9dSBarry Smith    computation of Jacobians.
444bbf0e169SBarry Smith 
445ef5ee4d1SLois Curfman McInnes    Collective on Mat
446ef5ee4d1SLois Curfman McInnes 
447639f9d9dSBarry Smith    Input Parameters:
448ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
4497086a01eSPeter Brune -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
450bbf0e169SBarry Smith 
451bbf0e169SBarry Smith     Output Parameter:
452639f9d9dSBarry Smith .   color - the new coloring context
453bbf0e169SBarry Smith 
45415091d37SBarry Smith     Level: intermediate
45515091d37SBarry Smith 
4568d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
457d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
4587086a01eSPeter Brune           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
459bbf0e169SBarry Smith @*/
4607087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
461bbf0e169SBarry Smith {
462639f9d9dSBarry Smith   MatFDColoring  c;
463639f9d9dSBarry Smith   MPI_Comm       comm;
464dfbe8321SBarry Smith   PetscErrorCode ierr;
46538baddfdSBarry Smith   PetscInt       M,N;
466639f9d9dSBarry Smith 
4673a40ed3dSBarry Smith   PetscFunctionBegin;
468c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
469f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
470d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
471639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
472ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
473f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
47473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
475639f9d9dSBarry Smith 
476b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
477b8f8c88eSHong Zhang 
47893dfae19SHong Zhang   if (mat->ops->fdcoloringcreate) {
47993dfae19SHong Zhang     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
48093dfae19SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
48193dfae19SHong Zhang 
4822a7a6963SBarry Smith   ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
483*e7e92044SBarry Smith   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
484*e7e92044SBarry Smith   ierr = VecPinToCPU(c->w1,PETSC_TRUE);CHKERRQ(ierr);
4853bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
486b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
487*e7e92044SBarry Smith   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
488*e7e92044SBarry Smith   ierr = VecPinToCPU(c->w2,PETSC_TRUE);CHKERRQ(ierr);
4893bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
490b8f8c88eSHong Zhang 
49177d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
49277d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
49349b058dcSBarry Smith   c->currentcolor = -1;
494efb30889SBarry Smith   c->htype        = "wp";
4954e269d77SPeter Brune   c->fset         = PETSC_FALSE;
496c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
497639f9d9dSBarry Smith 
498639f9d9dSBarry Smith   *color = c;
4994e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
500d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
5013a40ed3dSBarry Smith   PetscFunctionReturn(0);
502bbf0e169SBarry Smith }
503bbf0e169SBarry Smith 
50405869f15SSatish Balay /*@
505639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
506639f9d9dSBarry Smith     via MatFDColoringCreate().
507bbf0e169SBarry Smith 
508ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
509ef5ee4d1SLois Curfman McInnes 
510b4fc646aSLois Curfman McInnes     Input Parameter:
511639f9d9dSBarry Smith .   c - coloring context
512bbf0e169SBarry Smith 
51315091d37SBarry Smith     Level: intermediate
51415091d37SBarry Smith 
515639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
516bbf0e169SBarry Smith @*/
5176bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
518bbf0e169SBarry Smith {
5196849ba73SBarry Smith   PetscErrorCode ierr;
52038baddfdSBarry Smith   PetscInt       i;
5210df34763SHong Zhang   MatFDColoring  color = *c;
522bbf0e169SBarry Smith 
5233a40ed3dSBarry Smith   PetscFunctionBegin;
5246bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
5250df34763SHong Zhang   if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);}
526d4bb536fSBarry Smith 
5270df34763SHong Zhang   for (i=0; i<color->ncolors; i++) {
5280df34763SHong Zhang     ierr = PetscFree(color->columns[i]);CHKERRQ(ierr);
529bbf0e169SBarry Smith   }
5300df34763SHong Zhang   ierr = PetscFree(color->ncolumns);CHKERRQ(ierr);
5310df34763SHong Zhang   ierr = PetscFree(color->columns);CHKERRQ(ierr);
5320df34763SHong Zhang   ierr = PetscFree(color->nrows);CHKERRQ(ierr);
5330df34763SHong Zhang   if (color->htype[0] == 'w') {
5340df34763SHong Zhang     ierr = PetscFree(color->matentry2);CHKERRQ(ierr);
5350df34763SHong Zhang   } else {
5360df34763SHong Zhang     ierr = PetscFree(color->matentry);CHKERRQ(ierr);
5370df34763SHong Zhang   }
5380df34763SHong Zhang   ierr = PetscFree(color->dy);CHKERRQ(ierr);
5390df34763SHong Zhang   if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);}
5400df34763SHong Zhang   ierr = VecDestroy(&color->w1);CHKERRQ(ierr);
5410df34763SHong Zhang   ierr = VecDestroy(&color->w2);CHKERRQ(ierr);
5420df34763SHong Zhang   ierr = VecDestroy(&color->w3);CHKERRQ(ierr);
543d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
545bbf0e169SBarry Smith }
54643a90d84SBarry Smith 
54749b058dcSBarry Smith /*@C
54849b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
54949b058dcSBarry Smith       that are currently being perturbed.
55049b058dcSBarry Smith 
55149b058dcSBarry Smith     Not Collective
55249b058dcSBarry Smith 
55349b058dcSBarry Smith     Input Parameters:
55449b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
55549b058dcSBarry Smith 
55649b058dcSBarry Smith     Output Parameters:
55749b058dcSBarry Smith +   n - the number of local columns being perturbed
55849b058dcSBarry Smith -   cols - the column indices, in global numbering
55949b058dcSBarry Smith 
56021fcc2ddSBarry Smith    Level: advanced
56149b058dcSBarry Smith 
562edaa7c33SBarry Smith    Fortran Note:
563edaa7c33SBarry Smith    This routine has a different interface for Fortran
56421fcc2ddSBarry Smith $     #include <petsc/finclude/petscmat.h>
56521fcc2ddSBarry Smith $          use petscmat
566edaa7c33SBarry Smith $          PetscInt, pointer :: array(:)
567edaa7c33SBarry Smith $          PetscErrorCode  ierr
568edaa7c33SBarry Smith $          MatFDColoring   i
569edaa7c33SBarry Smith $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
570edaa7c33SBarry Smith $      use the entries of array ...
571edaa7c33SBarry Smith $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
572edaa7c33SBarry Smith 
57349b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
57449b058dcSBarry Smith 
57549b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
57649b058dcSBarry Smith @*/
577edaa7c33SBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[])
57849b058dcSBarry Smith {
57949b058dcSBarry Smith   PetscFunctionBegin;
58049b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
58149b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
58249b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
58349b058dcSBarry Smith   } else {
58449b058dcSBarry Smith     *n = 0;
58549b058dcSBarry Smith   }
58649b058dcSBarry Smith   PetscFunctionReturn(0);
58749b058dcSBarry Smith }
58849b058dcSBarry Smith 
58943a90d84SBarry Smith /*@
590e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
591e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
59243a90d84SBarry Smith 
593fee21e36SBarry Smith     Collective on MatFDColoring
594fee21e36SBarry Smith 
595ef5ee4d1SLois Curfman McInnes     Input Parameters:
596ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
597ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
598ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5997850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
600ef5ee4d1SLois Curfman McInnes 
6018bba8e72SBarry Smith     Options Database Keys:
602ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
603b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
604e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
605e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
6068bba8e72SBarry Smith 
60715091d37SBarry Smith     Level: intermediate
60898d79826SSatish Balay 
6097850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
61043a90d84SBarry Smith 
61143a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
61243a90d84SBarry Smith @*/
613d1e9a80fSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
61443a90d84SBarry Smith {
6153acb8795SBarry Smith   PetscErrorCode ierr;
616684f2004SHong Zhang   PetscBool      flg = PETSC_FALSE;
6173acb8795SBarry Smith 
6183acb8795SBarry Smith   PetscFunctionBegin;
6190700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
6200700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
6210700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
622ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
623e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
624c8a9c622SHong Zhang   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
625684f2004SHong Zhang 
626684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
627c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
628684f2004SHong Zhang   if (flg) {
629684f2004SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
630684f2004SHong Zhang   } else {
631684f2004SHong Zhang     PetscBool assembled;
632684f2004SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
633684f2004SHong Zhang     if (assembled) {
634684f2004SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
635684f2004SHong Zhang     }
636684f2004SHong Zhang   }
637684f2004SHong Zhang 
6385922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
639d1e9a80fSBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr);
6405922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6415bdb020cSBarry Smith   if (!coloring->viewed) {
6425bdb020cSBarry Smith     ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr);
6435bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6445bdb020cSBarry Smith   }
6453acb8795SBarry Smith   PetscFunctionReturn(0);
6463acb8795SBarry Smith }
647