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