xref: /petsc/src/mat/matfd/fdmatrix.c (revision 071fcb05f2a6aea8aef2c2090580530b9e9df77e)
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 
96bbf0e169SBarry Smith @*/
977087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
98bbf0e169SBarry Smith {
996849ba73SBarry Smith   PetscErrorCode    ierr;
10038baddfdSBarry Smith   PetscInt          i,j;
101ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
102f3ef73ceSBarry Smith   PetscViewerFormat format;
103bbf0e169SBarry Smith 
1043a40ed3dSBarry Smith   PetscFunctionBegin;
1050700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1063050cee2SBarry Smith   if (!viewer) {
107ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
1083050cee2SBarry Smith   }
1090700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
110c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
111bbf0e169SBarry Smith 
112251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
113251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1140f5bd95cSBarry Smith   if (isdraw) {
115b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11632077d6dSBarry Smith   } else if (iascii) {
117dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr);
11857622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr);
11957622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%g\n",(double)c->umin);CHKERRQ(ierr);
12077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
121ae09f205SBarry Smith 
122b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
123fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
124b312708bSHong Zhang       PetscInt row,col,nz;
125b312708bSHong Zhang       nz = 0;
126b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
129b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
13077431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
131639f9d9dSBarry Smith         }
13277431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
1335bdb020cSBarry Smith         if (c->matentry) {
134b4fc646aSLois Curfman McInnes           for (j=0; j<c->nrows[i]; j++) {
135b312708bSHong Zhang             row  = c->matentry[nz].row;
136b312708bSHong Zhang             col  = c->matentry[nz++].col;
137b312708bSHong Zhang             ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",row,col);CHKERRQ(ierr);
138b4fc646aSLois Curfman McInnes           }
139bbf0e169SBarry Smith         }
140bbf0e169SBarry Smith       }
1415bdb020cSBarry Smith     }
142b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
143bbf0e169SBarry Smith   }
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
145639f9d9dSBarry Smith }
146639f9d9dSBarry Smith 
147639f9d9dSBarry Smith /*@
148b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
149b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
150639f9d9dSBarry Smith 
1513f9fe445SBarry Smith    Logically Collective on MatFDColoring
152ef5ee4d1SLois Curfman McInnes 
153ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
154ef5ee4d1SLois Curfman McInnes .vb
15565f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
156d461c19dSHong Zhang        htype = 'ds':
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)
160d461c19dSHong Zhang 
161d461c19dSHong Zhang        htype = 'wp':
162d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
163ef5ee4d1SLois Curfman McInnes .ve
164639f9d9dSBarry Smith 
165639f9d9dSBarry Smith    Input Parameters:
166ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
167639f9d9dSBarry Smith .  error_rel - relative error
168f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
169fee21e36SBarry Smith 
17015091d37SBarry Smith    Level: advanced
17115091d37SBarry Smith 
17295b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
17395b89fc3SBarry Smith 
174639f9d9dSBarry Smith @*/
1757087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
176639f9d9dSBarry Smith {
1773a40ed3dSBarry Smith   PetscFunctionBegin;
1780700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
179c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
180c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
181639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
182639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1833a40ed3dSBarry Smith   PetscFunctionReturn(0);
184639f9d9dSBarry Smith }
185639f9d9dSBarry Smith 
186f86b9fbaSHong Zhang /*@
187c8a9c622SHong Zhang    MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
188005c665bSBarry Smith 
189f86b9fbaSHong Zhang    Logically Collective on MatFDColoring
190f86b9fbaSHong Zhang 
191f86b9fbaSHong Zhang    Input Parameters:
192f86b9fbaSHong Zhang +  coloring - the coloring context
193f86b9fbaSHong Zhang .  brows - number of rows in the block
194f86b9fbaSHong Zhang -  bcols - number of columns in the block
195f86b9fbaSHong Zhang 
196f86b9fbaSHong Zhang    Level: intermediate
197f86b9fbaSHong Zhang 
198f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
199f86b9fbaSHong Zhang 
200f86b9fbaSHong Zhang @*/
201f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)
202f86b9fbaSHong Zhang {
203f86b9fbaSHong Zhang   PetscFunctionBegin;
204f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
205f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,brows,2);
206f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd,bcols,3);
207f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
208f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
209f86b9fbaSHong Zhang   PetscFunctionReturn(0);
210f86b9fbaSHong Zhang }
211f86b9fbaSHong Zhang 
212f86b9fbaSHong Zhang /*@
213f86b9fbaSHong Zhang    MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
214f86b9fbaSHong Zhang 
215f86b9fbaSHong Zhang    Collective on Mat
216f86b9fbaSHong Zhang 
217f86b9fbaSHong Zhang    Input Parameters:
218f86b9fbaSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
219f86b9fbaSHong Zhang .  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
220f86b9fbaSHong Zhang -  color - the matrix coloring context
221f86b9fbaSHong Zhang 
222f86b9fbaSHong Zhang    Level: beginner
223f86b9fbaSHong Zhang 
224f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy()
225f86b9fbaSHong Zhang @*/
226f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)
227f86b9fbaSHong Zhang {
228f86b9fbaSHong Zhang   PetscErrorCode ierr;
229f86b9fbaSHong Zhang 
230f86b9fbaSHong Zhang   PetscFunctionBegin;
231c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
232c8a9c622SHong Zhang   PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3);
233c8a9c622SHong Zhang   if (color->setupcalled) PetscFunctionReturn(0);
234c8a9c622SHong Zhang 
2350df34763SHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
236f86b9fbaSHong Zhang   if (mat->ops->fdcoloringsetup) {
237f86b9fbaSHong Zhang     ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr);
238f86b9fbaSHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
239c8a9c622SHong Zhang 
240c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2410df34763SHong Zhang    ierr   = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr);
242f86b9fbaSHong Zhang   PetscFunctionReturn(0);
243f86b9fbaSHong Zhang }
244ff0cfa39SBarry Smith 
2454a9d489dSBarry Smith /*@C
2464a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2474a9d489dSBarry Smith 
2483f9fe445SBarry Smith    Not Collective
2494a9d489dSBarry Smith 
2504a9d489dSBarry Smith    Input Parameters:
2514a9d489dSBarry Smith .  coloring - the coloring context
2524a9d489dSBarry Smith 
2534a9d489dSBarry Smith    Output Parameters:
2544a9d489dSBarry Smith +  f - the function
2554a9d489dSBarry Smith -  fctx - the optional user-defined function context
2564a9d489dSBarry Smith 
2574a9d489dSBarry Smith    Level: intermediate
2584a9d489dSBarry Smith 
25995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
26095b89fc3SBarry Smith 
2614a9d489dSBarry Smith @*/
2627087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2634a9d489dSBarry Smith {
2644a9d489dSBarry Smith   PetscFunctionBegin;
2650700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2664a9d489dSBarry Smith   if (f) *f = matfd->f;
2674a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2684a9d489dSBarry Smith   PetscFunctionReturn(0);
2694a9d489dSBarry Smith }
2704a9d489dSBarry Smith 
271d64ed03dSBarry Smith /*@C
272005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
273005c665bSBarry Smith 
2743f9fe445SBarry Smith    Logically Collective on MatFDColoring
275fee21e36SBarry Smith 
276ef5ee4d1SLois Curfman McInnes    Input Parameters:
277ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
278ef5ee4d1SLois Curfman McInnes .  f - the function
279ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
280ef5ee4d1SLois Curfman McInnes 
2817850c7c0SBarry Smith    Calling sequence of (*f) function:
2827850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
283ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
28415091d37SBarry Smith 
2857850c7c0SBarry Smith    Level: advanced
2867850c7c0SBarry Smith 
28795452b02SPatrick Sanan    Notes:
28895452b02SPatrick Sanan     This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
2898d359177SBarry Smith      SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by
290ab637aeaSJed Brown      calling MatFDColoringApply()
2917850c7c0SBarry Smith 
2927850c7c0SBarry Smith    Fortran Notes:
2937850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
294ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
295f881d145SBarry Smith 
29695b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
29795b89fc3SBarry Smith 
298005c665bSBarry Smith @*/
2997087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
300005c665bSBarry Smith {
3013a40ed3dSBarry Smith   PetscFunctionBegin;
3020700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
303005c665bSBarry Smith   matfd->f    = f;
304005c665bSBarry Smith   matfd->fctx = fctx;
3053a40ed3dSBarry Smith   PetscFunctionReturn(0);
306005c665bSBarry Smith }
307005c665bSBarry Smith 
308639f9d9dSBarry Smith /*@
309b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
310639f9d9dSBarry Smith    the options database.
311639f9d9dSBarry Smith 
312fee21e36SBarry Smith    Collective on MatFDColoring
313fee21e36SBarry Smith 
31465f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
315ef5ee4d1SLois Curfman McInnes .vb
31665f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
317f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
318f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
319ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
320ef5ee4d1SLois Curfman McInnes .ve
321ef5ee4d1SLois Curfman McInnes 
322ef5ee4d1SLois Curfman McInnes    Input Parameter:
323ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
324ef5ee4d1SLois Curfman McInnes 
325b4fc646aSLois Curfman McInnes    Options Database Keys:
3265620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
327f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3283ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
329ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
330e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
331e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
332639f9d9dSBarry Smith 
33315091d37SBarry Smith     Level: intermediate
33415091d37SBarry Smith 
335d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
336d1c39f55SBarry Smith 
337639f9d9dSBarry Smith @*/
3387087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
339639f9d9dSBarry Smith {
340dfbe8321SBarry Smith   PetscErrorCode ierr;
341ace3abfcSBarry Smith   PetscBool      flg;
342efb30889SBarry Smith   char           value[3];
3433a40ed3dSBarry Smith 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
3450700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
346639f9d9dSBarry Smith 
3473194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
34887828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
34987828ca2SBarry Smith   ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
3501bc75a8dSBarry Smith   ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
351efb30889SBarry Smith   if (flg) {
352efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
353efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
354e32f2f54SBarry Smith     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
355efb30889SBarry Smith   }
356f86b9fbaSHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr);
35793dfae19SHong Zhang   ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr);
35893dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
35993dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
36093dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
36193dfae19SHong Zhang   }
362f86b9fbaSHong Zhang 
3635d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
3640633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr);
365b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3663a40ed3dSBarry Smith   PetscFunctionReturn(0);
367005c665bSBarry Smith }
368005c665bSBarry Smith 
36961ab5df0SBarry Smith /*@C
37061ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
37161ab5df0SBarry Smith 
37261ab5df0SBarry Smith    Collective on MatFDColoring
37361ab5df0SBarry Smith 
37461ab5df0SBarry Smith    Input Parameters:
37561ab5df0SBarry Smith +  coloring - the coloring context
37661ab5df0SBarry Smith -  type - either MATMFFD_WP or MATMFFD_DS
37761ab5df0SBarry Smith 
37861ab5df0SBarry Smith    Options Database Keys:
37961ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
38061ab5df0SBarry Smith 
38161ab5df0SBarry Smith    Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries
38261ab5df0SBarry 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
38361ab5df0SBarry Smith          introducing another one.
38461ab5df0SBarry Smith 
38561ab5df0SBarry Smith    Level: intermediate
38661ab5df0SBarry Smith 
38761ab5df0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
38861ab5df0SBarry Smith 
38961ab5df0SBarry Smith @*/
39061ab5df0SBarry Smith PetscErrorCode  MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)
39161ab5df0SBarry Smith {
39261ab5df0SBarry Smith   PetscFunctionBegin;
39361ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
39461ab5df0SBarry Smith   /*
39561ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
39661ab5df0SBarry Smith      and this function is being provided as patch for a release so it shouldn't change the implementaton
39761ab5df0SBarry Smith   */
39861ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
39961ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
40061ab5df0SBarry Smith   else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type);
40161ab5df0SBarry Smith   PetscFunctionReturn(0);
40261ab5df0SBarry Smith }
40361ab5df0SBarry Smith 
404146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])
405005c665bSBarry Smith {
406dfbe8321SBarry Smith   PetscErrorCode    ierr;
407e350db43SBarry Smith   PetscBool         flg;
4083050cee2SBarry Smith   PetscViewer       viewer;
409cffb1e40SBarry Smith   PetscViewerFormat format;
410005c665bSBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
412146574abSBarry Smith   if (prefix) {
41316413a6aSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
414146574abSBarry Smith   } else {
41516413a6aSBarry Smith     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr);
416146574abSBarry Smith   }
417005c665bSBarry Smith   if (flg) {
418cffb1e40SBarry Smith     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
4193050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
420cffb1e40SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
421cffb1e40SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
422005c665bSBarry Smith   }
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
424bbf0e169SBarry Smith }
425bbf0e169SBarry Smith 
42605869f15SSatish Balay /*@
427639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
428639f9d9dSBarry Smith    computation of Jacobians.
429bbf0e169SBarry Smith 
430ef5ee4d1SLois Curfman McInnes    Collective on Mat
431ef5ee4d1SLois Curfman McInnes 
432639f9d9dSBarry Smith    Input Parameters:
433ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
4347086a01eSPeter Brune -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
435bbf0e169SBarry Smith 
436bbf0e169SBarry Smith     Output Parameter:
437639f9d9dSBarry Smith .   color - the new coloring context
438bbf0e169SBarry Smith 
43915091d37SBarry Smith     Level: intermediate
44015091d37SBarry Smith 
4418d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(),
442d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
4437086a01eSPeter Brune           MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring()
444bbf0e169SBarry Smith @*/
4457087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
446bbf0e169SBarry Smith {
447639f9d9dSBarry Smith   MatFDColoring  c;
448639f9d9dSBarry Smith   MPI_Comm       comm;
449dfbe8321SBarry Smith   PetscErrorCode ierr;
45038baddfdSBarry Smith   PetscInt       M,N;
451639f9d9dSBarry Smith 
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
454f141eedfSHong Zhang   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
455d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
456639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
457ce94432eSBarry Smith   if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices");
458f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
45973107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
460639f9d9dSBarry Smith 
461b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
462b8f8c88eSHong Zhang 
46393dfae19SHong Zhang   if (mat->ops->fdcoloringcreate) {
46493dfae19SHong Zhang     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
46593dfae19SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);
46693dfae19SHong Zhang 
4672a7a6963SBarry Smith   ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr);
468e7e92044SBarry Smith   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
469e7e92044SBarry Smith   ierr = VecPinToCPU(c->w1,PETSC_TRUE);CHKERRQ(ierr);
4703bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr);
471b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
472e7e92044SBarry Smith   /* Vec is used instensively in particular piece of scalar CPU code; won't benifit from bouncing back and forth to the GPU */
473e7e92044SBarry Smith   ierr = VecPinToCPU(c->w2,PETSC_TRUE);CHKERRQ(ierr);
4743bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr);
475b8f8c88eSHong Zhang 
47677d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
47777d8c4bbSBarry Smith   c->umin         = 100.0*PETSC_SQRT_MACHINE_EPSILON;
47849b058dcSBarry Smith   c->currentcolor = -1;
479efb30889SBarry Smith   c->htype        = "wp";
4804e269d77SPeter Brune   c->fset         = PETSC_FALSE;
481c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
482639f9d9dSBarry Smith 
483639f9d9dSBarry Smith   *color = c;
4844e269d77SPeter Brune   ierr   = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr);
485d5ba7fb7SMatthew Knepley   ierr   = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4863a40ed3dSBarry Smith   PetscFunctionReturn(0);
487bbf0e169SBarry Smith }
488bbf0e169SBarry Smith 
48905869f15SSatish Balay /*@
490639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
491639f9d9dSBarry Smith     via MatFDColoringCreate().
492bbf0e169SBarry Smith 
493ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
494ef5ee4d1SLois Curfman McInnes 
495b4fc646aSLois Curfman McInnes     Input Parameter:
496639f9d9dSBarry Smith .   c - coloring context
497bbf0e169SBarry Smith 
49815091d37SBarry Smith     Level: intermediate
49915091d37SBarry Smith 
500639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
501bbf0e169SBarry Smith @*/
5026bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
503bbf0e169SBarry Smith {
5046849ba73SBarry Smith   PetscErrorCode ierr;
50538baddfdSBarry Smith   PetscInt       i;
5060df34763SHong Zhang   MatFDColoring  color = *c;
507bbf0e169SBarry Smith 
5083a40ed3dSBarry Smith   PetscFunctionBegin;
5096bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
5100df34763SHong Zhang   if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);}
511d4bb536fSBarry Smith 
512*071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
5130df34763SHong Zhang   for (i=0; i<color->ncolors; i++) {
514*071fcb05SBarry Smith     ierr = ISDestroy(&color->isa[i]);CHKERRQ(ierr);
515bbf0e169SBarry Smith   }
516*071fcb05SBarry Smith   ierr = PetscFree(color->isa);CHKERRQ(ierr);
517*071fcb05SBarry Smith   ierr = PetscFree2(color->ncolumns,color->columns);CHKERRQ(ierr);
5180df34763SHong Zhang   ierr = PetscFree(color->nrows);CHKERRQ(ierr);
5190df34763SHong Zhang   if (color->htype[0] == 'w') {
5200df34763SHong Zhang     ierr = PetscFree(color->matentry2);CHKERRQ(ierr);
5210df34763SHong Zhang   } else {
5220df34763SHong Zhang     ierr = PetscFree(color->matentry);CHKERRQ(ierr);
5230df34763SHong Zhang   }
5240df34763SHong Zhang   ierr = PetscFree(color->dy);CHKERRQ(ierr);
5250df34763SHong Zhang   if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);}
5260df34763SHong Zhang   ierr = VecDestroy(&color->w1);CHKERRQ(ierr);
5270df34763SHong Zhang   ierr = VecDestroy(&color->w2);CHKERRQ(ierr);
5280df34763SHong Zhang   ierr = VecDestroy(&color->w3);CHKERRQ(ierr);
529d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
5303a40ed3dSBarry Smith   PetscFunctionReturn(0);
531bbf0e169SBarry Smith }
53243a90d84SBarry Smith 
53349b058dcSBarry Smith /*@C
53449b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
53549b058dcSBarry Smith       that are currently being perturbed.
53649b058dcSBarry Smith 
53749b058dcSBarry Smith     Not Collective
53849b058dcSBarry Smith 
53949b058dcSBarry Smith     Input Parameters:
54049b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
54149b058dcSBarry Smith 
54249b058dcSBarry Smith     Output Parameters:
54349b058dcSBarry Smith +   n - the number of local columns being perturbed
54449b058dcSBarry Smith -   cols - the column indices, in global numbering
54549b058dcSBarry Smith 
54621fcc2ddSBarry Smith    Level: advanced
54749b058dcSBarry Smith 
548edaa7c33SBarry Smith    Fortran Note:
549edaa7c33SBarry Smith    This routine has a different interface for Fortran
55021fcc2ddSBarry Smith $     #include <petsc/finclude/petscmat.h>
55121fcc2ddSBarry Smith $          use petscmat
552edaa7c33SBarry Smith $          PetscInt, pointer :: array(:)
553edaa7c33SBarry Smith $          PetscErrorCode  ierr
554edaa7c33SBarry Smith $          MatFDColoring   i
555edaa7c33SBarry Smith $          call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
556edaa7c33SBarry Smith $      use the entries of array ...
557edaa7c33SBarry Smith $          call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
558edaa7c33SBarry Smith 
55949b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
56049b058dcSBarry Smith 
56149b058dcSBarry Smith @*/
562edaa7c33SBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[])
56349b058dcSBarry Smith {
56449b058dcSBarry Smith   PetscFunctionBegin;
56549b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
56649b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
56749b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
56849b058dcSBarry Smith   } else {
56949b058dcSBarry Smith     *n = 0;
57049b058dcSBarry Smith   }
57149b058dcSBarry Smith   PetscFunctionReturn(0);
57249b058dcSBarry Smith }
57349b058dcSBarry Smith 
57443a90d84SBarry Smith /*@
575e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
576e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
57743a90d84SBarry Smith 
578fee21e36SBarry Smith     Collective on MatFDColoring
579fee21e36SBarry Smith 
580ef5ee4d1SLois Curfman McInnes     Input Parameters:
581ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
582ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
583ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
5847850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
585ef5ee4d1SLois Curfman McInnes 
5868bba8e72SBarry Smith     Options Database Keys:
587ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
588b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
589e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
590e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5918bba8e72SBarry Smith 
59215091d37SBarry Smith     Level: intermediate
59398d79826SSatish Balay 
5947850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
59543a90d84SBarry Smith 
59643a90d84SBarry Smith @*/
597d1e9a80fSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx)
59843a90d84SBarry Smith {
5993acb8795SBarry Smith   PetscErrorCode ierr;
6003acb8795SBarry Smith 
6013acb8795SBarry Smith   PetscFunctionBegin;
6020700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
6030700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
6040700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
605ce94432eSBarry Smith   if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
606e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
607c8a9c622SHong Zhang   if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()");
608684f2004SHong Zhang 
609684f2004SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
6105922145eSHong Zhang   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
611d1e9a80fSBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr);
6125922145eSHong Zhang   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
6135bdb020cSBarry Smith   if (!coloring->viewed) {
6145bdb020cSBarry Smith     ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr);
6155bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6165bdb020cSBarry Smith   }
6173acb8795SBarry Smith   PetscFunctionReturn(0);
6183acb8795SBarry Smith }
619