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 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatFDColoringSetF" 11 PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F) 12 { 13 PetscErrorCode ierr; 14 15 PetscFunctionBegin; 16 if (F) { 17 ierr = VecCopy(F,fd->w1);CHKERRQ(ierr); 18 fd->fset = PETSC_TRUE; 19 } else { 20 fd->fset = PETSC_FALSE; 21 } 22 PetscFunctionReturn(0); 23 } 24 25 #include <petscdraw.h> 26 #undef __FUNCT__ 27 #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 28 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 29 { 30 MatFDColoring fd = (MatFDColoring)Aa; 31 PetscErrorCode ierr; 32 PetscInt i,j,nz,row; 33 PetscReal x,y; 34 MatEntry *Jentry=fd->matentry; 35 36 PetscFunctionBegin; 37 /* loop over colors */ 38 nz = 0; 39 for (i=0; i<fd->ncolors; i++) { 40 for (j=0; j<fd->nrows[i]; j++) { 41 row = Jentry[nz].row; 42 y = fd->M - row - fd->rstart; 43 x = (PetscReal)Jentry[nz++].col; 44 ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 45 } 46 } 47 PetscFunctionReturn(0); 48 } 49 50 #undef __FUNCT__ 51 #define __FUNCT__ "MatFDColoringView_Draw" 52 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 53 { 54 PetscErrorCode ierr; 55 PetscBool isnull; 56 PetscDraw draw; 57 PetscReal xr,yr,xl,yl,h,w; 58 59 PetscFunctionBegin; 60 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 61 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 62 63 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 64 65 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 66 xr += w; yr += h; xl = -w; yl = -h; 67 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 68 ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 69 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr); 70 PetscFunctionReturn(0); 71 } 72 73 #undef __FUNCT__ 74 #define __FUNCT__ "MatFDColoringView" 75 /*@C 76 MatFDColoringView - Views a finite difference coloring context. 77 78 Collective on MatFDColoring 79 80 Input Parameters: 81 + c - the coloring context 82 - viewer - visualization context 83 84 Level: intermediate 85 86 Notes: 87 The available visualization contexts include 88 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 89 . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 90 output where only the first processor opens 91 the file. All other processors send their 92 data to the first processor to print. 93 - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 94 95 Notes: 96 Since PETSc uses only a small number of basic colors (currently 33), if the coloring 97 involves more than 33 then some seemingly identical colors are displayed making it look 98 like an illegal coloring. This is just a graphical artifact. 99 100 .seealso: MatFDColoringCreate() 101 102 .keywords: Mat, finite differences, coloring, view 103 @*/ 104 PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 105 { 106 PetscErrorCode ierr; 107 PetscInt i,j; 108 PetscBool isdraw,iascii; 109 PetscViewerFormat format; 110 111 PetscFunctionBegin; 112 PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 113 if (!viewer) { 114 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr); 115 } 116 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 117 PetscCheckSameComm(c,1,viewer,2); 118 119 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 120 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 121 if (isdraw) { 122 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 123 } else if (iascii) { 124 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr); 125 ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr); 126 ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",(double)c->umin);CHKERRQ(ierr); 127 ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 128 129 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 130 if (format != PETSC_VIEWER_ASCII_INFO) { 131 PetscInt row,col,nz; 132 nz = 0; 133 for (i=0; i<c->ncolors; i++) { 134 ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 135 ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 136 for (j=0; j<c->ncolumns[i]; j++) { 137 ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 138 } 139 ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 140 for (j=0; j<c->nrows[i]; j++) { 141 row = c->matentry[nz].row; 142 col = c->matentry[nz++].col; 143 ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",row,col);CHKERRQ(ierr); 144 } 145 } 146 } 147 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 148 } 149 PetscFunctionReturn(0); 150 } 151 152 #undef __FUNCT__ 153 #define __FUNCT__ "MatFDColoringSetParameters" 154 /*@ 155 MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 156 a Jacobian matrix using finite differences. 157 158 Logically Collective on MatFDColoring 159 160 The Jacobian is estimated with the differencing approximation 161 .vb 162 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 163 h = error_rel*u[i] if abs(u[i]) > umin 164 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 165 dx_{i} = (0, ... 1, .... 0) 166 .ve 167 168 Input Parameters: 169 + coloring - the coloring context 170 . error_rel - relative error 171 - umin - minimum allowable u-value magnitude 172 173 Level: advanced 174 175 .keywords: Mat, finite differences, coloring, set, parameters 176 177 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 178 179 @*/ 180 PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 181 { 182 PetscFunctionBegin; 183 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 184 PetscValidLogicalCollectiveReal(matfd,error,2); 185 PetscValidLogicalCollectiveReal(matfd,umin,3); 186 if (error != PETSC_DEFAULT) matfd->error_rel = error; 187 if (umin != PETSC_DEFAULT) matfd->umin = umin; 188 PetscFunctionReturn(0); 189 } 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "MatFDColoringSetBlockSize" 193 /*@ 194 MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix. 195 196 Logically Collective on MatFDColoring 197 198 Input Parameters: 199 + coloring - the coloring context 200 . brows - number of rows in the block 201 - bcols - number of columns in the block 202 203 Level: intermediate 204 205 .keywords: Mat, coloring 206 207 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 208 209 @*/ 210 PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols) 211 { 212 PetscFunctionBegin; 213 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 214 PetscValidLogicalCollectiveInt(matfd,brows,2); 215 PetscValidLogicalCollectiveInt(matfd,bcols,3); 216 if (brows != PETSC_DEFAULT) matfd->brows = brows; 217 if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "MatFDColoringSetUp" 223 /*@ 224 MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 225 226 Collective on Mat 227 228 Input Parameters: 229 + mat - the matrix containing the nonzero structure of the Jacobian 230 . iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 231 - color - the matrix coloring context 232 233 Level: beginner 234 235 .keywords: MatFDColoring, setup 236 237 .seealso: MatFDColoringCreate(), MatFDColoringDestroy() 238 @*/ 239 PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color) 240 { 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 245 PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3); 246 if (color->setupcalled) PetscFunctionReturn(0); 247 248 ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 249 if (mat->ops->fdcoloringsetup) { 250 ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr); 251 } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 252 253 color->setupcalled = PETSC_TRUE; 254 ierr = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 255 PetscFunctionReturn(0); 256 } 257 258 #undef __FUNCT__ 259 #define __FUNCT__ "MatFDColoringGetFunction" 260 /*@C 261 MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 262 263 Not Collective 264 265 Input Parameters: 266 . coloring - the coloring context 267 268 Output Parameters: 269 + f - the function 270 - fctx - the optional user-defined function context 271 272 Level: intermediate 273 274 .keywords: Mat, Jacobian, finite differences, set, function 275 276 .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 277 278 @*/ 279 PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 280 { 281 PetscFunctionBegin; 282 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 283 if (f) *f = matfd->f; 284 if (fctx) *fctx = matfd->fctx; 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "MatFDColoringSetFunction" 290 /*@C 291 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 292 293 Logically Collective on MatFDColoring 294 295 Input Parameters: 296 + coloring - the coloring context 297 . f - the function 298 - fctx - the optional user-defined function context 299 300 Calling sequence of (*f) function: 301 For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 302 If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 303 304 Level: advanced 305 306 Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 307 SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by 308 calling MatFDColoringApply() 309 310 Fortran Notes: 311 In Fortran you must call MatFDColoringSetFunction() for a coloring object to 312 be used without SNES or within the SNES solvers. 313 314 .keywords: Mat, Jacobian, finite differences, set, function 315 316 .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 317 318 @*/ 319 PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 320 { 321 PetscFunctionBegin; 322 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 323 matfd->f = f; 324 matfd->fctx = fctx; 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "MatFDColoringSetFromOptions" 330 /*@ 331 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 332 the options database. 333 334 Collective on MatFDColoring 335 336 The Jacobian, F'(u), is estimated with the differencing approximation 337 .vb 338 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 339 h = error_rel*u[i] if abs(u[i]) > umin 340 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 341 dx_{i} = (0, ... 1, .... 0) 342 .ve 343 344 Input Parameter: 345 . coloring - the coloring context 346 347 Options Database Keys: 348 + -mat_fd_coloring_err <err> - Sets <err> (square root 349 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,_p_MatFDColoring,int,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 = MatGetVecs(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,MatStructure *flag,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,flag,sctx);CHKERRQ(ierr); 610 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 611 PetscFunctionReturn(0); 612 } 613