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 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> 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 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 @*/ 97 PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer) 98 { 99 PetscInt i, j; 100 PetscBool isdraw, iascii; 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, &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)); 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 @*/ 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 @*/ 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 @*/ 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 252 - fctx - the optional user-defined function context 253 254 Level: intermediate 255 256 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()` 257 @*/ 258 PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), 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 275 - fctx - the optional user-defined function context 276 277 Level: advanced 278 279 Note: 280 `f` has two possible calling configurations\: 281 $ PetscErrorCode f(SNES snes, Vec in, Vec out, void *fctx) 282 + snes - the nonlinear solver `SNES` object 283 . in - the location where the Jacobian is to be computed 284 . out - the location to put the computed function value 285 - fctx - the function context 286 287 and 288 $ PetscErrorCode f(void *dummy, Vec in, Vec out, void *fctx) 289 + dummy - an unused parameter 290 . in - the location where the Jacobian is to be computed 291 . out - the location to put the computed function value 292 - fctx - the function context 293 294 This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument 295 `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by 296 calling `MatFDColoringApply()` 297 298 Fortran Notes: 299 In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to 300 be used without `SNES` or within the `SNES` solvers. 301 302 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()` 303 @*/ 304 PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx) 305 { 306 PetscFunctionBegin; 307 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 308 matfd->f = f; 309 matfd->fctx = fctx; 310 PetscFunctionReturn(PETSC_SUCCESS); 311 } 312 313 /*@ 314 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 315 the options database. 316 317 Collective 318 319 The Jacobian, F'(u), is estimated with the differencing approximation 320 .vb 321 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 322 h = error_rel*u[i] if abs(u[i]) > umin 323 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 324 dx_{i} = (0, ... 1, .... 0) 325 .ve 326 327 Input Parameter: 328 . matfd - the coloring context 329 330 Options Database Keys: 331 + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 332 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 333 . -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`) 334 . -mat_fd_coloring_view - Activates basic viewing 335 . -mat_fd_coloring_view ::ascii_info - Activates viewing info 336 - -mat_fd_coloring_view draw - Activates drawing 337 338 Level: intermediate 339 340 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 341 @*/ 342 PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 343 { 344 PetscBool flg; 345 char value[3]; 346 347 PetscFunctionBegin; 348 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 349 350 PetscObjectOptionsBegin((PetscObject)matfd); 351 PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL)); 352 PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL)); 353 PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg)); 354 if (flg) { 355 if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 356 else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 357 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value); 358 } 359 PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL)); 360 PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg)); 361 if (flg && matfd->bcols > matfd->ncolors) { 362 /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */ 363 matfd->bcols = matfd->ncolors; 364 } 365 366 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 367 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject)); 368 PetscOptionsEnd(); 369 PetscFunctionReturn(PETSC_SUCCESS); 370 } 371 372 /*@ 373 MatFDColoringSetType - Sets the approach for computing the finite difference parameter 374 375 Collective 376 377 Input Parameters: 378 + matfd - the coloring context 379 - type - either `MATMFFD_WP` or `MATMFFD_DS` 380 381 Options Database Key: 382 . -mat_fd_type - "wp" or "ds" 383 384 Level: intermediate 385 386 Note: 387 It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries 388 but the process of computing the entries is the same as with the `MATMFFD` operation so we should reuse the names instead of 389 introducing another one. 390 391 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 392 @*/ 393 PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type) 394 { 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 397 /* 398 It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype 399 and this function is being provided as patch for a release so it shouldn't change the implementation 400 */ 401 if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp"; 402 else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds"; 403 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type); 404 PetscCall(PetscObjectChangeTypeName((PetscObject)matfd, type)); 405 PetscFunctionReturn(PETSC_SUCCESS); 406 } 407 408 static PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[]) 409 { 410 PetscBool flg; 411 PetscViewer viewer; 412 PetscViewerFormat format; 413 414 PetscFunctionBegin; 415 if (prefix) { 416 PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg)); 417 } else { 418 PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg)); 419 } 420 if (flg) { 421 PetscCall(PetscViewerPushFormat(viewer, format)); 422 PetscCall(MatFDColoringView(fd, viewer)); 423 PetscCall(PetscViewerPopFormat(viewer)); 424 PetscCall(PetscViewerDestroy(&viewer)); 425 } 426 PetscFunctionReturn(PETSC_SUCCESS); 427 } 428 429 /*@ 430 MatFDColoringCreate - Creates a matrix coloring context for finite difference 431 computation of Jacobians. 432 433 Collective 434 435 Input Parameters: 436 + mat - the matrix containing the nonzero structure of the Jacobian 437 - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()` 438 439 Output Parameter: 440 . color - the new coloring context 441 442 Level: intermediate 443 444 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`, 445 `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`, 446 `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()` 447 @*/ 448 PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color) 449 { 450 MatFDColoring c; 451 MPI_Comm comm; 452 PetscInt M, N; 453 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 456 PetscAssertPointer(color, 3); 457 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();"); 458 PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0)); 459 PetscCall(MatGetSize(mat, &M, &N)); 460 PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices"); 461 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 462 PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView)); 463 464 c->ctype = iscoloring->ctype; 465 PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid)); 466 467 PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c); 468 469 PetscCall(MatCreateVecs(mat, NULL, &c->w1)); 470 /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 471 PetscCall(VecBindToCPU(c->w1, PETSC_TRUE)); 472 PetscCall(VecDuplicate(c->w1, &c->w2)); 473 /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 474 PetscCall(VecBindToCPU(c->w2, PETSC_TRUE)); 475 476 c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 477 c->umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 478 c->currentcolor = -1; 479 c->htype = "wp"; 480 c->fset = PETSC_FALSE; 481 c->setupcalled = PETSC_FALSE; 482 483 *color = c; 484 PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c)); 485 PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0)); 486 PetscFunctionReturn(PETSC_SUCCESS); 487 } 488 489 /*@ 490 MatFDColoringDestroy - Destroys a matrix coloring context that was created 491 via `MatFDColoringCreate()`. 492 493 Collective 494 495 Input Parameter: 496 . c - coloring context 497 498 Level: intermediate 499 500 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()` 501 @*/ 502 PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 503 { 504 PetscInt i; 505 MatFDColoring color = *c; 506 507 PetscFunctionBegin; 508 if (!*c) PetscFunctionReturn(PETSC_SUCCESS); 509 if (--((PetscObject)color)->refct > 0) { 510 *c = NULL; 511 PetscFunctionReturn(PETSC_SUCCESS); 512 } 513 514 /* we do not free the column arrays since their entries are owned by the ISs in color->isa */ 515 for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i])); 516 PetscCall(PetscFree(color->isa)); 517 PetscCall(PetscFree2(color->ncolumns, color->columns)); 518 PetscCall(PetscFree(color->nrows)); 519 if (color->htype[0] == 'w') { 520 PetscCall(PetscFree(color->matentry2)); 521 } else { 522 PetscCall(PetscFree(color->matentry)); 523 } 524 PetscCall(PetscFree(color->dy)); 525 if (color->vscale) PetscCall(VecDestroy(&color->vscale)); 526 PetscCall(VecDestroy(&color->w1)); 527 PetscCall(VecDestroy(&color->w2)); 528 PetscCall(VecDestroy(&color->w3)); 529 PetscCall(PetscHeaderDestroy(c)); 530 PetscFunctionReturn(PETSC_SUCCESS); 531 } 532 533 /*@C 534 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 535 that are currently being perturbed. 536 537 Not Collective 538 539 Input Parameter: 540 . coloring - coloring context created with `MatFDColoringCreate()` 541 542 Output Parameters: 543 + n - the number of local columns being perturbed 544 - cols - the column indices, in global numbering 545 546 Level: advanced 547 548 Note: 549 IF the matrix type is `MATBAIJ`, then the block column indices are returned 550 551 Fortran Note: 552 .vb 553 PetscInt, pointer :: cols(:) 554 .ve 555 Use `PETSC_NULL_INT` if `n` is not needed 556 557 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()` 558 @*/ 559 PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[]) 560 { 561 PetscFunctionBegin; 562 if (coloring->currentcolor >= 0) { 563 *n = coloring->ncolumns[coloring->currentcolor]; 564 *cols = coloring->columns[coloring->currentcolor]; 565 } else { 566 *n = 0; 567 } 568 PetscFunctionReturn(PETSC_SUCCESS); 569 } 570 571 /*@ 572 MatFDColoringApply - Given a matrix for which a `MatFDColoring` context 573 has been created, computes the Jacobian for a function via finite differences. 574 575 Collective 576 577 Input Parameters: 578 + J - matrix to store Jacobian entries into 579 . coloring - coloring context created with `MatFDColoringCreate()` 580 . x1 - location at which Jacobian is to be computed 581 - sctx - context required by function, if this is being used with the `SNES` solver then it is `SNES` object, otherwise it is `NULL` 582 583 Options Database Keys: 584 + -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`) 585 . -mat_fd_coloring_view - Activates basic viewing or coloring 586 . -mat_fd_coloring_view draw - Activates drawing of coloring 587 - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 588 589 Level: intermediate 590 591 .seealso: `Mat`, `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()` 592 @*/ 593 PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx) 594 { 595 PetscBool eq; 596 597 PetscFunctionBegin; 598 PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 599 PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2); 600 PetscValidHeaderSpecific(x1, VEC_CLASSID, 3); 601 PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq)); 602 PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()"); 603 PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()"); 604 PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()"); 605 606 PetscCall(MatSetUnfactored(J)); 607 PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0)); 608 PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx); 609 PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0)); 610 if (!coloring->viewed) { 611 PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view")); 612 coloring->viewed = PETSC_TRUE; 613 } 614 PetscFunctionReturn(PETSC_SUCCESS); 615 } 616