Actual source code: fdmatrix.c
1: /*
2: This is where the abstract matrix operations are defined that are
3: used for finite difference computations of Jacobians using coloring.
4: */
6: #include <petsc/private/matimpl.h>
7: #include <petsc/private/isimpl.h>
9: PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F)
10: {
11: PetscFunctionBegin;
12: if (F) {
13: PetscCall(VecCopy(F, fd->w1));
14: fd->fset = PETSC_TRUE;
15: } else {
16: fd->fset = PETSC_FALSE;
17: }
18: PetscFunctionReturn(PETSC_SUCCESS);
19: }
21: #include <petscdraw.h>
22: static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa)
23: {
24: MatFDColoring fd = (MatFDColoring)Aa;
25: PetscMPIInt i, j, nz;
26: PetscInt row;
27: PetscReal x, y;
28: MatEntry *Jentry = fd->matentry;
30: PetscFunctionBegin;
31: /* loop over colors */
32: nz = 0;
33: for (i = 0; i < fd->ncolors; i++) {
34: for (j = 0; j < fd->nrows[i]; j++) {
35: row = Jentry[nz].row;
36: y = fd->M - row - fd->rstart;
37: x = (PetscReal)Jentry[nz++].col;
38: PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1));
39: }
40: }
41: PetscFunctionReturn(PETSC_SUCCESS);
42: }
44: static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer)
45: {
46: PetscBool isnull;
47: PetscDraw draw;
48: PetscReal xr, yr, xl, yl, h, w;
50: PetscFunctionBegin;
51: PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
52: PetscCall(PetscDrawIsNull(draw, &isnull));
53: if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
55: xr = fd->N;
56: yr = fd->M;
57: h = yr / 10.0;
58: w = xr / 10.0;
59: xr += w;
60: yr += h;
61: xl = -w;
62: yl = -h;
63: PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
64: PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer));
65: PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd));
66: PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL));
67: PetscCall(PetscDrawSave(draw));
68: PetscFunctionReturn(PETSC_SUCCESS);
69: }
71: /*@
72: MatFDColoringView - Views a finite difference coloring context.
74: Collective
76: Input Parameters:
77: + c - the coloring context
78: - viewer - visualization context
80: Level: intermediate
82: Notes:
83: The available visualization contexts include
84: + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
85: . `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
86: output where only the first processor opens
87: the file. All other processors send their
88: data to the first processor to print.
89: - `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure
91: Since PETSc uses only a small number of basic colors (currently 33), if the coloring
92: involves more than 33 then some seemingly identical colors are displayed making it look
93: like an illegal coloring. This is just a graphical artifact.
95: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
96: @*/
97: PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
98: {
99: PetscBool isdraw, isascii;
100: PetscViewerFormat format;
102: PetscFunctionBegin;
104: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
106: PetscCheckSameComm(c, 1, viewer, 2);
108: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
109: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
110: if (isdraw) {
111: PetscCall(MatFDColoringView_Draw(c, viewer));
112: } else if (isascii) {
113: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
114: PetscCall(PetscViewerASCIIPrintf(viewer, " Error tolerance=%g\n", (double)c->error_rel));
115: PetscCall(PetscViewerASCIIPrintf(viewer, " Umin=%g\n", (double)c->umin));
116: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of colors=%" PetscInt_FMT "\n", c->ncolors));
118: PetscCall(PetscViewerGetFormat(viewer, &format));
119: if (format != PETSC_VIEWER_ASCII_INFO) {
120: PetscInt row, col, nz;
121: nz = 0;
122: for (PetscInt i = 0; i < c->ncolors; i++) {
123: PetscCall(PetscViewerASCIIPrintf(viewer, " Information for color %" PetscInt_FMT "\n", i));
124: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
125: for (PetscInt j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", c->columns[i][j]));
126: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
127: if (c->matentry) {
128: for (PetscInt j = 0; j < c->nrows[i]; j++) {
129: row = c->matentry[nz].row;
130: col = c->matentry[nz++].col;
131: PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
132: }
133: }
134: }
135: }
136: PetscCall(PetscViewerFlush(viewer));
137: }
138: PetscFunctionReturn(PETSC_SUCCESS);
139: }
141: /*@
142: MatFDColoringSetParameters - Sets the parameters for the approximation of
143: a sparse Jacobian matrix using finite differences and matrix coloring
145: Logically Collective
147: Input Parameters:
148: + matfd - the coloring context
149: . error - relative error
150: - umin - minimum allowable u-value magnitude
152: Level: advanced
154: Note:
155: The Jacobian is estimated with the differencing approximation
156: .vb
157: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
158: htype = 'ds':
159: h = error_rel*u[i] if abs(u[i]) > umin
160: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
161: dx_{i} = (0, ... 1, .... 0)
163: htype = 'wp':
164: h = error_rel * sqrt(1 + ||u||)
165: .ve
167: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
168: @*/
169: PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
170: {
171: PetscFunctionBegin;
175: if (error != (PetscReal)PETSC_DEFAULT) matfd->error_rel = error;
176: if (umin != (PetscReal)PETSC_DEFAULT) matfd->umin = umin;
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: /*@
181: MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
183: Logically Collective
185: Input Parameters:
186: + matfd - the coloring context
187: . brows - number of rows in the block
188: - bcols - number of columns in the block
190: Level: intermediate
192: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
193: @*/
194: PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
195: {
196: PetscFunctionBegin;
200: if (brows != PETSC_DEFAULT) matfd->brows = brows;
201: if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
202: PetscFunctionReturn(PETSC_SUCCESS);
203: }
205: /*@
206: MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
208: Collective
210: Input Parameters:
211: + mat - the matrix containing the nonzero structure of the Jacobian
212: . iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()`
213: - color - the matrix coloring context
215: Level: beginner
217: Notes:
218: When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
220: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
221: @*/
222: PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
223: {
224: PetscBool eq;
226: PetscFunctionBegin;
229: if (color->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
230: PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
231: PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
233: PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
234: PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
236: color->setupcalled = PETSC_TRUE;
237: PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*@C
242: MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
244: Not Collective
246: Input Parameter:
247: . matfd - the coloring context
249: Output Parameters:
250: + f - the function, see `MatFDColoringFn` for the calling sequence
251: - fctx - the optional user-defined function context
253: Level: intermediate
255: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
256: @*/
257: PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, MatFDColoringFn **f, void **fctx)
258: {
259: PetscFunctionBegin;
261: if (f) *f = matfd->f;
262: if (fctx) *fctx = matfd->fctx;
263: PetscFunctionReturn(PETSC_SUCCESS);
264: }
266: /*@C
267: MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
269: Logically Collective
271: Input Parameters:
272: + matfd - the coloring context
273: . f - the function, see `MatFDColoringFn` for the calling sequence
274: - fctx - the optional user-defined function context
276: Level: advanced
278: Note:
279: This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
280: `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
281: calling `MatFDColoringApply()`
283: Fortran Note:
284: In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
285: be used without `SNES` or within the `SNES` solvers.
287: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
288: @*/
289: PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, MatFDColoringFn *f, void *fctx)
290: {
291: PetscFunctionBegin;
293: matfd->f = f;
294: matfd->fctx = fctx;
295: PetscFunctionReturn(PETSC_SUCCESS);
296: }
298: /*@
299: MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
300: the options database.
302: Collective
304: The Jacobian, F'(u), is estimated with the differencing approximation
305: .vb
306: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
307: h = error_rel*u[i] if abs(u[i]) > umin
308: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
309: dx_{i} = (0, ... 1, .... 0)
310: .ve
312: Input Parameter:
313: . matfd - the coloring context
315: Options Database Keys:
316: + -mat_fd_coloring_err err - Sets err (square root of relative error in the function)
317: . -mat_fd_coloring_umin umin - Sets umin, the minimum allowable u-value magnitude
318: . -mat_fd_type (wp|ds) - See `MATMFFD_WP` and `MATMFFD_DS`
319: . -mat_fd_coloring_view - Activates basic viewing
320: . -mat_fd_coloring_view ::ascii_info - Activates viewing info
321: - -mat_fd_coloring_view draw - Activates drawing
323: Level: intermediate
325: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
326: @*/
327: PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
328: {
329: PetscBool flg;
330: char value[3];
332: PetscFunctionBegin;
335: PetscObjectOptionsBegin((PetscObject)matfd);
336: PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
337: PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
338: PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
339: if (flg) {
340: if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
341: else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
342: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
343: }
344: PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
345: PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
346: if (flg && matfd->bcols > matfd->ncolors) {
347: /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
348: matfd->bcols = matfd->ncolors;
349: }
351: /* process any options handlers added with PetscObjectAddOptionsHandler() */
352: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
353: PetscOptionsEnd();
354: PetscFunctionReturn(PETSC_SUCCESS);
355: }
357: /*@
358: MatFDColoringSetType - Sets the approach for computing the finite difference parameter
360: Collective
362: Input Parameters:
363: + matfd - the coloring context
364: - type - either `MATMFFD_WP` or `MATMFFD_DS`
366: Options Database Key:
367: . -mat_fd_type - "wp" or "ds"
369: Level: intermediate
371: Note:
372: It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
373: but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of
374: introducing another one.
376: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
377: @*/
378: PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
379: {
380: PetscFunctionBegin;
382: /*
383: It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
384: and this function is being provided as patch for a release so it shouldn't change the implementation
385: */
386: if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
387: else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
388: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
389: PetscCall(PetscObjectChangeTypeName((PetscObject)matfd, type));
390: PetscFunctionReturn(PETSC_SUCCESS);
391: }
393: static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
394: {
395: PetscBool flg;
396: PetscViewer viewer;
397: PetscViewerFormat format;
399: PetscFunctionBegin;
400: if (prefix) {
401: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
402: } else {
403: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
404: }
405: if (flg) {
406: PetscCall(PetscViewerPushFormat(viewer, format));
407: PetscCall(MatFDColoringView(fd, viewer));
408: PetscCall(PetscViewerPopFormat(viewer));
409: PetscCall(PetscViewerDestroy(&viewer));
410: }
411: PetscFunctionReturn(PETSC_SUCCESS);
412: }
414: /*@
415: MatFDColoringCreate - Creates a matrix coloring context for finite difference
416: computation of Jacobians.
418: Collective
420: Input Parameters:
421: + mat - the matrix containing the nonzero structure of the Jacobian
422: - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
424: Output Parameter:
425: . color - the new coloring context
427: Level: intermediate
429: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
430: `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
431: `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
432: @*/
433: PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
434: {
435: MatFDColoring c;
436: MPI_Comm comm;
437: PetscInt M, N;
439: PetscFunctionBegin;
441: PetscAssertPointer(color, 3);
442: PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
443: PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
444: PetscCall(MatGetSize(mat, &M, &N));
445: PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
446: PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
447: PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
449: c->ctype = iscoloring->ctype;
450: PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
452: PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
454: PetscCall(MatCreateVecs(mat, NULL, &c->w1));
455: /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
456: PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
457: PetscCall(VecDuplicate(c->w1, &c->w2));
458: /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
459: PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
461: c->error_rel = PETSC_SQRT_MACHINE_EPSILON;
462: c->umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
463: c->currentcolor = -1;
464: c->htype = "wp";
465: c->fset = PETSC_FALSE;
466: c->setupcalled = PETSC_FALSE;
468: *color = c;
469: PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
470: PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
471: PetscFunctionReturn(PETSC_SUCCESS);
472: }
474: /*@
475: MatFDColoringDestroy - Destroys a matrix coloring context that was created
476: via `MatFDColoringCreate()`.
478: Collective
480: Input Parameter:
481: . c - coloring context
483: Level: intermediate
485: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
486: @*/
487: PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
488: {
489: MatFDColoring color = *c;
491: PetscFunctionBegin;
492: if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
493: if (--((PetscObject)color)->refct > 0) {
494: *c = NULL;
495: PetscFunctionReturn(PETSC_SUCCESS);
496: }
498: /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
499: for (PetscInt i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
500: PetscCall(PetscFree(color->isa));
501: PetscCall(PetscFree2(color->ncolumns, color->columns));
502: PetscCall(PetscFree(color->nrows));
503: if (color->htype[0] == 'w') {
504: PetscCall(PetscFree(color->matentry2));
505: } else {
506: PetscCall(PetscFree(color->matentry));
507: }
508: PetscCall(PetscFree(color->dy));
509: PetscCall(VecDestroy(&color->vscale));
510: PetscCall(VecDestroy(&color->w1));
511: PetscCall(VecDestroy(&color->w2));
512: PetscCall(VecDestroy(&color->w3));
513: PetscCall(PetscHeaderDestroy(c));
514: PetscFunctionReturn(PETSC_SUCCESS);
515: }
517: /*@C
518: MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
519: that are currently being perturbed.
521: Not Collective
523: Input Parameter:
524: . coloring - coloring context created with `MatFDColoringCreate()`
526: Output Parameters:
527: + n - the number of local columns being perturbed
528: - cols - the column indices, in global numbering
530: Level: advanced
532: Note:
533: IF the matrix type is `MATBAIJ`, then the block column indices are returned
535: Fortran Note:
536: .vb
537: PetscInt, pointer :: cols(:)
538: .ve
539: Use `PETSC_NULL_INTEGER` if `n` is not needed
541: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
542: @*/
543: PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
544: {
545: PetscFunctionBegin;
546: if (coloring->currentcolor >= 0) {
547: *n = coloring->ncolumns[coloring->currentcolor];
548: *cols = coloring->columns[coloring->currentcolor];
549: } else {
550: *n = 0;
551: }
552: PetscFunctionReturn(PETSC_SUCCESS);
553: }
555: /*@
556: MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
557: has been created, computes the Jacobian for a function via finite differences.
559: Collective
561: Input Parameters:
562: + J - matrix to store Jacobian entries into
563: . coloring - coloring context created with `MatFDColoringCreate()`
564: . x1 - location at which Jacobian is to be computed
565: - sctx - context required by function, if this is being used with the `SNES` solver then it is `SNES` object, otherwise it is `NULL`
567: Options Database Keys:
568: + -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
569: . -mat_fd_coloring_view - Activates basic viewing or coloring
570: . -mat_fd_coloring_view draw - Activates drawing of coloring
571: - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
573: Level: intermediate
575: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
576: @*/
577: PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
578: {
579: PetscBool eq;
581: PetscFunctionBegin;
585: PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
586: PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
587: PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
588: PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
590: PetscCall(MatSetUnfactored(J));
591: PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
592: PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
593: PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
594: if (!coloring->viewed) {
595: PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
596: coloring->viewed = PETSC_TRUE;
597: }
598: PetscFunctionReturn(PETSC_SUCCESS);
599: }