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