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