xref: /petsc/src/mat/matfd/fdmatrix.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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
157f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
158f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
159ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
160ef5ee4d1SLois Curfman McInnes .ve
161639f9d9dSBarry Smith 
162639f9d9dSBarry Smith    Input Parameters:
163ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
164639f9d9dSBarry Smith .  error_rel - relative error
165f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
166fee21e36SBarry Smith 
16715091d37SBarry Smith    Level: advanced
16815091d37SBarry Smith 
169b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
170b4fc646aSLois Curfman McInnes 
17195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17295b89fc3SBarry Smith 
173639f9d9dSBarry Smith @*/
1747087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
175639f9d9dSBarry Smith {
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1770700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
178c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
179c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
180639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
181639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
183639f9d9dSBarry Smith }
184639f9d9dSBarry Smith 
185f86b9fbaSHong Zhang /*@
186c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
187005c665bSBarry Smith 
188f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
189f86b9fbaSHong Zhang 
190f86b9fbaSHong Zhang    Input Parameters:
191f86b9fbaSHong Zhang +  coloring - the coloring context
192f86b9fbaSHong Zhang .  brows - number of rows in the block
193f86b9fbaSHong Zhang -  bcols - number of columns in the block
194f86b9fbaSHong Zhang 
195f86b9fbaSHong Zhang    Level: intermediate
196f86b9fbaSHong Zhang 
197f86b9fbaSHong Zhang .keywords: Mat, coloring
198f86b9fbaSHong Zhang 
199f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
200f86b9fbaSHong Zhang 
201f86b9fbaSHong Zhang @*/
202f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
203f86b9fbaSHong Zhang {
204f86b9fbaSHong Zhang   PetscFunctionBegin;
205f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
206f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,brows,2);
207f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,bcols,3);
208f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
209f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
210f86b9fbaSHong Zhang   PetscFunctionReturn(0);
211f86b9fbaSHong Zhang }
212f86b9fbaSHong Zhang 
213f86b9fbaSHong Zhang /*@
214f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
215f86b9fbaSHong Zhang 
216f86b9fbaSHong Zhang    Collective on Mat
217f86b9fbaSHong Zhang 
218f86b9fbaSHong Zhang    Input Parameters:
219f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
220f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
221f86b9fbaSHong Zhang -  color - the matrix coloring context
222f86b9fbaSHong Zhang 
223f86b9fbaSHong Zhang    Level: beginner
224f86b9fbaSHong Zhang 
225f86b9fbaSHong Zhang .keywords: MatFDColoring, setup
226f86b9fbaSHong Zhang 
227f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
228f86b9fbaSHong Zhang @*/
229f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
230f86b9fbaSHong Zhang {
231f86b9fbaSHong Zhang   PetscErrorCode ierr;
232f86b9fbaSHong Zhang 
233f86b9fbaSHong Zhang   PetscFunctionBegin;
234c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
235c8a9c622SHong Zhang   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
236c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
237c8a9c622SHong Zhang 
2380df34763SHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
239f86b9fbaSHong Zhang   if (mat->ops->fdcoloringsetup) {
240f86b9fbaSHong Zhang     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
241f86b9fbaSHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
242c8a9c622SHong Zhang 
243c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2440df34763SHong Zhang    ierr   = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
245f86b9fbaSHong Zhang   PetscFunctionReturn(0);
246f86b9fbaSHong Zhang }
247ff0cfa39SBarry Smith 
2484a9d489dSBarry Smith /*@C
2494a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2504a9d489dSBarry Smith 
2513f9fe445SBarry Smith    Not Collective
2524a9d489dSBarry Smith 
2534a9d489dSBarry Smith    Input Parameters:
2544a9d489dSBarry Smith .  coloring - the coloring context
2554a9d489dSBarry Smith 
2564a9d489dSBarry Smith    Output Parameters:
2574a9d489dSBarry Smith +  f - the function
2584a9d489dSBarry Smith -  fctx - the optional user-defined function context
2594a9d489dSBarry Smith 
2604a9d489dSBarry Smith    Level: intermediate
2614a9d489dSBarry Smith 
2624a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
26395b89fc3SBarry Smith 
26495b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
26595b89fc3SBarry Smith 
2664a9d489dSBarry Smith @*/
2677087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2684a9d489dSBarry Smith {
2694a9d489dSBarry Smith   PetscFunctionBegin;
2700700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2714a9d489dSBarry Smith   if (f) *f = matfd->f;
2724a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2734a9d489dSBarry Smith   PetscFunctionReturn(0);
2744a9d489dSBarry Smith }
2754a9d489dSBarry Smith 
276d64ed03dSBarry Smith /*@C
277005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
278005c665bSBarry Smith 
2793f9fe445SBarry Smith    Logically Collective on MatFDColoring
280fee21e36SBarry Smith 
281ef5ee4d1SLois Curfman McInnes    Input Parameters:
282ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
283ef5ee4d1SLois Curfman McInnes .  f - the function
284ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
285ef5ee4d1SLois Curfman McInnes 
2867850c7c0SBarry Smith    Calling sequence of (*f) function:
2877850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
288ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
28915091d37SBarry Smith 
2907850c7c0SBarry Smith    Level: advanced
2917850c7c0SBarry Smith 
292*95452b02SPatrick Sanan    Notes:
293*95452b02SPatrick Sanan     This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2948d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
295ab637aeaSJed Brown      calling MatFDColoringApply()
2967850c7c0SBarry Smith 
2977850c7c0SBarry Smith    Fortran Notes:
2987850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
299ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
300f881d145SBarry Smith 
301b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
30295b89fc3SBarry Smith 
30395b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
30495b89fc3SBarry Smith 
305005c665bSBarry Smith @*/
3067087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
307005c665bSBarry Smith {
3083a40ed3dSBarry Smith   PetscFunctionBegin;
3090700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
310005c665bSBarry Smith   matfd->f    = f;
311005c665bSBarry Smith   matfd->fctx = fctx;
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
313005c665bSBarry Smith }
314005c665bSBarry Smith 
315639f9d9dSBarry Smith /*@
316b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
317639f9d9dSBarry Smith    the options database.
318639f9d9dSBarry Smith 
319fee21e36SBarry Smith    Collective on MatFDColoring
320fee21e36SBarry Smith 
32165f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
322ef5ee4d1SLois Curfman McInnes .vb
32365f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
324f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
325f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
326ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
327ef5ee4d1SLois Curfman McInnes .ve
328ef5ee4d1SLois Curfman McInnes 
329ef5ee4d1SLois Curfman McInnes    Input Parameter:
330ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
331ef5ee4d1SLois Curfman McInnes 
332b4fc646aSLois Curfman McInnes    Options Database Keys:
3335620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
334f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3353ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
336ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
337e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
338e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
339639f9d9dSBarry Smith 
34015091d37SBarry Smith     Level: intermediate
34115091d37SBarry Smith 
342b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
343d1c39f55SBarry Smith 
344d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
345d1c39f55SBarry Smith 
346639f9d9dSBarry Smith @*/
3477087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
348639f9d9dSBarry Smith {
349dfbe8321SBarry Smith   PetscErrorCode ierr;
350ace3abfcSBarry Smith   PetscBool      flg;
351efb30889SBarry Smith   char           value[3];
3523a40ed3dSBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
3540700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
355639f9d9dSBarry Smith 
3563194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
35787828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
35887828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3591bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
360efb30889SBarry Smith   if (flg) {
361efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
362efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
363e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
364efb30889SBarry Smith   }
365f86b9fbaSHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
36693dfae19SHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
36793dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
36893dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
36993dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
37093dfae19SHong Zhang   }
371f86b9fbaSHong Zhang 
3725d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3730633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr);
374b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3753a40ed3dSBarry Smith   PetscFunctionReturn(0);
376005c665bSBarry Smith }
377005c665bSBarry Smith 
37861ab5df0SBarry Smith /*@C
37961ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
38061ab5df0SBarry Smith 
38161ab5df0SBarry Smith    Collective on MatFDColoring
38261ab5df0SBarry Smith 
38361ab5df0SBarry Smith    Input Parameters:
38461ab5df0SBarry Smith +  coloring - the coloring context
38561ab5df0SBarry Smith -  type - either MATMFFD_WP or MATMFFD_DS
38661ab5df0SBarry Smith 
38761ab5df0SBarry Smith    Options Database Keys:
38861ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
38961ab5df0SBarry Smith 
39061ab5df0SBarry Smith    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
39161ab5df0SBarry 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
39261ab5df0SBarry Smith          introducing another one.
39361ab5df0SBarry Smith 
39461ab5df0SBarry Smith    Level: intermediate
39561ab5df0SBarry Smith 
39661ab5df0SBarry Smith .keywords: Mat, finite differences, parameters
39761ab5df0SBarry Smith 
39861ab5df0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
39961ab5df0SBarry Smith 
40061ab5df0SBarry Smith @*/
40161ab5df0SBarry Smith PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
40261ab5df0SBarry Smith {
40361ab5df0SBarry Smith   PetscFunctionBegin;
40461ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
40561ab5df0SBarry Smith   /*
40661ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
40761ab5df0SBarry Smith      and this function is being provided as patch for a release so it shouldn't change the implementaton
40861ab5df0SBarry Smith   */
40961ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
41061ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
41161ab5df0SBarry Smith   else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
41261ab5df0SBarry Smith   PetscFunctionReturn(0);
41361ab5df0SBarry Smith }
41461ab5df0SBarry Smith 
415146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
416005c665bSBarry Smith {
417dfbe8321SBarry Smith   PetscErrorCode    ierr;
418e350db43SBarry Smith   PetscBool         flg;
4193050cee2SBarry Smith   PetscViewer       viewer;
420cffb1e40SBarry Smith   PetscViewerFormat format;
421005c665bSBarry Smith 
4223a40ed3dSBarry Smith   PetscFunctionBegin;
423146574abSBarry Smith   if (prefix) {
424146574abSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
425146574abSBarry Smith   } else {
426ce94432eSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
427146574abSBarry Smith   }
428005c665bSBarry Smith   if (flg) {
429cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
4303050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
431cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
432cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
433005c665bSBarry Smith   }
4343a40ed3dSBarry Smith   PetscFunctionReturn(0);
435bbf0e169SBarry Smith }
436bbf0e169SBarry Smith 
43705869f15SSatish Balay /*@
438639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
439639f9d9dSBarry Smith    computation of Jacobians.
440bbf0e169SBarry Smith 
441ef5ee4d1SLois Curfman McInnes    Collective on Mat
442ef5ee4d1SLois Curfman McInnes 
443639f9d9dSBarry Smith    Input Parameters:
444ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
4457086a01eSPeter Brune -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
446bbf0e169SBarry Smith 
447bbf0e169SBarry Smith     Output Parameter:
448639f9d9dSBarry Smith .   color - the new coloring context
449bbf0e169SBarry Smith 
45015091d37SBarry Smith     Level: intermediate
45115091d37SBarry Smith 
4528d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
453d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
4547086a01eSPeter Brune           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
455bbf0e169SBarry Smith @*/
4567087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
457bbf0e169SBarry Smith {
458639f9d9dSBarry Smith   MatFDColoring  c;
459639f9d9dSBarry Smith   MPI_Comm       comm;
460dfbe8321SBarry Smith   PetscErrorCode ierr;
46138baddfdSBarry Smith   PetscInt       M,N;
462639f9d9dSBarry Smith 
4633a40ed3dSBarry Smith   PetscFunctionBegin;
464c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
465f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
466d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
467639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
468ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
469f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
47073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
471639f9d9dSBarry Smith 
472b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
473b8f8c88eSHong Zhang 
47493dfae19SHong Zhang   if (mat->ops->fdcoloringcreate) {
47593dfae19SHong Zhang     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
47693dfae19SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
47793dfae19SHong Zhang 
4782a7a6963SBarry Smith   ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
4793bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
480b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
4813bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
482b8f8c88eSHong Zhang 
48377d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
48477d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
48549b058dcSBarry Smith   c->currentcolor = -1;
486efb30889SBarry Smith   c->htype        = "wp";
4874e269d77SPeter Brune   c->fset         = PETSC_FALSE;
488c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
489639f9d9dSBarry Smith 
490639f9d9dSBarry Smith   *color = c;
4914e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
492d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
494bbf0e169SBarry Smith }
495bbf0e169SBarry Smith 
49605869f15SSatish Balay /*@
497639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
498639f9d9dSBarry Smith     via MatFDColoringCreate().
499bbf0e169SBarry Smith 
500ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
501ef5ee4d1SLois Curfman McInnes 
502b4fc646aSLois Curfman McInnes     Input Parameter:
503639f9d9dSBarry Smith .   c - coloring context
504bbf0e169SBarry Smith 
50515091d37SBarry Smith     Level: intermediate
50615091d37SBarry Smith 
507639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
508bbf0e169SBarry Smith @*/
5096bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
510bbf0e169SBarry Smith {
5116849ba73SBarry Smith   PetscErrorCode ierr;
51238baddfdSBarry Smith   PetscInt       i;
5130df34763SHong Zhang   MatFDColoring  color = *c;
514bbf0e169SBarry Smith 
5153a40ed3dSBarry Smith   PetscFunctionBegin;
5166bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
5170df34763SHong Zhang   if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);}
518d4bb536fSBarry Smith 
5190df34763SHong Zhang   for (i=0; i<color->ncolors; i++) {
5200df34763SHong Zhang     ierr = PetscFree(color->columns[i]);CHKERRQ(ierr);
521bbf0e169SBarry Smith   }
5220df34763SHong Zhang   ierr = PetscFree(color->ncolumns);CHKERRQ(ierr);
5230df34763SHong Zhang   ierr = PetscFree(color->columns);CHKERRQ(ierr);
5240df34763SHong Zhang   ierr = PetscFree(color->nrows);CHKERRQ(ierr);
5250df34763SHong Zhang   if (color->htype[0] == 'w') {
5260df34763SHong Zhang     ierr = PetscFree(color->matentry2);CHKERRQ(ierr);
5270df34763SHong Zhang   } else {
5280df34763SHong Zhang     ierr = PetscFree(color->matentry);CHKERRQ(ierr);
5290df34763SHong Zhang   }
5300df34763SHong Zhang   ierr = PetscFree(color->dy);CHKERRQ(ierr);
5310df34763SHong Zhang   if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);}
5320df34763SHong Zhang   ierr = VecDestroy(&color->w1);CHKERRQ(ierr);
5330df34763SHong Zhang   ierr = VecDestroy(&color->w2);CHKERRQ(ierr);
5340df34763SHong Zhang   ierr = VecDestroy(&color->w3);CHKERRQ(ierr);
535d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5363a40ed3dSBarry Smith   PetscFunctionReturn(0);
537bbf0e169SBarry Smith }
53843a90d84SBarry Smith 
53949b058dcSBarry Smith /*@C
54049b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
54149b058dcSBarry Smith       that are currently being perturbed.
54249b058dcSBarry Smith 
54349b058dcSBarry Smith     Not Collective
54449b058dcSBarry Smith 
54549b058dcSBarry Smith     Input Parameters:
54649b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
54749b058dcSBarry Smith 
54849b058dcSBarry Smith     Output Parameters:
54949b058dcSBarry Smith +   n - the number of local columns being perturbed
55049b058dcSBarry Smith -   cols - the column indices, in global numbering
55149b058dcSBarry Smith 
55221fcc2ddSBarry Smith    Level: advanced
55349b058dcSBarry Smith 
554edaa7c33SBarry Smith    Fortran Note:
555edaa7c33SBarry Smith    This routine has a different interface for Fortran
55621fcc2ddSBarry Smith $     #include <petsc/finclude/petscmat.h>
55721fcc2ddSBarry Smith $          use petscmat
558edaa7c33SBarry Smith $          PetscInt, pointer :: array(:)
559edaa7c33SBarry Smith $          PetscErrorCode  ierr
560edaa7c33SBarry Smith $          MatFDColoring   i
561edaa7c33SBarry Smith $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
562edaa7c33SBarry Smith $      use the entries of array ...
563edaa7c33SBarry Smith $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
564edaa7c33SBarry Smith 
56549b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
56649b058dcSBarry Smith 
56749b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
56849b058dcSBarry Smith @*/
569edaa7c33SBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[])
57049b058dcSBarry Smith {
57149b058dcSBarry Smith   PetscFunctionBegin;
57249b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
57349b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
57449b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
57549b058dcSBarry Smith   } else {
57649b058dcSBarry Smith     *n = 0;
57749b058dcSBarry Smith   }
57849b058dcSBarry Smith   PetscFunctionReturn(0);
57949b058dcSBarry Smith }
58049b058dcSBarry Smith 
58143a90d84SBarry Smith /*@
582e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
583e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
58443a90d84SBarry Smith 
585fee21e36SBarry Smith     Collective on MatFDColoring
586fee21e36SBarry Smith 
587ef5ee4d1SLois Curfman McInnes     Input Parameters:
588ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
589ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
590ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5917850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
592ef5ee4d1SLois Curfman McInnes 
5938bba8e72SBarry Smith     Options Database Keys:
594ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
595b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
596e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
597e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5988bba8e72SBarry Smith 
59915091d37SBarry Smith     Level: intermediate
60098d79826SSatish Balay 
6017850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
60243a90d84SBarry Smith 
60343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
60443a90d84SBarry Smith @*/
605d1e9a80fSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
60643a90d84SBarry Smith {
6073acb8795SBarry Smith   PetscErrorCode ierr;
608684f2004SHong Zhang   PetscBool      flg = PETSC_FALSE;
6093acb8795SBarry Smith 
6103acb8795SBarry Smith   PetscFunctionBegin;
6110700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
6120700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
6130700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
614ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
615e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
616c8a9c622SHong Zhang   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
617684f2004SHong Zhang 
618684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
619c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
620684f2004SHong Zhang   if (flg) {
621684f2004SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
622684f2004SHong Zhang   } else {
623684f2004SHong Zhang     PetscBool assembled;
624684f2004SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
625684f2004SHong Zhang     if (assembled) {
626684f2004SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
627684f2004SHong Zhang     }
628684f2004SHong Zhang   }
629684f2004SHong Zhang 
6305922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
631d1e9a80fSBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr);
6325922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6335bdb020cSBarry Smith   if (!coloring->viewed) {
6345bdb020cSBarry Smith     ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr);
6355bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6365bdb020cSBarry Smith   }
6373acb8795SBarry Smith   PetscFunctionReturn(0);
6383acb8795SBarry Smith }
639