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