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