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