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