xref: /petsc/src/mat/matfd/fdmatrix.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1bbf0e169SBarry Smith /*
2639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
3639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
4bbf0e169SBarry Smith */
5bbf0e169SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
7af0996ceSBarry Smith #include <petsc/private/isimpl.h>
8bbf0e169SBarry Smith 
MatFDColoringSetF(MatFDColoring fd,Vec F)9d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F)
10d71ae5a4SJacob Faibussowitsch {
113a7fca6bSBarry Smith   PetscFunctionBegin;
124e269d77SPeter Brune   if (F) {
139566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, fd->w1));
144e269d77SPeter Brune     fd->fset = PETSC_TRUE;
154e269d77SPeter Brune   } else {
164e269d77SPeter Brune     fd->fset = PETSC_FALSE;
174e269d77SPeter Brune   }
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
193a7fca6bSBarry Smith }
203a7fca6bSBarry Smith 
219804daf3SBarry Smith #include <petscdraw.h>
MatFDColoringView_Draw_Zoom(PetscDraw draw,void * Aa)22d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa)
23d71ae5a4SJacob Faibussowitsch {
249194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
256497c311SBarry Smith   PetscMPIInt   i, j, nz;
266497c311SBarry Smith   PetscInt      row;
279194eea9SBarry Smith   PetscReal     x, y;
28b312708bSHong Zhang   MatEntry     *Jentry = fd->matentry;
299194eea9SBarry Smith 
309194eea9SBarry Smith   PetscFunctionBegin;
319194eea9SBarry Smith   /* loop over colors  */
32b312708bSHong Zhang   nz = 0;
339194eea9SBarry Smith   for (i = 0; i < fd->ncolors; i++) {
349194eea9SBarry Smith     for (j = 0; j < fd->nrows[i]; j++) {
35b312708bSHong Zhang       row = Jentry[nz].row;
36b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
37b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
389566063dSJacob Faibussowitsch       PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1));
399194eea9SBarry Smith     }
409194eea9SBarry Smith   }
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429194eea9SBarry Smith }
439194eea9SBarry Smith 
MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)44d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer)
45d71ae5a4SJacob Faibussowitsch {
46ace3abfcSBarry Smith   PetscBool isnull;
47b0a32e0cSBarry Smith   PetscDraw draw;
4854d96fa3SBarry Smith   PetscReal xr, yr, xl, yl, h, w;
49005c665bSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
529566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
533ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
54005c665bSBarry Smith 
559371c9d4SSatish Balay   xr = fd->N;
569371c9d4SSatish Balay   yr = fd->M;
579371c9d4SSatish Balay   h  = yr / 10.0;
589371c9d4SSatish Balay   w  = xr / 10.0;
599371c9d4SSatish Balay   xr += w;
609371c9d4SSatish Balay   yr += h;
619371c9d4SSatish Balay   xl = -w;
629371c9d4SSatish Balay   yl = -h;
639566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer));
659566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd));
669566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69005c665bSBarry Smith }
70005c665bSBarry Smith 
71ffeef943SBarry Smith /*@
72639f9d9dSBarry Smith   MatFDColoringView - Views a finite difference coloring context.
73bbf0e169SBarry Smith 
742ef1f0ffSBarry Smith   Collective
75fee21e36SBarry Smith 
76ef5ee4d1SLois Curfman McInnes   Input Parameters:
77ef5ee4d1SLois Curfman McInnes + c      - the coloring context
78ef5ee4d1SLois Curfman McInnes - viewer - visualization context
79ef5ee4d1SLois Curfman McInnes 
8015091d37SBarry Smith   Level: intermediate
8115091d37SBarry Smith 
82b4fc646aSLois Curfman McInnes   Notes:
83b4fc646aSLois Curfman McInnes   The available visualization contexts include
8411a5261eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
8511a5261eSBarry Smith .     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
86ef5ee4d1SLois Curfman McInnes   output where only the first processor opens
87ef5ee4d1SLois Curfman McInnes   the file.  All other processors send their
88ef5ee4d1SLois Curfman McInnes   data to the first processor to print.
8911a5261eSBarry Smith -     `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure
90bbf0e169SBarry Smith 
919194eea9SBarry Smith   Since PETSc uses only a small number of basic colors (currently 33), if the coloring
92b0a32e0cSBarry Smith   involves more than 33 then some seemingly identical colors are displayed making it look
939194eea9SBarry Smith   like an illegal coloring. This is just a graphical artifact.
949194eea9SBarry Smith 
9567be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
96bbf0e169SBarry Smith @*/
MatFDColoringView(MatFDColoring c,PetscViewer viewer)97d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
98d71ae5a4SJacob Faibussowitsch {
9938baddfdSBarry Smith   PetscInt          i, j;
100*9f196a02SMartin Diehl   PetscBool         isdraw, isascii;
101f3ef73ceSBarry Smith   PetscViewerFormat format;
102bbf0e169SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
1040700a824SBarry Smith   PetscValidHeaderSpecific(c, MAT_FDCOLORING_CLASSID, 1);
10548a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
1060700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
107c9780b6fSBarry Smith   PetscCheckSameComm(c, 1, viewer, 2);
108bbf0e169SBarry Smith 
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
110*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1110f5bd95cSBarry Smith   if (isdraw) {
1129566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView_Draw(c, viewer));
113*9f196a02SMartin Diehl   } else if (isascii) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Error tolerance=%g\n", (double)c->error_rel));
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Umin=%g\n", (double)c->umin));
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of colors=%" PetscInt_FMT "\n", c->ncolors));
118ae09f205SBarry Smith 
1199566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
120fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
121b312708bSHong Zhang       PetscInt row, col, nz;
122b312708bSHong Zhang       nz = 0;
123b4fc646aSLois Curfman McInnes       for (i = 0; i < c->ncolors; i++) {
1249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Information for color %" PetscInt_FMT "\n", i));
1259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
12648a46eb9SPierre Jolivet         for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT "\n", c->columns[i][j]));
1279566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
1285bdb020cSBarry Smith         if (c->matentry) {
129b4fc646aSLois Curfman McInnes           for (j = 0; j < c->nrows[i]; j++) {
130b312708bSHong Zhang             row = c->matentry[nz].row;
131b312708bSHong Zhang             col = c->matentry[nz++].col;
1329566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
133b4fc646aSLois Curfman McInnes           }
134bbf0e169SBarry Smith         }
135bbf0e169SBarry Smith       }
1365bdb020cSBarry Smith     }
1379566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
138bbf0e169SBarry Smith   }
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140639f9d9dSBarry Smith }
141639f9d9dSBarry Smith 
142639f9d9dSBarry Smith /*@
14326a11704SBarry Smith   MatFDColoringSetParameters - Sets the parameters for the approximation of
14426a11704SBarry Smith   a sparse Jacobian matrix using finite differences and matrix coloring
145639f9d9dSBarry Smith 
146c3339decSBarry Smith   Logically Collective
147ef5ee4d1SLois Curfman McInnes 
1482fe279fdSBarry Smith   Input Parameters:
1492fe279fdSBarry Smith + matfd - the coloring context
1502fe279fdSBarry Smith . error - relative error
1512fe279fdSBarry Smith - umin  - minimum allowable u-value magnitude
1522fe279fdSBarry Smith 
1532fe279fdSBarry Smith   Level: advanced
1542fe279fdSBarry Smith 
1552fe279fdSBarry Smith   Note:
156ef5ee4d1SLois Curfman McInnes   The Jacobian is estimated with the differencing approximation
157ef5ee4d1SLois Curfman McInnes .vb
15865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
159d461c19dSHong Zhang        htype = 'ds':
160f23b5b22SLois Curfman McInnes          h = error_rel*u[i]                 if  abs(u[i]) > umin
161f23b5b22SLois Curfman McInnes            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
162ef5ee4d1SLois Curfman McInnes          dx_{i} = (0, ... 1, .... 0)
163d461c19dSHong Zhang 
164d461c19dSHong Zhang        htype = 'wp':
165d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
166ef5ee4d1SLois Curfman McInnes .ve
167639f9d9dSBarry Smith 
16867be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
169639f9d9dSBarry Smith @*/
MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)170d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
171d71ae5a4SJacob Faibussowitsch {
1723a40ed3dSBarry Smith   PetscFunctionBegin;
1730700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
174c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, error, 2);
175c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, umin, 3);
17613bcc0bdSJacob Faibussowitsch   if (error != (PetscReal)PETSC_DEFAULT) matfd->error_rel = error;
17713bcc0bdSJacob Faibussowitsch   if (umin != (PetscReal)PETSC_DEFAULT) matfd->umin = umin;
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179639f9d9dSBarry Smith }
180639f9d9dSBarry Smith 
181f86b9fbaSHong Zhang /*@
182c8a9c622SHong Zhang   MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
183005c665bSBarry Smith 
184c3339decSBarry Smith   Logically Collective
185f86b9fbaSHong Zhang 
186f86b9fbaSHong Zhang   Input Parameters:
187fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
188f86b9fbaSHong Zhang . brows - number of rows in the block
189f86b9fbaSHong Zhang - bcols - number of columns in the block
190f86b9fbaSHong Zhang 
191f86b9fbaSHong Zhang   Level: intermediate
192f86b9fbaSHong Zhang 
19367be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
194f86b9fbaSHong Zhang @*/
MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)195d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
196d71ae5a4SJacob Faibussowitsch {
197f86b9fbaSHong Zhang   PetscFunctionBegin;
198f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
199f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, brows, 2);
200f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, bcols, 3);
201f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
202f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204f86b9fbaSHong Zhang }
205f86b9fbaSHong Zhang 
206f86b9fbaSHong Zhang /*@
207f86b9fbaSHong Zhang   MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
208f86b9fbaSHong Zhang 
209c3339decSBarry Smith   Collective
210f86b9fbaSHong Zhang 
211f86b9fbaSHong Zhang   Input Parameters:
212f86b9fbaSHong Zhang + mat        - the matrix containing the nonzero structure of the Jacobian
21311a5261eSBarry Smith . iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()`
214f86b9fbaSHong Zhang - color      - the matrix coloring context
215f86b9fbaSHong Zhang 
216f86b9fbaSHong Zhang   Level: beginner
217f86b9fbaSHong Zhang 
21811a5261eSBarry Smith   Notes:
21911a5261eSBarry Smith   When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
220bdaf1daeSBarry Smith 
22167be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
222f86b9fbaSHong Zhang @*/
MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)223d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
224d71ae5a4SJacob Faibussowitsch {
225bdaf1daeSBarry Smith   PetscBool eq;
226f86b9fbaSHong Zhang 
227f86b9fbaSHong Zhang   PetscFunctionBegin;
228c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
229c8a9c622SHong Zhang   PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3);
2303ba16761SJacob Faibussowitsch   if (color->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
2319566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
23228b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
233c8a9c622SHong Zhang 
2349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
235dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
236c8a9c622SHong Zhang 
237c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240f86b9fbaSHong Zhang }
241ff0cfa39SBarry Smith 
2424a9d489dSBarry Smith /*@C
2434a9d489dSBarry Smith   MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2444a9d489dSBarry Smith 
2453f9fe445SBarry Smith   Not Collective
2464a9d489dSBarry Smith 
247f899ff85SJose E. Roman   Input Parameter:
248fe59aa6dSJacob Faibussowitsch . matfd - the coloring context
2494a9d489dSBarry Smith 
2504a9d489dSBarry Smith   Output Parameters:
2512ba42892SBarry Smith + f    - the function, see `MatFDColoringFn` for the calling sequence
2524a9d489dSBarry Smith - fctx - the optional user-defined function context
2534a9d489dSBarry Smith 
2544a9d489dSBarry Smith   Level: intermediate
2554a9d489dSBarry Smith 
2562ba42892SBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
2574a9d489dSBarry Smith @*/
MatFDColoringGetFunction(MatFDColoring matfd,MatFDColoringFn ** f,void ** fctx)2582ba42892SBarry Smith PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, MatFDColoringFn **f, void **fctx)
259d71ae5a4SJacob Faibussowitsch {
2604a9d489dSBarry Smith   PetscFunctionBegin;
2610700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
2624a9d489dSBarry Smith   if (f) *f = matfd->f;
2634a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2654a9d489dSBarry Smith }
2664a9d489dSBarry Smith 
267d64ed03dSBarry Smith /*@C
268005c665bSBarry Smith   MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
269005c665bSBarry Smith 
270c3339decSBarry Smith   Logically Collective
271fee21e36SBarry Smith 
272ef5ee4d1SLois Curfman McInnes   Input Parameters:
273fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
2742ba42892SBarry Smith . f     - the function, see `MatFDColoringFn` for the calling sequence
275ef5ee4d1SLois Curfman McInnes - fctx  - the optional user-defined function context
276ef5ee4d1SLois Curfman McInnes 
2772920cce0SJacob Faibussowitsch   Level: advanced
2782920cce0SJacob Faibussowitsch 
2792920cce0SJacob Faibussowitsch   Note:
28011a5261eSBarry Smith   This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
28111a5261eSBarry Smith   `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
28211a5261eSBarry Smith   calling `MatFDColoringApply()`
2837850c7c0SBarry Smith 
2842ba42892SBarry Smith   Fortran Note:
28511a5261eSBarry Smith   In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
28611a5261eSBarry Smith   be used without `SNES` or within the `SNES` solvers.
287f881d145SBarry Smith 
2882ba42892SBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
289005c665bSBarry Smith @*/
MatFDColoringSetFunction(MatFDColoring matfd,MatFDColoringFn * f,void * fctx)2902ba42892SBarry Smith PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, MatFDColoringFn *f, void *fctx)
291d71ae5a4SJacob Faibussowitsch {
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2930700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
294005c665bSBarry Smith   matfd->f    = f;
295005c665bSBarry Smith   matfd->fctx = fctx;
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297005c665bSBarry Smith }
298005c665bSBarry Smith 
299639f9d9dSBarry Smith /*@
300b4fc646aSLois Curfman McInnes   MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
301639f9d9dSBarry Smith   the options database.
302639f9d9dSBarry Smith 
303c3339decSBarry Smith   Collective
304fee21e36SBarry Smith 
30565f2ba5bSLois Curfman McInnes   The Jacobian, F'(u), is estimated with the differencing approximation
306ef5ee4d1SLois Curfman McInnes .vb
30765f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
308f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
309f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
310ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
311ef5ee4d1SLois Curfman McInnes .ve
312ef5ee4d1SLois Curfman McInnes 
313ef5ee4d1SLois Curfman McInnes   Input Parameter:
314fe59aa6dSJacob Faibussowitsch . matfd - the coloring context
315ef5ee4d1SLois Curfman McInnes 
316b4fc646aSLois Curfman McInnes   Options Database Keys:
3175620d6dcSBarry Smith + -mat_fd_coloring_err <err>         - Sets <err> (square root of relative error in the function)
318f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin>       - Sets umin, the minimum allowable u-value magnitude
31926a11704SBarry Smith . -mat_fd_type                       - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
320ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view              - Activates basic viewing
321e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info
322e350db43SBarry Smith - -mat_fd_coloring_view draw         - Activates drawing
323639f9d9dSBarry Smith 
32415091d37SBarry Smith   Level: intermediate
32515091d37SBarry Smith 
32667be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
327639f9d9dSBarry Smith @*/
MatFDColoringSetFromOptions(MatFDColoring matfd)328d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
329d71ae5a4SJacob Faibussowitsch {
330ace3abfcSBarry Smith   PetscBool flg;
331efb30889SBarry Smith   char      value[3];
3323a40ed3dSBarry Smith 
3333a40ed3dSBarry Smith   PetscFunctionBegin;
3340700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
335639f9d9dSBarry Smith 
336d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)matfd);
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
3399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
340efb30889SBarry Smith   if (flg) {
341efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
342efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
34398921bdaSJacob Faibussowitsch     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
344efb30889SBarry Smith   }
3459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
3469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
34793dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
34893dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
34993dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
35093dfae19SHong Zhang   }
351f86b9fbaSHong Zhang 
3525d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
353dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
354d0609cedSBarry Smith   PetscOptionsEnd();
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356005c665bSBarry Smith }
357005c665bSBarry Smith 
358cc4c1da9SBarry Smith /*@
35961ab5df0SBarry Smith   MatFDColoringSetType - Sets the approach for computing the finite difference parameter
36061ab5df0SBarry Smith 
361c3339decSBarry Smith   Collective
36261ab5df0SBarry Smith 
36361ab5df0SBarry Smith   Input Parameters:
364fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
36511a5261eSBarry Smith - type  - either `MATMFFD_WP` or `MATMFFD_DS`
36661ab5df0SBarry Smith 
3672fe279fdSBarry Smith   Options Database Key:
36861ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds"
36961ab5df0SBarry Smith 
37067be906fSBarry Smith   Level: intermediate
37167be906fSBarry Smith 
37211a5261eSBarry Smith   Note:
37311a5261eSBarry Smith   It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
37415229ffcSPierre Jolivet   but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of
37561ab5df0SBarry Smith   introducing another one.
37661ab5df0SBarry Smith 
37767be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
37861ab5df0SBarry Smith @*/
MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)379d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
380d71ae5a4SJacob Faibussowitsch {
38161ab5df0SBarry Smith   PetscFunctionBegin;
38261ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
38361ab5df0SBarry Smith   /*
38461ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
385da81f932SPierre Jolivet      and this function is being provided as patch for a release so it shouldn't change the implementation
38661ab5df0SBarry Smith   */
38761ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
38861ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
38998921bdaSJacob Faibussowitsch   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
390758219f1SBarry Smith   PetscCall(PetscObjectChangeTypeName((PetscObject)matfd, type));
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39261ab5df0SBarry Smith }
39361ab5df0SBarry Smith 
MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])39466976f2fSJacob Faibussowitsch static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
395d71ae5a4SJacob Faibussowitsch {
396e350db43SBarry Smith   PetscBool         flg;
3973050cee2SBarry Smith   PetscViewer       viewer;
398cffb1e40SBarry Smith   PetscViewerFormat format;
399005c665bSBarry Smith 
4003a40ed3dSBarry Smith   PetscFunctionBegin;
401146574abSBarry Smith   if (prefix) {
402648c30bcSBarry Smith     PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
403146574abSBarry Smith   } else {
404648c30bcSBarry Smith     PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
405146574abSBarry Smith   }
406005c665bSBarry Smith   if (flg) {
4079566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
4089566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView(fd, viewer));
4099566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
410648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
411005c665bSBarry Smith   }
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413bbf0e169SBarry Smith }
414bbf0e169SBarry Smith 
41505869f15SSatish Balay /*@
416639f9d9dSBarry Smith   MatFDColoringCreate - Creates a matrix coloring context for finite difference
417639f9d9dSBarry Smith   computation of Jacobians.
418bbf0e169SBarry Smith 
419c3339decSBarry Smith   Collective
420ef5ee4d1SLois Curfman McInnes 
421639f9d9dSBarry Smith   Input Parameters:
422ef5ee4d1SLois Curfman McInnes + mat        - the matrix containing the nonzero structure of the Jacobian
42311a5261eSBarry Smith - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
424bbf0e169SBarry Smith 
425bbf0e169SBarry Smith   Output Parameter:
426639f9d9dSBarry Smith . color - the new coloring context
427bbf0e169SBarry Smith 
42815091d37SBarry Smith   Level: intermediate
42915091d37SBarry Smith 
43067be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
431db781477SPatrick Sanan           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
432db781477SPatrick Sanan           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
433bbf0e169SBarry Smith @*/
MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring * color)434d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
435d71ae5a4SJacob Faibussowitsch {
436639f9d9dSBarry Smith   MatFDColoring c;
437639f9d9dSBarry Smith   MPI_Comm      comm;
43838baddfdSBarry Smith   PetscInt      M, N;
439639f9d9dSBarry Smith 
4403a40ed3dSBarry Smith   PetscFunctionBegin;
441c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
442377f809aSBarry Smith   PetscAssertPointer(color, 3);
44328b400f6SJacob Faibussowitsch   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
4449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
4459566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
44608401ef6SPierre Jolivet   PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4489566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
449639f9d9dSBarry Smith 
450b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
4519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
452b8f8c88eSHong Zhang 
453dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
45493dfae19SHong Zhang 
4559566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &c->w1));
456ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4579566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
4589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c->w1, &c->w2));
459ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4609566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
461b8f8c88eSHong Zhang 
46277d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
46377d8c4bbSBarry Smith   c->umin         = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
46449b058dcSBarry Smith   c->currentcolor = -1;
465efb30889SBarry Smith   c->htype        = "wp";
4664e269d77SPeter Brune   c->fset         = PETSC_FALSE;
467c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
468639f9d9dSBarry Smith 
469639f9d9dSBarry Smith   *color = c;
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
4719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473bbf0e169SBarry Smith }
474bbf0e169SBarry Smith 
47505869f15SSatish Balay /*@
476639f9d9dSBarry Smith   MatFDColoringDestroy - Destroys a matrix coloring context that was created
47711a5261eSBarry Smith   via `MatFDColoringCreate()`.
478bbf0e169SBarry Smith 
479c3339decSBarry Smith   Collective
480ef5ee4d1SLois Curfman McInnes 
481b4fc646aSLois Curfman McInnes   Input Parameter:
482639f9d9dSBarry Smith . c - coloring context
483bbf0e169SBarry Smith 
48415091d37SBarry Smith   Level: intermediate
48515091d37SBarry Smith 
48667be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
487bbf0e169SBarry Smith @*/
MatFDColoringDestroy(MatFDColoring * c)488d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
489d71ae5a4SJacob Faibussowitsch {
49038baddfdSBarry Smith   PetscInt      i;
4910df34763SHong Zhang   MatFDColoring color = *c;
492bbf0e169SBarry Smith 
4933a40ed3dSBarry Smith   PetscFunctionBegin;
4943ba16761SJacob Faibussowitsch   if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
4959371c9d4SSatish Balay   if (--((PetscObject)color)->refct > 0) {
4969371c9d4SSatish Balay     *c = NULL;
4973ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4989371c9d4SSatish Balay   }
499d4bb536fSBarry Smith 
500071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
50148a46eb9SPierre Jolivet   for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
5029566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->isa));
5039566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color->ncolumns, color->columns));
5049566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->nrows));
5050df34763SHong Zhang   if (color->htype[0] == 'w') {
5069566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry2));
5070df34763SHong Zhang   } else {
5089566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry));
5090df34763SHong Zhang   }
5109566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->dy));
5119566063dSJacob Faibussowitsch   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
5129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w1));
5139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w2));
5149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w3));
5159566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
517bbf0e169SBarry Smith }
51843a90d84SBarry Smith 
51949b058dcSBarry Smith /*@C
52049b058dcSBarry Smith   MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
52149b058dcSBarry Smith   that are currently being perturbed.
52249b058dcSBarry Smith 
52349b058dcSBarry Smith   Not Collective
52449b058dcSBarry Smith 
525f899ff85SJose E. Roman   Input Parameter:
52611a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()`
52749b058dcSBarry Smith 
52849b058dcSBarry Smith   Output Parameters:
52949b058dcSBarry Smith + n    - the number of local columns being perturbed
53049b058dcSBarry Smith - cols - the column indices, in global numbering
53149b058dcSBarry Smith 
53221fcc2ddSBarry Smith   Level: advanced
53349b058dcSBarry Smith 
53411a5261eSBarry Smith   Note:
53511a5261eSBarry Smith   IF the matrix type is `MATBAIJ`, then the block column indices are returned
536e66d0d93SMatthew Knepley 
537ce78bad3SBarry Smith   Fortran Note:
53867be906fSBarry Smith .vb
539687234bfSJose E. Roman   PetscInt, pointer :: cols(:)
54067be906fSBarry Smith .ve
541feaf08eaSBarry Smith   Use `PETSC_NULL_INTEGER` if `n` is not needed
542edaa7c33SBarry Smith 
54367be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
54449b058dcSBarry Smith @*/
MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt * n,const PetscInt * cols[])545d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
546d71ae5a4SJacob Faibussowitsch {
54749b058dcSBarry Smith   PetscFunctionBegin;
54849b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
54949b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
55049b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
55149b058dcSBarry Smith   } else {
55249b058dcSBarry Smith     *n = 0;
55349b058dcSBarry Smith   }
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55549b058dcSBarry Smith }
55649b058dcSBarry Smith 
55743a90d84SBarry Smith /*@
55811a5261eSBarry Smith   MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
559e0907662SLois Curfman McInnes   has been created, computes the Jacobian for a function via finite differences.
56043a90d84SBarry Smith 
561c3339decSBarry Smith   Collective
562fee21e36SBarry Smith 
563ef5ee4d1SLois Curfman McInnes   Input Parameters:
56426a11704SBarry Smith + J        - matrix to store Jacobian entries into
56511a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()`
566ef5ee4d1SLois Curfman McInnes . x1       - location at which Jacobian is to be computed
56726a11704SBarry Smith - sctx     - context required by function, if this is being used with the `SNES` solver then it is `SNES` object, otherwise it is `NULL`
568ef5ee4d1SLois Curfman McInnes 
5698bba8e72SBarry Smith   Options Database Keys:
57011a5261eSBarry Smith + -mat_fd_type                       - "wp" or "ds"  (see `MATMFFD_WP` or `MATMFFD_DS`)
571b9722508SLois Curfman McInnes . -mat_fd_coloring_view              - Activates basic viewing or coloring
572e350db43SBarry Smith . -mat_fd_coloring_view draw         - Activates drawing of coloring
573e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5748bba8e72SBarry Smith 
57515091d37SBarry Smith   Level: intermediate
57698d79826SSatish Balay 
57767be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
57843a90d84SBarry Smith @*/
MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void * sctx)579d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
580d71ae5a4SJacob Faibussowitsch {
581bdaf1daeSBarry Smith   PetscBool eq;
5823acb8795SBarry Smith 
5833acb8795SBarry Smith   PetscFunctionBegin;
5840700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
5850700a824SBarry Smith   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
5860700a824SBarry Smith   PetscValidHeaderSpecific(x1, VEC_CLASSID, 3);
5879566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
58828b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
58928b400f6SJacob Faibussowitsch   PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
59028b400f6SJacob Faibussowitsch   PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
591684f2004SHong Zhang 
5929566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(J));
5939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
594dbbe0bcdSBarry Smith   PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
5959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
5965bdb020cSBarry Smith   if (!coloring->viewed) {
5979566063dSJacob Faibussowitsch     PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
5985bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
5995bdb020cSBarry Smith   }
6003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6013acb8795SBarry Smith }
602