1 /*
2 This is where the abstract matrix operations are defined that are
3 used for finite difference computations of Jacobians using coloring.
4 */
5
6 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
7 #include <petsc/private/isimpl.h>
8
MatFDColoringSetF(MatFDColoring fd,Vec F)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 }
20
21 #include <petscdraw.h>
MatFDColoringView_Draw_Zoom(PetscDraw draw,void * Aa)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;
29
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 }
43
MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)44 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer)
45 {
46 PetscBool isnull;
47 PetscDraw draw;
48 PetscReal xr, yr, xl, yl, h, w;
49
50 PetscFunctionBegin;
51 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
52 PetscCall(PetscDrawIsNull(draw, &isnull));
53 if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
54
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 }
70
71 /*@
72 MatFDColoringView - Views a finite difference coloring context.
73
74 Collective
75
76 Input Parameters:
77 + c - the coloring context
78 - viewer - visualization context
79
80 Level: intermediate
81
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
90
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.
94
95 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
96 @*/
MatFDColoringView(MatFDColoring c,PetscViewer viewer)97 PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer)
98 {
99 PetscInt i, j;
100 PetscBool isdraw, isascii;
101 PetscViewerFormat format;
102
103 PetscFunctionBegin;
104 PetscValidHeaderSpecific(c, MAT_FDCOLORING_CLASSID, 1);
105 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
106 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
107 PetscCheckSameComm(c, 1, viewer, 2);
108
109 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
110 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
111 if (isdraw) {
112 PetscCall(MatFDColoringView_Draw(c, viewer));
113 } else if (isascii) {
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));
118
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 }
141
142 /*@
143 MatFDColoringSetParameters - Sets the parameters for the approximation of
144 a sparse Jacobian matrix using finite differences and matrix coloring
145
146 Logically Collective
147
148 Input Parameters:
149 + matfd - the coloring context
150 . error - relative error
151 - umin - minimum allowable u-value magnitude
152
153 Level: advanced
154
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)
163
164 htype = 'wp':
165 h = error_rel * sqrt(1 + ||u||)
166 .ve
167
168 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
169 @*/
MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)170 PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin)
171 {
172 PetscFunctionBegin;
173 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
174 PetscValidLogicalCollectiveReal(matfd, error, 2);
175 PetscValidLogicalCollectiveReal(matfd, umin, 3);
176 if (error != (PetscReal)PETSC_DEFAULT) matfd->error_rel = error;
177 if (umin != (PetscReal)PETSC_DEFAULT) matfd->umin = umin;
178 PetscFunctionReturn(PETSC_SUCCESS);
179 }
180
181 /*@
182 MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix.
183
184 Logically Collective
185
186 Input Parameters:
187 + matfd - the coloring context
188 . brows - number of rows in the block
189 - bcols - number of columns in the block
190
191 Level: intermediate
192
193 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()`
194 @*/
MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols)195 PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols)
196 {
197 PetscFunctionBegin;
198 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
199 PetscValidLogicalCollectiveInt(matfd, brows, 2);
200 PetscValidLogicalCollectiveInt(matfd, bcols, 3);
201 if (brows != PETSC_DEFAULT) matfd->brows = brows;
202 if (bcols != PETSC_DEFAULT) matfd->bcols = bcols;
203 PetscFunctionReturn(PETSC_SUCCESS);
204 }
205
206 /*@
207 MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use.
208
209 Collective
210
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
215
216 Level: beginner
217
218 Notes:
219 When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns.
220
221 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`
222 @*/
MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color)223 PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color)
224 {
225 PetscBool eq;
226
227 PetscFunctionBegin;
228 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
229 PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3);
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()");
233
234 PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0));
235 PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color);
236
237 color->setupcalled = PETSC_TRUE;
238 PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0));
239 PetscFunctionReturn(PETSC_SUCCESS);
240 }
241
242 /*@C
243 MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
244
245 Not Collective
246
247 Input Parameter:
248 . matfd - the coloring context
249
250 Output Parameters:
251 + f - the function, see `MatFDColoringFn` for the calling sequence
252 - fctx - the optional user-defined function context
253
254 Level: intermediate
255
256 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
257 @*/
MatFDColoringGetFunction(MatFDColoring matfd,MatFDColoringFn ** f,void ** fctx)258 PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, MatFDColoringFn **f, void **fctx)
259 {
260 PetscFunctionBegin;
261 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
262 if (f) *f = matfd->f;
263 if (fctx) *fctx = matfd->fctx;
264 PetscFunctionReturn(PETSC_SUCCESS);
265 }
266
267 /*@C
268 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
269
270 Logically Collective
271
272 Input Parameters:
273 + matfd - the coloring context
274 . f - the function, see `MatFDColoringFn` for the calling sequence
275 - fctx - the optional user-defined function context
276
277 Level: advanced
278
279 Note:
280 This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument
281 `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by
282 calling `MatFDColoringApply()`
283
284 Fortran Note:
285 In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to
286 be used without `SNES` or within the `SNES` solvers.
287
288 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringFn`
289 @*/
MatFDColoringSetFunction(MatFDColoring matfd,MatFDColoringFn * f,void * fctx)290 PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, MatFDColoringFn *f, void *fctx)
291 {
292 PetscFunctionBegin;
293 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
294 matfd->f = f;
295 matfd->fctx = fctx;
296 PetscFunctionReturn(PETSC_SUCCESS);
297 }
298
299 /*@
300 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
301 the options database.
302
303 Collective
304
305 The Jacobian, F'(u), is estimated with the differencing approximation
306 .vb
307 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
308 h = error_rel*u[i] if abs(u[i]) > umin
309 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i]
310 dx_{i} = (0, ... 1, .... 0)
311 .ve
312
313 Input Parameter:
314 . matfd - the coloring context
315
316 Options Database Keys:
317 + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function)
318 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
319 . -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
320 . -mat_fd_coloring_view - Activates basic viewing
321 . -mat_fd_coloring_view ::ascii_info - Activates viewing info
322 - -mat_fd_coloring_view draw - Activates drawing
323
324 Level: intermediate
325
326 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
327 @*/
MatFDColoringSetFromOptions(MatFDColoring matfd)328 PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd)
329 {
330 PetscBool flg;
331 char value[3];
332
333 PetscFunctionBegin;
334 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
335
336 PetscObjectOptionsBegin((PetscObject)matfd);
337 PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL));
338 PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL));
339 PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg));
340 if (flg) {
341 if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
342 else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
343 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value);
344 }
345 PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL));
346 PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg));
347 if (flg && matfd->bcols > matfd->ncolors) {
348 /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */
349 matfd->bcols = matfd->ncolors;
350 }
351
352 /* process any options handlers added with PetscObjectAddOptionsHandler() */
353 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject));
354 PetscOptionsEnd();
355 PetscFunctionReturn(PETSC_SUCCESS);
356 }
357
358 /*@
359 MatFDColoringSetType - Sets the approach for computing the finite difference parameter
360
361 Collective
362
363 Input Parameters:
364 + matfd - the coloring context
365 - type - either `MATMFFD_WP` or `MATMFFD_DS`
366
367 Options Database Key:
368 . -mat_fd_type - "wp" or "ds"
369
370 Level: intermediate
371
372 Note:
373 It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries
374 but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of
375 introducing another one.
376
377 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()`
378 @*/
MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type)379 PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type)
380 {
381 PetscFunctionBegin;
382 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1);
383 /*
384 It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype
385 and this function is being provided as patch for a release so it shouldn't change the implementation
386 */
387 if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp";
388 else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds";
389 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type);
390 PetscCall(PetscObjectChangeTypeName((PetscObject)matfd, type));
391 PetscFunctionReturn(PETSC_SUCCESS);
392 }
393
MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[])394 static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[])
395 {
396 PetscBool flg;
397 PetscViewer viewer;
398 PetscViewerFormat format;
399
400 PetscFunctionBegin;
401 if (prefix) {
402 PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg));
403 } else {
404 PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg));
405 }
406 if (flg) {
407 PetscCall(PetscViewerPushFormat(viewer, format));
408 PetscCall(MatFDColoringView(fd, viewer));
409 PetscCall(PetscViewerPopFormat(viewer));
410 PetscCall(PetscViewerDestroy(&viewer));
411 }
412 PetscFunctionReturn(PETSC_SUCCESS);
413 }
414
415 /*@
416 MatFDColoringCreate - Creates a matrix coloring context for finite difference
417 computation of Jacobians.
418
419 Collective
420
421 Input Parameters:
422 + mat - the matrix containing the nonzero structure of the Jacobian
423 - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`
424
425 Output Parameter:
426 . color - the new coloring context
427
428 Level: intermediate
429
430 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`,
431 `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`,
432 `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()`
433 @*/
MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring * color)434 PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color)
435 {
436 MatFDColoring c;
437 MPI_Comm comm;
438 PetscInt M, N;
439
440 PetscFunctionBegin;
441 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
442 PetscAssertPointer(color, 3);
443 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();");
444 PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0));
445 PetscCall(MatGetSize(mat, &M, &N));
446 PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices");
447 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
448 PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView));
449
450 c->ctype = iscoloring->ctype;
451 PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid));
452
453 PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c);
454
455 PetscCall(MatCreateVecs(mat, NULL, &c->w1));
456 /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
457 PetscCall(VecBindToCPU(c->w1, PETSC_TRUE));
458 PetscCall(VecDuplicate(c->w1, &c->w2));
459 /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
460 PetscCall(VecBindToCPU(c->w2, PETSC_TRUE));
461
462 c->error_rel = PETSC_SQRT_MACHINE_EPSILON;
463 c->umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
464 c->currentcolor = -1;
465 c->htype = "wp";
466 c->fset = PETSC_FALSE;
467 c->setupcalled = PETSC_FALSE;
468
469 *color = c;
470 PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c));
471 PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0));
472 PetscFunctionReturn(PETSC_SUCCESS);
473 }
474
475 /*@
476 MatFDColoringDestroy - Destroys a matrix coloring context that was created
477 via `MatFDColoringCreate()`.
478
479 Collective
480
481 Input Parameter:
482 . c - coloring context
483
484 Level: intermediate
485
486 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`
487 @*/
MatFDColoringDestroy(MatFDColoring * c)488 PetscErrorCode MatFDColoringDestroy(MatFDColoring *c)
489 {
490 PetscInt i;
491 MatFDColoring color = *c;
492
493 PetscFunctionBegin;
494 if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
495 if (--((PetscObject)color)->refct > 0) {
496 *c = NULL;
497 PetscFunctionReturn(PETSC_SUCCESS);
498 }
499
500 /* we do not free the column arrays since their entries are owned by the ISs in color->isa */
501 for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i]));
502 PetscCall(PetscFree(color->isa));
503 PetscCall(PetscFree2(color->ncolumns, color->columns));
504 PetscCall(PetscFree(color->nrows));
505 if (color->htype[0] == 'w') {
506 PetscCall(PetscFree(color->matentry2));
507 } else {
508 PetscCall(PetscFree(color->matentry));
509 }
510 PetscCall(PetscFree(color->dy));
511 if (color->vscale) PetscCall(VecDestroy(&color->vscale));
512 PetscCall(VecDestroy(&color->w1));
513 PetscCall(VecDestroy(&color->w2));
514 PetscCall(VecDestroy(&color->w3));
515 PetscCall(PetscHeaderDestroy(c));
516 PetscFunctionReturn(PETSC_SUCCESS);
517 }
518
519 /*@C
520 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
521 that are currently being perturbed.
522
523 Not Collective
524
525 Input Parameter:
526 . coloring - coloring context created with `MatFDColoringCreate()`
527
528 Output Parameters:
529 + n - the number of local columns being perturbed
530 - cols - the column indices, in global numbering
531
532 Level: advanced
533
534 Note:
535 IF the matrix type is `MATBAIJ`, then the block column indices are returned
536
537 Fortran Note:
538 .vb
539 PetscInt, pointer :: cols(:)
540 .ve
541 Use `PETSC_NULL_INTEGER` if `n` is not needed
542
543 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()`
544 @*/
MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt * n,const PetscInt * cols[])545 PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[])
546 {
547 PetscFunctionBegin;
548 if (coloring->currentcolor >= 0) {
549 *n = coloring->ncolumns[coloring->currentcolor];
550 *cols = coloring->columns[coloring->currentcolor];
551 } else {
552 *n = 0;
553 }
554 PetscFunctionReturn(PETSC_SUCCESS);
555 }
556
557 /*@
558 MatFDColoringApply - Given a matrix for which a `MatFDColoring` context
559 has been created, computes the Jacobian for a function via finite differences.
560
561 Collective
562
563 Input Parameters:
564 + J - matrix to store Jacobian entries into
565 . coloring - coloring context created with `MatFDColoringCreate()`
566 . x1 - location at which Jacobian is to be computed
567 - sctx - context required by function, if this is being used with the `SNES` solver then it is `SNES` object, otherwise it is `NULL`
568
569 Options Database Keys:
570 + -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`)
571 . -mat_fd_coloring_view - Activates basic viewing or coloring
572 . -mat_fd_coloring_view draw - Activates drawing of coloring
573 - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info
574
575 Level: intermediate
576
577 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()`
578 @*/
MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void * sctx)579 PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx)
580 {
581 PetscBool eq;
582
583 PetscFunctionBegin;
584 PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
585 PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2);
586 PetscValidHeaderSpecific(x1, VEC_CLASSID, 3);
587 PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq));
588 PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()");
589 PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()");
590 PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()");
591
592 PetscCall(MatSetUnfactored(J));
593 PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0));
594 PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx);
595 PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0));
596 if (!coloring->viewed) {
597 PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view"));
598 coloring->viewed = PETSC_TRUE;
599 }
600 PetscFunctionReturn(PETSC_SUCCESS);
601 }
602