xref: /petsc/src/mat/matfd/fdmatrix.c (revision ffeef943c8ee50edff320d8a3135bb0c94853e4c)
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 
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>
22d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa)
23d71ae5a4SJacob Faibussowitsch {
249194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
25b312708bSHong Zhang   PetscInt      i, j, nz, row;
269194eea9SBarry Smith   PetscReal     x, y;
27b312708bSHong Zhang   MatEntry     *Jentry = fd->matentry;
289194eea9SBarry Smith 
299194eea9SBarry Smith   PetscFunctionBegin;
309194eea9SBarry Smith   /* loop over colors  */
31b312708bSHong Zhang   nz = 0;
329194eea9SBarry Smith   for (i = 0; i < fd->ncolors; i++) {
339194eea9SBarry Smith     for (j = 0; j < fd->nrows[i]; j++) {
34b312708bSHong Zhang       row = Jentry[nz].row;
35b312708bSHong Zhang       y   = fd->M - row - fd->rstart;
36b312708bSHong Zhang       x   = (PetscReal)Jentry[nz++].col;
379566063dSJacob Faibussowitsch       PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1));
389194eea9SBarry Smith     }
399194eea9SBarry Smith   }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
419194eea9SBarry Smith }
429194eea9SBarry Smith 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer)
44d71ae5a4SJacob Faibussowitsch {
45ace3abfcSBarry Smith   PetscBool isnull;
46b0a32e0cSBarry Smith   PetscDraw draw;
4754d96fa3SBarry Smith   PetscReal xr, yr, xl, yl, h, w;
48005c665bSBarry Smith 
493a40ed3dSBarry Smith   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
519566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
523ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
53005c665bSBarry Smith 
549371c9d4SSatish Balay   xr = fd->N;
559371c9d4SSatish Balay   yr = fd->M;
569371c9d4SSatish Balay   h  = yr / 10.0;
579371c9d4SSatish Balay   w  = xr / 10.0;
589371c9d4SSatish Balay   xr += w;
599371c9d4SSatish Balay   yr += h;
609371c9d4SSatish Balay   xl = -w;
619371c9d4SSatish Balay   yl = -h;
629566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer));
649566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd));
659566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68005c665bSBarry Smith }
69005c665bSBarry Smith 
70*ffeef943SBarry Smith /*@
71639f9d9dSBarry Smith   MatFDColoringView - Views a finite difference coloring context.
72bbf0e169SBarry Smith 
732ef1f0ffSBarry Smith   Collective
74fee21e36SBarry Smith 
75ef5ee4d1SLois Curfman McInnes   Input Parameters:
76ef5ee4d1SLois Curfman McInnes + c      - the coloring context
77ef5ee4d1SLois Curfman McInnes - viewer - visualization context
78ef5ee4d1SLois Curfman McInnes 
7915091d37SBarry Smith   Level: intermediate
8015091d37SBarry Smith 
81b4fc646aSLois Curfman McInnes   Notes:
82b4fc646aSLois Curfman McInnes   The available visualization contexts include
8311a5261eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
8411a5261eSBarry Smith .     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
85ef5ee4d1SLois Curfman McInnes   output where only the first processor opens
86ef5ee4d1SLois Curfman McInnes   the file.  All other processors send their
87ef5ee4d1SLois Curfman McInnes   data to the first processor to print.
8811a5261eSBarry Smith -     `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure
89bbf0e169SBarry Smith 
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 
9467be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
95bbf0e169SBarry Smith @*/
96d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
97d71ae5a4SJacob Faibussowitsch {
9838baddfdSBarry Smith   PetscInt          i, j;
99ace3abfcSBarry Smith   PetscBool         isdraw, iascii;
100f3ef73ceSBarry Smith   PetscViewerFormat format;
101bbf0e169SBarry Smith 
1023a40ed3dSBarry Smith   PetscFunctionBegin;
1030700a824SBarry Smith   PetscValidHeaderSpecific(c, MAT_FDCOLORING_CLASSID, 1);
10448a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
1050700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
106c9780b6fSBarry Smith   PetscCheckSameComm(c, 1, viewer, 2);
107bbf0e169SBarry Smith 
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1100f5bd95cSBarry Smith   if (isdraw) {
1119566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView_Draw(c, viewer));
11232077d6dSBarry Smith   } else if (iascii) {
1139566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
1149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Error tolerance=%g\n", (double)c->error_rel));
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Umin=%g\n", (double)c->umin));
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of colors=%" PetscInt_FMT "\n", c->ncolors));
117ae09f205SBarry Smith 
1189566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
119fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
120b312708bSHong Zhang       PetscInt row, col, nz;
121b312708bSHong Zhang       nz = 0;
122b4fc646aSLois Curfman McInnes       for (i = 0; i < c->ncolors; i++) {
1239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Information for color %" PetscInt_FMT "\n", i));
1249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
12548a46eb9SPierre Jolivet         for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT "\n", c->columns[i][j]));
1269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
1275bdb020cSBarry Smith         if (c->matentry) {
128b4fc646aSLois Curfman McInnes           for (j = 0; j < c->nrows[i]; j++) {
129b312708bSHong Zhang             row = c->matentry[nz].row;
130b312708bSHong Zhang             col = c->matentry[nz++].col;
1319566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "      %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
132b4fc646aSLois Curfman McInnes           }
133bbf0e169SBarry Smith         }
134bbf0e169SBarry Smith       }
1355bdb020cSBarry Smith     }
1369566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
137bbf0e169SBarry Smith   }
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139639f9d9dSBarry Smith }
140639f9d9dSBarry Smith 
141639f9d9dSBarry Smith /*@
142b4fc646aSLois Curfman McInnes   MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
143b4fc646aSLois Curfman McInnes   a Jacobian matrix using finite differences.
144639f9d9dSBarry Smith 
145c3339decSBarry Smith   Logically Collective
146ef5ee4d1SLois Curfman McInnes 
1472fe279fdSBarry Smith   Input Parameters:
1482fe279fdSBarry Smith + matfd - the coloring context
1492fe279fdSBarry Smith . error - relative error
1502fe279fdSBarry Smith - umin  - minimum allowable u-value magnitude
1512fe279fdSBarry Smith 
1522fe279fdSBarry Smith   Level: advanced
1532fe279fdSBarry Smith 
1542fe279fdSBarry Smith   Note:
155ef5ee4d1SLois Curfman McInnes   The Jacobian is estimated with the differencing approximation
156ef5ee4d1SLois Curfman McInnes .vb
15765f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
158d461c19dSHong Zhang        htype = 'ds':
159f23b5b22SLois Curfman McInnes          h = error_rel*u[i]                 if  abs(u[i]) > umin
160f23b5b22SLois Curfman McInnes            = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
161ef5ee4d1SLois Curfman McInnes          dx_{i} = (0, ... 1, .... 0)
162d461c19dSHong Zhang 
163d461c19dSHong Zhang        htype = 'wp':
164d461c19dSHong Zhang          h = error_rel * sqrt(1 + ||u||)
165ef5ee4d1SLois Curfman McInnes .ve
166639f9d9dSBarry Smith 
16767be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
168639f9d9dSBarry Smith @*/
169d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
170d71ae5a4SJacob Faibussowitsch {
1713a40ed3dSBarry Smith   PetscFunctionBegin;
1720700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
173c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, error, 2);
174c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd, umin, 3);
17513bcc0bdSJacob Faibussowitsch   if (error != (PetscReal)PETSC_DEFAULT) matfd->error_rel = error;
17613bcc0bdSJacob Faibussowitsch   if (umin != (PetscReal)PETSC_DEFAULT) matfd->umin = umin;
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178639f9d9dSBarry Smith }
179639f9d9dSBarry Smith 
180f86b9fbaSHong Zhang /*@
181c8a9c622SHong Zhang   MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
182005c665bSBarry Smith 
183c3339decSBarry Smith   Logically Collective
184f86b9fbaSHong Zhang 
185f86b9fbaSHong Zhang   Input Parameters:
186fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
187f86b9fbaSHong Zhang . brows - number of rows in the block
188f86b9fbaSHong Zhang - bcols - number of columns in the block
189f86b9fbaSHong Zhang 
190f86b9fbaSHong Zhang   Level: intermediate
191f86b9fbaSHong Zhang 
19267be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
193f86b9fbaSHong Zhang @*/
194d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
195d71ae5a4SJacob Faibussowitsch {
196f86b9fbaSHong Zhang   PetscFunctionBegin;
197f86b9fbaSHong Zhang   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
198f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, brows, 2);
199f86b9fbaSHong Zhang   PetscValidLogicalCollectiveInt(matfd, bcols, 3);
200f86b9fbaSHong Zhang   if (brows != PETSC_DEFAULT) matfd->brows = brows;
201f86b9fbaSHong Zhang   if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203f86b9fbaSHong Zhang }
204f86b9fbaSHong Zhang 
205f86b9fbaSHong Zhang /*@
206f86b9fbaSHong Zhang   MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
207f86b9fbaSHong Zhang 
208c3339decSBarry Smith   Collective
209f86b9fbaSHong Zhang 
210f86b9fbaSHong Zhang   Input Parameters:
211f86b9fbaSHong Zhang + mat        - the matrix containing the nonzero structure of the Jacobian
21211a5261eSBarry Smith . iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()`
213f86b9fbaSHong Zhang - color      - the matrix coloring context
214f86b9fbaSHong Zhang 
215f86b9fbaSHong Zhang   Level: beginner
216f86b9fbaSHong Zhang 
21711a5261eSBarry Smith   Notes:
21811a5261eSBarry Smith   When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
219bdaf1daeSBarry Smith 
22067be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
221f86b9fbaSHong Zhang @*/
222d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
223d71ae5a4SJacob Faibussowitsch {
224bdaf1daeSBarry Smith   PetscBool eq;
225f86b9fbaSHong Zhang 
226f86b9fbaSHong Zhang   PetscFunctionBegin;
227c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
228c8a9c622SHong Zhang   PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3);
2293ba16761SJacob Faibussowitsch   if (color->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
23128b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
232c8a9c622SHong Zhang 
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
234dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
235c8a9c622SHong Zhang 
236c8a9c622SHong Zhang   color->setupcalled = PETSC_TRUE;
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239f86b9fbaSHong Zhang }
240ff0cfa39SBarry Smith 
2414a9d489dSBarry Smith /*@C
2424a9d489dSBarry Smith   MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
2434a9d489dSBarry Smith 
2443f9fe445SBarry Smith   Not Collective
2454a9d489dSBarry Smith 
246f899ff85SJose E. Roman   Input Parameter:
247fe59aa6dSJacob Faibussowitsch . matfd - the coloring context
2484a9d489dSBarry Smith 
2494a9d489dSBarry Smith   Output Parameters:
2504a9d489dSBarry Smith + f    - the function
2514a9d489dSBarry Smith - fctx - the optional user-defined function context
2524a9d489dSBarry Smith 
2534a9d489dSBarry Smith   Level: intermediate
2544a9d489dSBarry Smith 
25567be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`
2564a9d489dSBarry Smith @*/
257d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), void **fctx)
258d71ae5a4SJacob Faibussowitsch {
2594a9d489dSBarry Smith   PetscFunctionBegin;
2600700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
2614a9d489dSBarry Smith   if (f) *f = matfd->f;
2624a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2644a9d489dSBarry Smith }
2654a9d489dSBarry Smith 
266d64ed03dSBarry Smith /*@C
267005c665bSBarry Smith   MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
268005c665bSBarry Smith 
269c3339decSBarry Smith   Logically Collective
270fee21e36SBarry Smith 
271ef5ee4d1SLois Curfman McInnes   Input Parameters:
272fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
273ef5ee4d1SLois Curfman McInnes . f     - the function
274ef5ee4d1SLois Curfman McInnes - fctx  - the optional user-defined function context
275ef5ee4d1SLois Curfman McInnes 
2762920cce0SJacob Faibussowitsch   Level: advanced
2772920cce0SJacob Faibussowitsch 
2782920cce0SJacob Faibussowitsch   Note:
2792920cce0SJacob Faibussowitsch   `f` has two possible calling configurations\:
28037fdd005SBarry Smith $ PetscErrorCode f(SNES snes, Vec in, Vec out, void *fctx)
28137fdd005SBarry Smith + snes - the nonlinear solver `SNES` object
28237fdd005SBarry Smith . in   - the location where the Jacobian is to be computed
28337fdd005SBarry Smith . out  - the location to put the computed function value
28437fdd005SBarry Smith - fctx - the function context
28520f4b53cSBarry Smith 
2862920cce0SJacob Faibussowitsch   and
28720f4b53cSBarry Smith $ PetscErrorCode f(void *dummy, Vec in, Vec out, void *fctx)
28837fdd005SBarry Smith + dummy - an unused parameter
28937fdd005SBarry Smith . in    - the location where the Jacobian is to be computed
29037fdd005SBarry Smith . out   - the location to put the computed function value
29137fdd005SBarry Smith - fctx  - the function context
29215091d37SBarry Smith 
29311a5261eSBarry Smith   This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
29411a5261eSBarry Smith   `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
29511a5261eSBarry Smith   calling `MatFDColoringApply()`
2967850c7c0SBarry Smith 
297fe59aa6dSJacob Faibussowitsch   Fortran Notes:
29811a5261eSBarry Smith   In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
29911a5261eSBarry Smith   be used without `SNES` or within the `SNES` solvers.
300f881d145SBarry Smith 
30167be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`
302005c665bSBarry Smith @*/
303d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx)
304d71ae5a4SJacob Faibussowitsch {
3053a40ed3dSBarry Smith   PetscFunctionBegin;
3060700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
307005c665bSBarry Smith   matfd->f    = f;
308005c665bSBarry Smith   matfd->fctx = fctx;
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
310005c665bSBarry Smith }
311005c665bSBarry Smith 
312639f9d9dSBarry Smith /*@
313b4fc646aSLois Curfman McInnes   MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
314639f9d9dSBarry Smith   the options database.
315639f9d9dSBarry Smith 
316c3339decSBarry Smith   Collective
317fee21e36SBarry Smith 
31865f2ba5bSLois Curfman McInnes   The Jacobian, F'(u), is estimated with the differencing approximation
319ef5ee4d1SLois Curfman McInnes .vb
32065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
321f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
322f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
323ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
324ef5ee4d1SLois Curfman McInnes .ve
325ef5ee4d1SLois Curfman McInnes 
326ef5ee4d1SLois Curfman McInnes   Input Parameter:
327fe59aa6dSJacob Faibussowitsch . matfd - the coloring context
328ef5ee4d1SLois Curfman McInnes 
329b4fc646aSLois Curfman McInnes   Options Database Keys:
3305620d6dcSBarry Smith + -mat_fd_coloring_err <err>         - Sets <err> (square root of relative error in the function)
331f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin>       - Sets umin, the minimum allowable u-value magnitude
3323ec795f1SBarry Smith . -mat_fd_type                       - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
333ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view              - Activates basic viewing
334e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info
335e350db43SBarry Smith - -mat_fd_coloring_view draw         - Activates drawing
336639f9d9dSBarry Smith 
33715091d37SBarry Smith   Level: intermediate
33815091d37SBarry Smith 
33967be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
340639f9d9dSBarry Smith @*/
341d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
342d71ae5a4SJacob Faibussowitsch {
343ace3abfcSBarry Smith   PetscBool flg;
344efb30889SBarry Smith   char      value[3];
3453a40ed3dSBarry Smith 
3463a40ed3dSBarry Smith   PetscFunctionBegin;
3470700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
348639f9d9dSBarry Smith 
349d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)matfd);
3509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
3519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
3529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
353efb30889SBarry Smith   if (flg) {
354efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
355efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
35698921bdaSJacob Faibussowitsch     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
357efb30889SBarry Smith   }
3589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
36093dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
36193dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
36293dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
36393dfae19SHong Zhang   }
364f86b9fbaSHong Zhang 
3655d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
366dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
367d0609cedSBarry Smith   PetscOptionsEnd();
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369005c665bSBarry Smith }
370005c665bSBarry Smith 
371cc4c1da9SBarry Smith /*@
37261ab5df0SBarry Smith   MatFDColoringSetType - Sets the approach for computing the finite difference parameter
37361ab5df0SBarry Smith 
374c3339decSBarry Smith   Collective
37561ab5df0SBarry Smith 
37661ab5df0SBarry Smith   Input Parameters:
377fe59aa6dSJacob Faibussowitsch + matfd - the coloring context
37811a5261eSBarry Smith - type  - either `MATMFFD_WP` or `MATMFFD_DS`
37961ab5df0SBarry Smith 
3802fe279fdSBarry Smith   Options Database Key:
38161ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds"
38261ab5df0SBarry Smith 
38367be906fSBarry Smith   Level: intermediate
38467be906fSBarry Smith 
38511a5261eSBarry Smith   Note:
38611a5261eSBarry Smith   It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
38715229ffcSPierre Jolivet   but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of
38861ab5df0SBarry Smith   introducing another one.
38961ab5df0SBarry Smith 
39067be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
39161ab5df0SBarry Smith @*/
392d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
393d71ae5a4SJacob Faibussowitsch {
39461ab5df0SBarry Smith   PetscFunctionBegin;
39561ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
39661ab5df0SBarry Smith   /*
39761ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
398da81f932SPierre Jolivet      and this function is being provided as patch for a release so it shouldn't change the implementation
39961ab5df0SBarry Smith   */
40061ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
40161ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
40298921bdaSJacob Faibussowitsch   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40461ab5df0SBarry Smith }
40561ab5df0SBarry Smith 
40666976f2fSJacob Faibussowitsch static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
407d71ae5a4SJacob Faibussowitsch {
408e350db43SBarry Smith   PetscBool         flg;
4093050cee2SBarry Smith   PetscViewer       viewer;
410cffb1e40SBarry Smith   PetscViewerFormat format;
411005c665bSBarry Smith 
4123a40ed3dSBarry Smith   PetscFunctionBegin;
413146574abSBarry Smith   if (prefix) {
4149566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
415146574abSBarry Smith   } else {
4169566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
417146574abSBarry Smith   }
418005c665bSBarry Smith   if (flg) {
4199566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
4209566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView(fd, viewer));
4219566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
422cd791dc2SBarry Smith     PetscCall(PetscOptionsRestoreViewer(&viewer));
423005c665bSBarry Smith   }
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425bbf0e169SBarry Smith }
426bbf0e169SBarry Smith 
42705869f15SSatish Balay /*@
428639f9d9dSBarry Smith   MatFDColoringCreate - Creates a matrix coloring context for finite difference
429639f9d9dSBarry Smith   computation of Jacobians.
430bbf0e169SBarry Smith 
431c3339decSBarry Smith   Collective
432ef5ee4d1SLois Curfman McInnes 
433639f9d9dSBarry Smith   Input Parameters:
434ef5ee4d1SLois Curfman McInnes + mat        - the matrix containing the nonzero structure of the Jacobian
43511a5261eSBarry Smith - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
436bbf0e169SBarry Smith 
437bbf0e169SBarry Smith   Output Parameter:
438639f9d9dSBarry Smith . color - the new coloring context
439bbf0e169SBarry Smith 
44015091d37SBarry Smith   Level: intermediate
44115091d37SBarry Smith 
44267be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
443db781477SPatrick Sanan           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
444db781477SPatrick Sanan           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
445bbf0e169SBarry Smith @*/
446d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
447d71ae5a4SJacob Faibussowitsch {
448639f9d9dSBarry Smith   MatFDColoring c;
449639f9d9dSBarry Smith   MPI_Comm      comm;
45038baddfdSBarry Smith   PetscInt      M, N;
451639f9d9dSBarry Smith 
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
45428b400f6SJacob Faibussowitsch   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
4559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
4569566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
45708401ef6SPierre Jolivet   PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
4589566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4599566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
460639f9d9dSBarry Smith 
461b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
463b8f8c88eSHong Zhang 
464dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
46593dfae19SHong Zhang 
4669566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &c->w1));
467ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4689566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
4699566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c->w1, &c->w2));
470ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4719566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
472b8f8c88eSHong Zhang 
47377d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
47477d8c4bbSBarry Smith   c->umin         = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
47549b058dcSBarry Smith   c->currentcolor = -1;
476efb30889SBarry Smith   c->htype        = "wp";
4774e269d77SPeter Brune   c->fset         = PETSC_FALSE;
478c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
479639f9d9dSBarry Smith 
480639f9d9dSBarry Smith   *color = c;
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
4829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
484bbf0e169SBarry Smith }
485bbf0e169SBarry Smith 
48605869f15SSatish Balay /*@
487639f9d9dSBarry Smith   MatFDColoringDestroy - Destroys a matrix coloring context that was created
48811a5261eSBarry Smith   via `MatFDColoringCreate()`.
489bbf0e169SBarry Smith 
490c3339decSBarry Smith   Collective
491ef5ee4d1SLois Curfman McInnes 
492b4fc646aSLois Curfman McInnes   Input Parameter:
493639f9d9dSBarry Smith . c - coloring context
494bbf0e169SBarry Smith 
49515091d37SBarry Smith   Level: intermediate
49615091d37SBarry Smith 
49767be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
498bbf0e169SBarry Smith @*/
499d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
500d71ae5a4SJacob Faibussowitsch {
50138baddfdSBarry Smith   PetscInt      i;
5020df34763SHong Zhang   MatFDColoring color = *c;
503bbf0e169SBarry Smith 
5043a40ed3dSBarry Smith   PetscFunctionBegin;
5053ba16761SJacob Faibussowitsch   if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
5069371c9d4SSatish Balay   if (--((PetscObject)color)->refct > 0) {
5079371c9d4SSatish Balay     *c = NULL;
5083ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5099371c9d4SSatish Balay   }
510d4bb536fSBarry Smith 
511071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
51248a46eb9SPierre Jolivet   for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->isa));
5149566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color->ncolumns, color->columns));
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->nrows));
5160df34763SHong Zhang   if (color->htype[0] == 'w') {
5179566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry2));
5180df34763SHong Zhang   } else {
5199566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry));
5200df34763SHong Zhang   }
5219566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->dy));
5229566063dSJacob Faibussowitsch   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
5239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w1));
5249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w2));
5259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w3));
5269566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
5273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
528bbf0e169SBarry Smith }
52943a90d84SBarry Smith 
53049b058dcSBarry Smith /*@C
53149b058dcSBarry Smith   MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
53249b058dcSBarry Smith   that are currently being perturbed.
53349b058dcSBarry Smith 
53449b058dcSBarry Smith   Not Collective
53549b058dcSBarry Smith 
536f899ff85SJose E. Roman   Input Parameter:
53711a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()`
53849b058dcSBarry Smith 
53949b058dcSBarry Smith   Output Parameters:
54049b058dcSBarry Smith + n    - the number of local columns being perturbed
54149b058dcSBarry Smith - cols - the column indices, in global numbering
54249b058dcSBarry Smith 
54321fcc2ddSBarry Smith   Level: advanced
54449b058dcSBarry Smith 
54511a5261eSBarry Smith   Note:
54611a5261eSBarry Smith   IF the matrix type is `MATBAIJ`, then the block column indices are returned
547e66d0d93SMatthew Knepley 
548fe59aa6dSJacob Faibussowitsch   Fortran Notes:
549edaa7c33SBarry Smith   This routine has a different interface for Fortran
55067be906fSBarry Smith .vb
55167be906fSBarry Smith      #include <petsc/finclude/petscmat.h>
55267be906fSBarry Smith           use petscmat
55367be906fSBarry Smith           PetscInt, pointer :: array(:)
55467be906fSBarry Smith           PetscErrorCode  ierr
55567be906fSBarry Smith           MatFDColoring   i
55667be906fSBarry Smith           call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
55767be906fSBarry Smith       use the entries of array ...
55867be906fSBarry Smith           call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
55967be906fSBarry Smith .ve
560edaa7c33SBarry Smith 
56167be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
56249b058dcSBarry Smith @*/
563d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
564d71ae5a4SJacob Faibussowitsch {
56549b058dcSBarry Smith   PetscFunctionBegin;
56649b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
56749b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
56849b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
56949b058dcSBarry Smith   } else {
57049b058dcSBarry Smith     *n = 0;
57149b058dcSBarry Smith   }
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57349b058dcSBarry Smith }
57449b058dcSBarry Smith 
57543a90d84SBarry Smith /*@
57611a5261eSBarry Smith   MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
577e0907662SLois Curfman McInnes   has been created, computes the Jacobian for a function via finite differences.
57843a90d84SBarry Smith 
579c3339decSBarry Smith   Collective
580fee21e36SBarry Smith 
581ef5ee4d1SLois Curfman McInnes   Input Parameters:
582fe59aa6dSJacob Faibussowitsch + J        - location to store Jacobian
58311a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()`
584ef5ee4d1SLois Curfman McInnes . x1       - location at which Jacobian is to be computed
58511a5261eSBarry Smith - sctx     - context required by function, if this is being used with the SNES solver then it is `SNES` object, otherwise it is null
586ef5ee4d1SLois Curfman McInnes 
5878bba8e72SBarry Smith   Options Database Keys:
58811a5261eSBarry Smith + -mat_fd_type                       - "wp" or "ds"  (see `MATMFFD_WP` or `MATMFFD_DS`)
589b9722508SLois Curfman McInnes . -mat_fd_coloring_view              - Activates basic viewing or coloring
590e350db43SBarry Smith . -mat_fd_coloring_view draw         - Activates drawing of coloring
591e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5928bba8e72SBarry Smith 
59315091d37SBarry Smith   Level: intermediate
59498d79826SSatish Balay 
59567be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
59643a90d84SBarry Smith @*/
597d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
598d71ae5a4SJacob Faibussowitsch {
599bdaf1daeSBarry Smith   PetscBool eq;
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);
6059566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
60628b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
60728b400f6SJacob Faibussowitsch   PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
60828b400f6SJacob Faibussowitsch   PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
609684f2004SHong Zhang 
6109566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(J));
6119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
612dbbe0bcdSBarry Smith   PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
6139566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
6145bdb020cSBarry Smith   if (!coloring->viewed) {
6159566063dSJacob Faibussowitsch     PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
6165bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6175bdb020cSBarry Smith   }
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6193acb8795SBarry Smith }
620