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