xref: /petsc/src/mat/matfd/fdmatrix.c (revision 37fdd0055130c7e84fc09516fb012d75c18be859)
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 
10d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F)
11d71ae5a4SJacob Faibussowitsch {
123a7fca6bSBarry Smith   PetscFunctionBegin;
134e269d77SPeter Brune   if (F) {
149566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, fd->w1));
154e269d77SPeter Brune     fd->fset = PETSC_TRUE;
164e269d77SPeter Brune   } else {
174e269d77SPeter Brune     fd->fset = PETSC_FALSE;
184e269d77SPeter Brune   }
193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203a7fca6bSBarry Smith }
213a7fca6bSBarry Smith 
229804daf3SBarry Smith #include <petscdraw.h>
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa)
24d71ae5a4SJacob Faibussowitsch {
259194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
26b312708bSHong Zhang   PetscInt      i, j, nz, 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 
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 
71bbf0e169SBarry Smith /*@C
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 @*/
97d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
98d71ae5a4SJacob Faibussowitsch {
9938baddfdSBarry Smith   PetscInt          i, j;
100ace3abfcSBarry Smith   PetscBool         isdraw, iascii;
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));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1110f5bd95cSBarry Smith   if (isdraw) {
1129566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView_Draw(c, viewer));
11332077d6dSBarry Smith   } else if (iascii) {
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 /*@
143b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
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 @*/
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:
187f86b9fbaSHong Zhang +  coloring - 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 @*/
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 @*/
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:
2484a9d489dSBarry Smith .  coloring - the coloring context
2494a9d489dSBarry Smith 
2504a9d489dSBarry Smith    Output Parameters:
2514a9d489dSBarry Smith +  f - the function
2524a9d489dSBarry Smith -  fctx - the optional user-defined function context
2534a9d489dSBarry Smith 
2544a9d489dSBarry Smith    Level: intermediate
2554a9d489dSBarry Smith 
25667be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`
2574a9d489dSBarry Smith @*/
258d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), 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:
273ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
274ef5ee4d1SLois Curfman McInnes .  f - the function
275ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
276ef5ee4d1SLois Curfman McInnes 
27720f4b53cSBarry Smith    Calling sequence with `SNES` of `f`:
278*37fdd005SBarry Smith $   PetscErrorCode f(SNES snes, Vec in, Vec out, void *fctx)
279*37fdd005SBarry Smith +  snes - the nonlinear solver `SNES` object
280*37fdd005SBarry Smith .  in - the location where the Jacobian is to be computed
281*37fdd005SBarry Smith .  out - the location to put the computed function value
282*37fdd005SBarry Smith -  fctx - the function context
28320f4b53cSBarry Smith 
28420f4b53cSBarry Smith    Calling sequence without `SNES` of `f`:
28520f4b53cSBarry Smith $   PetscErrorCode f(void *dummy, Vec in, Vec out, void *fctx)
286*37fdd005SBarry Smith +  dummy - an unused parameter
287*37fdd005SBarry Smith .  in - the location where the Jacobian is to be computed
288*37fdd005SBarry Smith .  out - the location to put the computed function value
289*37fdd005SBarry Smith -  fctx - the function context
29015091d37SBarry Smith 
2917850c7c0SBarry Smith    Level: advanced
2927850c7c0SBarry Smith 
29311a5261eSBarry Smith    Note:
29411a5261eSBarry Smith     This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
29511a5261eSBarry Smith      `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
29611a5261eSBarry Smith      calling `MatFDColoringApply()`
2977850c7c0SBarry Smith 
29811a5261eSBarry Smith    Fortran Note:
29911a5261eSBarry Smith     In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
30011a5261eSBarry Smith   be used without `SNES` or within the `SNES` solvers.
301f881d145SBarry Smith 
30267be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`
303005c665bSBarry Smith @*/
304d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx)
305d71ae5a4SJacob Faibussowitsch {
3063a40ed3dSBarry Smith   PetscFunctionBegin;
3070700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
308005c665bSBarry Smith   matfd->f    = f;
309005c665bSBarry Smith   matfd->fctx = fctx;
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311005c665bSBarry Smith }
312005c665bSBarry Smith 
313639f9d9dSBarry Smith /*@
314b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
315639f9d9dSBarry Smith    the options database.
316639f9d9dSBarry Smith 
317c3339decSBarry Smith    Collective
318fee21e36SBarry Smith 
31965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
320ef5ee4d1SLois Curfman McInnes .vb
32165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
322f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
323f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
324ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
325ef5ee4d1SLois Curfman McInnes .ve
326ef5ee4d1SLois Curfman McInnes 
327ef5ee4d1SLois Curfman McInnes    Input Parameter:
328ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
329ef5ee4d1SLois Curfman McInnes 
330b4fc646aSLois Curfman McInnes    Options Database Keys:
3315620d6dcSBarry Smith +  -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
332f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
3333ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
334ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
335e350db43SBarry Smith .  -mat_fd_coloring_view ::ascii_info - Activates viewing info
336e350db43SBarry Smith -  -mat_fd_coloring_view draw - Activates drawing
337639f9d9dSBarry Smith 
33815091d37SBarry Smith     Level: intermediate
33915091d37SBarry Smith 
34067be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
341639f9d9dSBarry Smith @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
343d71ae5a4SJacob Faibussowitsch {
344ace3abfcSBarry Smith   PetscBool flg;
345efb30889SBarry Smith   char      value[3];
3463a40ed3dSBarry Smith 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
3480700a824SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
349639f9d9dSBarry Smith 
350d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)matfd);
3519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
3529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
3539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
354efb30889SBarry Smith   if (flg) {
355efb30889SBarry Smith     if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
356efb30889SBarry Smith     else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
35798921bdaSJacob Faibussowitsch     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
358efb30889SBarry Smith   }
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
36193dfae19SHong Zhang   if (flg && matfd->bcols > matfd->ncolors) {
36293dfae19SHong Zhang     /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
36393dfae19SHong Zhang     matfd->bcols = matfd->ncolors;
36493dfae19SHong Zhang   }
365f86b9fbaSHong Zhang 
3665d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
367dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
368d0609cedSBarry Smith   PetscOptionsEnd();
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
370005c665bSBarry Smith }
371005c665bSBarry Smith 
37261ab5df0SBarry Smith /*@C
37361ab5df0SBarry Smith    MatFDColoringSetType - Sets the approach for computing the finite difference parameter
37461ab5df0SBarry Smith 
375c3339decSBarry Smith    Collective
37661ab5df0SBarry Smith 
37761ab5df0SBarry Smith    Input Parameters:
37861ab5df0SBarry Smith +  coloring - the coloring context
37911a5261eSBarry Smith -  type - either `MATMFFD_WP` or `MATMFFD_DS`
38061ab5df0SBarry Smith 
3812fe279fdSBarry Smith    Options Database Key:
38261ab5df0SBarry Smith .  -mat_fd_type - "wp" or "ds"
38361ab5df0SBarry Smith 
38467be906fSBarry Smith    Level: intermediate
38567be906fSBarry Smith 
38611a5261eSBarry Smith    Note:
38711a5261eSBarry Smith    It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
38811a5261eSBarry 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
38961ab5df0SBarry Smith          introducing another one.
39061ab5df0SBarry Smith 
39167be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
39261ab5df0SBarry Smith @*/
393d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
394d71ae5a4SJacob Faibussowitsch {
39561ab5df0SBarry Smith   PetscFunctionBegin;
39661ab5df0SBarry Smith   PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
39761ab5df0SBarry Smith   /*
39861ab5df0SBarry Smith      It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
399da81f932SPierre Jolivet      and this function is being provided as patch for a release so it shouldn't change the implementation
40061ab5df0SBarry Smith   */
40161ab5df0SBarry Smith   if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
40261ab5df0SBarry Smith   else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
40398921bdaSJacob Faibussowitsch   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
4043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40561ab5df0SBarry Smith }
40661ab5df0SBarry Smith 
407d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
408d71ae5a4SJacob Faibussowitsch {
409e350db43SBarry Smith   PetscBool         flg;
4103050cee2SBarry Smith   PetscViewer       viewer;
411cffb1e40SBarry Smith   PetscViewerFormat format;
412005c665bSBarry Smith 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
414146574abSBarry Smith   if (prefix) {
4159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
416146574abSBarry Smith   } else {
4179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
418146574abSBarry Smith   }
419005c665bSBarry Smith   if (flg) {
4209566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
4219566063dSJacob Faibussowitsch     PetscCall(MatFDColoringView(fd, viewer));
4229566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
4239566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
424005c665bSBarry Smith   }
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
426bbf0e169SBarry Smith }
427bbf0e169SBarry Smith 
42805869f15SSatish Balay /*@
429639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
430639f9d9dSBarry Smith    computation of Jacobians.
431bbf0e169SBarry Smith 
432c3339decSBarry Smith    Collective
433ef5ee4d1SLois Curfman McInnes 
434639f9d9dSBarry Smith    Input Parameters:
435ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
43611a5261eSBarry Smith -  iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
437bbf0e169SBarry Smith 
438bbf0e169SBarry Smith     Output Parameter:
439639f9d9dSBarry Smith .   color - the new coloring context
440bbf0e169SBarry Smith 
44115091d37SBarry Smith     Level: intermediate
44215091d37SBarry Smith 
44367be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
444db781477SPatrick Sanan           `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
445db781477SPatrick Sanan           `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
446bbf0e169SBarry Smith @*/
447d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
448d71ae5a4SJacob Faibussowitsch {
449639f9d9dSBarry Smith   MatFDColoring c;
450639f9d9dSBarry Smith   MPI_Comm      comm;
45138baddfdSBarry Smith   PetscInt      M, N;
452639f9d9dSBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454c8a9c622SHong Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
45528b400f6SJacob Faibussowitsch   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
4569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
4579566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
45808401ef6SPierre Jolivet   PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4609566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
461639f9d9dSBarry Smith 
462b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
464b8f8c88eSHong Zhang 
465dbbe0bcdSBarry Smith   PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
46693dfae19SHong Zhang 
4679566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &c->w1));
468ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4699566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
4709566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c->w1, &c->w2));
471ce911d59SRichard Tran Mills   /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
4729566063dSJacob Faibussowitsch   PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
473b8f8c88eSHong Zhang 
47477d8c4bbSBarry Smith   c->error_rel    = PETSC_SQRT_MACHINE_EPSILON;
47577d8c4bbSBarry Smith   c->umin         = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
47649b058dcSBarry Smith   c->currentcolor = -1;
477efb30889SBarry Smith   c->htype        = "wp";
4784e269d77SPeter Brune   c->fset         = PETSC_FALSE;
479c8a9c622SHong Zhang   c->setupcalled  = PETSC_FALSE;
480639f9d9dSBarry Smith 
481639f9d9dSBarry Smith   *color = c;
4829566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
4839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485bbf0e169SBarry Smith }
486bbf0e169SBarry Smith 
48705869f15SSatish Balay /*@
488639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
48911a5261eSBarry Smith     via `MatFDColoringCreate()`.
490bbf0e169SBarry Smith 
491c3339decSBarry Smith     Collective
492ef5ee4d1SLois Curfman McInnes 
493b4fc646aSLois Curfman McInnes     Input Parameter:
494639f9d9dSBarry Smith .   c - coloring context
495bbf0e169SBarry Smith 
49615091d37SBarry Smith     Level: intermediate
49715091d37SBarry Smith 
49867be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
499bbf0e169SBarry Smith @*/
500d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
501d71ae5a4SJacob Faibussowitsch {
50238baddfdSBarry Smith   PetscInt      i;
5030df34763SHong Zhang   MatFDColoring color = *c;
504bbf0e169SBarry Smith 
5053a40ed3dSBarry Smith   PetscFunctionBegin;
5063ba16761SJacob Faibussowitsch   if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
5079371c9d4SSatish Balay   if (--((PetscObject)color)->refct > 0) {
5089371c9d4SSatish Balay     *c = NULL;
5093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5109371c9d4SSatish Balay   }
511d4bb536fSBarry Smith 
512071fcb05SBarry Smith   /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
51348a46eb9SPierre Jolivet   for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
5149566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->isa));
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color->ncolumns, color->columns));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->nrows));
5170df34763SHong Zhang   if (color->htype[0] == 'w') {
5189566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry2));
5190df34763SHong Zhang   } else {
5209566063dSJacob Faibussowitsch     PetscCall(PetscFree(color->matentry));
5210df34763SHong Zhang   }
5229566063dSJacob Faibussowitsch   PetscCall(PetscFree(color->dy));
5239566063dSJacob Faibussowitsch   if (color->vscale) PetscCall(VecDestroy(&color->vscale));
5249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w1));
5259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w2));
5269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&color->w3));
5279566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
529bbf0e169SBarry Smith }
53043a90d84SBarry Smith 
53149b058dcSBarry Smith /*@C
53249b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
53349b058dcSBarry Smith       that are currently being perturbed.
53449b058dcSBarry Smith 
53549b058dcSBarry Smith     Not Collective
53649b058dcSBarry Smith 
537f899ff85SJose E. Roman     Input Parameter:
53811a5261eSBarry Smith .   coloring - coloring context created with `MatFDColoringCreate()`
53949b058dcSBarry Smith 
54049b058dcSBarry Smith     Output Parameters:
54149b058dcSBarry Smith +   n - the number of local columns being perturbed
54249b058dcSBarry Smith -   cols - the column indices, in global numbering
54349b058dcSBarry Smith 
54421fcc2ddSBarry Smith    Level: advanced
54549b058dcSBarry Smith 
54611a5261eSBarry Smith    Note:
54711a5261eSBarry Smith    IF the matrix type is `MATBAIJ`, then the block column indices are returned
548e66d0d93SMatthew Knepley 
549edaa7c33SBarry Smith    Fortran Note:
550edaa7c33SBarry Smith    This routine has a different interface for Fortran
55167be906fSBarry Smith .vb
55267be906fSBarry Smith      #include <petsc/finclude/petscmat.h>
55367be906fSBarry Smith           use petscmat
55467be906fSBarry Smith           PetscInt, pointer :: array(:)
55567be906fSBarry Smith           PetscErrorCode  ierr
55667be906fSBarry Smith           MatFDColoring   i
55767be906fSBarry Smith           call MatFDColoringGetPerturbedColumnsF90(i,array,ierr)
55867be906fSBarry Smith       use the entries of array ...
55967be906fSBarry Smith           call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr)
56067be906fSBarry Smith .ve
561edaa7c33SBarry Smith 
56267be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
56349b058dcSBarry Smith @*/
564d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
565d71ae5a4SJacob Faibussowitsch {
56649b058dcSBarry Smith   PetscFunctionBegin;
56749b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
56849b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
56949b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
57049b058dcSBarry Smith   } else {
57149b058dcSBarry Smith     *n = 0;
57249b058dcSBarry Smith   }
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57449b058dcSBarry Smith }
57549b058dcSBarry Smith 
57643a90d84SBarry Smith /*@
57711a5261eSBarry Smith     MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
578e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
57943a90d84SBarry Smith 
580c3339decSBarry Smith     Collective
581fee21e36SBarry Smith 
582ef5ee4d1SLois Curfman McInnes     Input Parameters:
583ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
58411a5261eSBarry Smith .   coloring - coloring context created with `MatFDColoringCreate()`
585ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
58611a5261eSBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is `SNES` object, otherwise it is null
587ef5ee4d1SLois Curfman McInnes 
5888bba8e72SBarry Smith     Options Database Keys:
58911a5261eSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see `MATMFFD_WP` or `MATMFFD_DS`)
590b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
591e350db43SBarry Smith .    -mat_fd_coloring_view draw - Activates drawing of coloring
592e350db43SBarry Smith -    -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
5938bba8e72SBarry Smith 
59415091d37SBarry Smith     Level: intermediate
59598d79826SSatish Balay 
59667be906fSBarry Smith .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
59743a90d84SBarry Smith @*/
598d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
599d71ae5a4SJacob Faibussowitsch {
600bdaf1daeSBarry Smith   PetscBool eq;
6013acb8795SBarry Smith 
6023acb8795SBarry Smith   PetscFunctionBegin;
6030700a824SBarry Smith   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
6040700a824SBarry Smith   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
6050700a824SBarry Smith   PetscValidHeaderSpecific(x1, VEC_CLASSID, 3);
6069566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
60728b400f6SJacob Faibussowitsch   PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
60828b400f6SJacob Faibussowitsch   PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
60928b400f6SJacob Faibussowitsch   PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
610684f2004SHong Zhang 
6119566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(J));
6129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
613dbbe0bcdSBarry Smith   PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
6149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
6155bdb020cSBarry Smith   if (!coloring->viewed) {
6169566063dSJacob Faibussowitsch     PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
6175bdb020cSBarry Smith     coloring->viewed = PETSC_TRUE;
6185bdb020cSBarry Smith   }
6193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6203acb8795SBarry Smith }
621