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: PetscInt i, j;
100: PetscBool isdraw, iascii;
101: PetscViewerFormat format;
103: PetscFunctionBegin;
105: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
107: PetscCheckSameComm(c, 1, viewer, 2);
109: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
110: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
111: if (isdraw) {
112: PetscCall(MatFDColoringView_Draw(c, viewer));
113: } else if (iascii) {
114: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer));
115: PetscCall(PetscViewerASCIIPrintf(viewer, " Error tolerance=%g\n", (double)c->error_rel));
116: PetscCall(PetscViewerASCIIPrintf(viewer, " Umin=%g\n", (double)c->umin));
117: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of colors=%" PetscInt_FMT "\n", c->ncolors));
119: PetscCall(PetscViewerGetFormat(viewer, &format));
120: if (format != PETSC_VIEWER_ASCII_INFO) {
121: PetscInt row, col, nz;
122: nz = 0;
123: for (i = 0; i < c->ncolors; i++) {
124: PetscCall(PetscViewerASCIIPrintf(viewer, " Information for color %" PetscInt_FMT "\n", i));
125: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of columns %" PetscInt_FMT "\n", c->ncolumns[i]));
126: for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", c->columns[i][j]));
127: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of rows %" PetscInt_FMT "\n", c->nrows[i]));
128: if (c->matentry) {
129: for (j = 0; j < c->nrows[i]; j++) {
130: row = c->matentry[nz].row;
131: col = c->matentry[nz++].col;
132: PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col));
133: }
134: }
135: }
136: }
137: PetscCall(PetscViewerFlush(viewer));
138: }
139: PetscFunctionReturn(PETSC_SUCCESS);
140: }
142: /*@
143: MatFDColoringSetParameters - Sets the parameters for the approximation of
144: a sparse Jacobian matrix using finite differences and matrix coloring
146: Logically Collective
148: Input Parameters:
149: + matfd - the coloring context
150: . error - relative error
151: - umin - minimum allowable u-value magnitude
153: Level: advanced
155: Note:
156: The Jacobian is estimated with the differencing approximation
157: .vb
158: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
159: htype = 'ds':
160: h = error_rel*u[i] if abs(u[i]) > umin
161: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
162: dx_{i} = (0, ... 1, .... 0)
164: htype = 'wp':
165: h = error_rel * sqrt(1 + ||u||)
166: .ve
168: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
169: @*/
170: PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
171: {
172: PetscFunctionBegin;
176: if (error != (PetscReal)PETSC_DEFAULT) matfd->error_rel = error;
177: if (umin != (PetscReal)PETSC_DEFAULT) matfd->umin = umin;
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@
182: MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
184: Logically Collective
186: Input Parameters:
187: + matfd - the coloring context
188: . brows - number of rows in the block
189: - bcols - number of columns in the block
191: Level: intermediate
193: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
194: @*/
195: PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
196: {
197: PetscFunctionBegin;
201: if (brows != PETSC_DEFAULT) matfd->brows = brows;
202: if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: /*@
207: MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
209: Collective
211: Input Parameters:
212: + mat - the matrix containing the nonzero structure of the Jacobian
213: . iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()`
214: - color - the matrix coloring context
216: Level: beginner
218: Notes:
219: When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
221: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
222: @*/
223: PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
224: {
225: PetscBool eq;
227: PetscFunctionBegin;
230: if (color->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
231: PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq));
232: PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()");
234: PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
235: PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
237: color->setupcalled = PETSC_TRUE;
238: PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: /*@C
243: MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
245: Not Collective
247: Input Parameter:
248: . matfd - the coloring context
250: Output Parameters:
251: + f - the function
252: - fctx - the optional user-defined function context
254: Level: intermediate
256: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`
257: @*/
258: PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), void **fctx)
259: {
260: PetscFunctionBegin;
262: if (f) *f = matfd->f;
263: if (fctx) *fctx = matfd->fctx;
264: PetscFunctionReturn(PETSC_SUCCESS);
265: }
267: /*@C
268: MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
270: Logically Collective
272: Input Parameters:
273: + matfd - the coloring context
274: . f - the function
275: - fctx - the optional user-defined function context
277: Level: advanced
279: Note:
280: `f` has two possible calling configurations\:
281: .vb
282: PetscErrorCode f(SNES snes, Vec in, Vec out, void *fctx)
283: .ve
284: + snes - the nonlinear solver `SNES` object
285: . in - the location where the Jacobian is to be computed
286: . out - the location to put the computed function value
287: - fctx - the function context
289: and
290: .vb
291: PetscErrorCode f(void *dummy, Vec in, Vec out, void *fctx)
292: .ve
293: + dummy - an unused parameter
294: . in - the location where the Jacobian is to be computed
295: . out - the location to put the computed function value
296: - fctx - the function context
298: This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
299: `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
300: calling `MatFDColoringApply()`
302: Fortran Notes:
303: In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
304: be used without `SNES` or within the `SNES` solvers.
306: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`
307: @*/
308: PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx)
309: {
310: PetscFunctionBegin;
312: matfd->f = f;
313: matfd->fctx = fctx;
314: PetscFunctionReturn(PETSC_SUCCESS);
315: }
317: /*@
318: MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
319: the options database.
321: Collective
323: The Jacobian, F'(u), is estimated with the differencing approximation
324: .vb
325: F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
326: h = error_rel*u[i] if abs(u[i]) > umin
327: = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
328: dx_{i} = (0, ... 1, .... 0)
329: .ve
331: Input Parameter:
332: . matfd - the coloring context
334: Options Database Keys:
335: + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
336: . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
337: . -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
338: . -mat_fd_coloring_view - Activates basic viewing
339: . -mat_fd_coloring_view ::ascii_info - Activates viewing info
340: - -mat_fd_coloring_view draw - Activates drawing
342: Level: intermediate
344: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
345: @*/
346: PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
347: {
348: PetscBool flg;
349: char value[3];
351: PetscFunctionBegin;
354: PetscObjectOptionsBegin((PetscObject)matfd);
355: PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
356: PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
357: PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
358: if (flg) {
359: if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
360: else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
361: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
362: }
363: PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
364: PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
365: if (flg && matfd->bcols > matfd->ncolors) {
366: /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
367: matfd->bcols = matfd->ncolors;
368: }
370: /* process any options handlers added with PetscObjectAddOptionsHandler() */
371: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
372: PetscOptionsEnd();
373: PetscFunctionReturn(PETSC_SUCCESS);
374: }
376: /*@
377: MatFDColoringSetType - Sets the approach for computing the finite difference parameter
379: Collective
381: Input Parameters:
382: + matfd - the coloring context
383: - type - either `MATMFFD_WP` or `MATMFFD_DS`
385: Options Database Key:
386: . -mat_fd_type - "wp" or "ds"
388: Level: intermediate
390: Note:
391: It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
392: but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of
393: introducing another one.
395: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
396: @*/
397: PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
398: {
399: PetscFunctionBegin;
401: /*
402: It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
403: and this function is being provided as patch for a release so it shouldn't change the implementation
404: */
405: if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
406: else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
407: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
408: PetscCall(PetscObjectChangeTypeName((PetscObject)matfd, type));
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
413: {
414: PetscBool flg;
415: PetscViewer viewer;
416: PetscViewerFormat format;
418: PetscFunctionBegin;
419: if (prefix) {
420: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
421: } else {
422: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
423: }
424: if (flg) {
425: PetscCall(PetscViewerPushFormat(viewer, format));
426: PetscCall(MatFDColoringView(fd, viewer));
427: PetscCall(PetscViewerPopFormat(viewer));
428: PetscCall(PetscViewerDestroy(&viewer));
429: }
430: PetscFunctionReturn(PETSC_SUCCESS);
431: }
433: /*@
434: MatFDColoringCreate - Creates a matrix coloring context for finite difference
435: computation of Jacobians.
437: Collective
439: Input Parameters:
440: + mat - the matrix containing the nonzero structure of the Jacobian
441: - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
443: Output Parameter:
444: . color - the new coloring context
446: Level: intermediate
448: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
449: `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
450: `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
451: @*/
452: PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
453: {
454: MatFDColoring c;
455: MPI_Comm comm;
456: PetscInt M, N;
458: PetscFunctionBegin;
460: PetscAssertPointer(color, 3);
461: PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
462: PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
463: PetscCall(MatGetSize(mat, &M, &N));
464: PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
465: PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
466: PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
468: c->ctype = iscoloring->ctype;
469: PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
471: PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
473: PetscCall(MatCreateVecs(mat, NULL, &c->w1));
474: /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
475: PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
476: PetscCall(VecDuplicate(c->w1, &c->w2));
477: /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
478: PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
480: c->error_rel = PETSC_SQRT_MACHINE_EPSILON;
481: c->umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
482: c->currentcolor = -1;
483: c->htype = "wp";
484: c->fset = PETSC_FALSE;
485: c->setupcalled = PETSC_FALSE;
487: *color = c;
488: PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
489: PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
490: PetscFunctionReturn(PETSC_SUCCESS);
491: }
493: /*@
494: MatFDColoringDestroy - Destroys a matrix coloring context that was created
495: via `MatFDColoringCreate()`.
497: Collective
499: Input Parameter:
500: . c - coloring context
502: Level: intermediate
504: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
505: @*/
506: PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
507: {
508: PetscInt i;
509: MatFDColoring color = *c;
511: PetscFunctionBegin;
512: if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
513: if (--((PetscObject)color)->refct > 0) {
514: *c = NULL;
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }
518: /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
519: for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
520: PetscCall(PetscFree(color->isa));
521: PetscCall(PetscFree2(color->ncolumns, color->columns));
522: PetscCall(PetscFree(color->nrows));
523: if (color->htype[0] == 'w') {
524: PetscCall(PetscFree(color->matentry2));
525: } else {
526: PetscCall(PetscFree(color->matentry));
527: }
528: PetscCall(PetscFree(color->dy));
529: if (color->vscale) PetscCall(VecDestroy(&color->vscale));
530: PetscCall(VecDestroy(&color->w1));
531: PetscCall(VecDestroy(&color->w2));
532: PetscCall(VecDestroy(&color->w3));
533: PetscCall(PetscHeaderDestroy(c));
534: PetscFunctionReturn(PETSC_SUCCESS);
535: }
537: /*@C
538: MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
539: that are currently being perturbed.
541: Not Collective
543: Input Parameter:
544: . coloring - coloring context created with `MatFDColoringCreate()`
546: Output Parameters:
547: + n - the number of local columns being perturbed
548: - cols - the column indices, in global numbering
550: Level: advanced
552: Note:
553: IF the matrix type is `MATBAIJ`, then the block column indices are returned
555: Fortran Note:
556: .vb
557: PetscInt, pointer :: cols(:)
558: .ve
559: Use `PETSC_NULL_INTEGER` if `n` is not needed
561: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
562: @*/
563: PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
564: {
565: PetscFunctionBegin;
566: if (coloring->currentcolor >= 0) {
567: *n = coloring->ncolumns[coloring->currentcolor];
568: *cols = coloring->columns[coloring->currentcolor];
569: } else {
570: *n = 0;
571: }
572: PetscFunctionReturn(PETSC_SUCCESS);
573: }
575: /*@
576: MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
577: has been created, computes the Jacobian for a function via finite differences.
579: Collective
581: Input Parameters:
582: + J - matrix to store Jacobian entries into
583: . coloring - coloring context created with `MatFDColoringCreate()`
584: . x1 - location at which Jacobian is to be computed
585: - sctx - context required by function, if this is being used with the `SNES` solver then it is `SNES` object, otherwise it is `NULL`
587: Options Database Keys:
588: + -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
589: . -mat_fd_coloring_view - Activates basic viewing or coloring
590: . -mat_fd_coloring_view draw - Activates drawing of coloring
591: - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
593: Level: intermediate
595: .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
596: @*/
597: PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
598: {
599: PetscBool eq;
601: PetscFunctionBegin;
605: PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
606: PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
607: PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
608: PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
610: PetscCall(MatSetUnfactored(J));
611: PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
612: PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
613: PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
614: if (!coloring->viewed) {
615: PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
616: coloring->viewed = PETSC_TRUE;
617: }
618: PetscFunctionReturn(PETSC_SUCCESS);
619: }