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