1 #define PETSCMAT_DLL 2 3 /* 4 This is where the abstract matrix operations are defined that are 5 used for finite difference computations of Jacobians using coloring. 6 */ 7 8 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 9 10 /* Logging support */ 11 PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE = 0; 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatFDColoringSetF" 15 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F) 16 { 17 PetscFunctionBegin; 18 fd->F = F; 19 PetscFunctionReturn(0); 20 } 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 24 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 25 { 26 MatFDColoring fd = (MatFDColoring)Aa; 27 PetscErrorCode ierr; 28 PetscInt i,j; 29 PetscReal x,y; 30 31 PetscFunctionBegin; 32 33 /* loop over colors */ 34 for (i=0; i<fd->ncolors; i++) { 35 for (j=0; j<fd->nrows[i]; j++) { 36 y = fd->M - fd->rows[i][j] - fd->rstart; 37 x = fd->columnsforrow[i][j]; 38 ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 39 } 40 } 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "MatFDColoringView_Draw" 46 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 47 { 48 PetscErrorCode ierr; 49 PetscTruth isnull; 50 PetscDraw draw; 51 PetscReal xr,yr,xl,yl,h,w; 52 53 PetscFunctionBegin; 54 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 55 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 56 57 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 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 = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 63 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 #undef __FUNCT__ 68 #define __FUNCT__ "MatFDColoringView" 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 PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer) 99 { 100 PetscErrorCode ierr; 101 PetscInt i,j; 102 PetscTruth isdraw,iascii; 103 PetscViewerFormat format; 104 105 PetscFunctionBegin; 106 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1); 107 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 108 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 109 PetscCheckSameComm(c,1,viewer,2); 110 111 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 112 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 113 if (isdraw) { 114 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 115 } else if (iascii) { 116 ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 117 ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 118 ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 119 ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 120 121 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 122 if (format != PETSC_VIEWER_ASCII_INFO) { 123 for (i=0; i<c->ncolors; i++) { 124 ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 125 ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 126 for (j=0; j<c->ncolumns[i]; j++) { 127 ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 128 } 129 ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 130 for (j=0; j<c->nrows[i]; j++) { 131 ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 132 } 133 } 134 } 135 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 136 } else { 137 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 138 } 139 PetscFunctionReturn(0); 140 } 141 142 #undef __FUNCT__ 143 #define __FUNCT__ "MatFDColoringSetParameters" 144 /*@ 145 MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 146 a Jacobian matrix using finite differences. 147 148 Collective on MatFDColoring 149 150 The Jacobian is estimated with the differencing approximation 151 .vb 152 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 153 h = error_rel*u[i] if abs(u[i]) > umin 154 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 155 dx_{i} = (0, ... 1, .... 0) 156 .ve 157 158 Input Parameters: 159 + coloring - the coloring context 160 . error_rel - relative error 161 - umin - minimum allowable u-value magnitude 162 163 Level: advanced 164 165 .keywords: Mat, finite differences, coloring, set, parameters 166 167 .seealso: MatFDColoringCreate() 168 @*/ 169 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 170 { 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 173 174 if (error != PETSC_DEFAULT) matfd->error_rel = error; 175 if (umin != PETSC_DEFAULT) matfd->umin = umin; 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "MatFDColoringSetFrequency" 181 /*@ 182 MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 183 matrices. 184 185 Collective on MatFDColoring 186 187 Input Parameters: 188 + coloring - the coloring context 189 - freq - frequency (default is 1) 190 191 Options Database Keys: 192 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 193 194 Level: advanced 195 196 Notes: 197 Using a modified Newton strategy, where the Jacobian remains fixed for several 198 iterations, can be cost effective in terms of overall nonlinear solution 199 efficiency. This parameter indicates that a new Jacobian will be computed every 200 <freq> nonlinear iterations. 201 202 .keywords: Mat, finite differences, coloring, set, frequency 203 204 .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute() 205 @*/ 206 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq) 207 { 208 PetscFunctionBegin; 209 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 210 211 matfd->freq = freq; 212 PetscFunctionReturn(0); 213 } 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "MatFDColoringGetFrequency" 217 /*@ 218 MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 219 matrices. 220 221 Not Collective 222 223 Input Parameters: 224 . coloring - the coloring context 225 226 Output Parameters: 227 . freq - frequency (default is 1) 228 229 Options Database Keys: 230 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 231 232 Level: advanced 233 234 Notes: 235 Using a modified Newton strategy, where the Jacobian remains fixed for several 236 iterations, can be cost effective in terms of overall nonlinear solution 237 efficiency. This parameter indicates that a new Jacobian will be computed every 238 <freq> nonlinear iterations. 239 240 .keywords: Mat, finite differences, coloring, get, frequency 241 242 .seealso: MatFDColoringSetFrequency() 243 @*/ 244 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq) 245 { 246 PetscFunctionBegin; 247 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 248 *freq = matfd->freq; 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "MatFDColoringGetFunction" 254 /*@C 255 MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 256 257 Collective on MatFDColoring 258 259 Input Parameters: 260 . coloring - the coloring context 261 262 Output Parameters: 263 + f - the function 264 - fctx - the optional user-defined function context 265 266 Level: intermediate 267 268 .keywords: Mat, Jacobian, finite differences, set, function 269 @*/ 270 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 271 { 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 274 if (f) *f = matfd->f; 275 if (fctx) *fctx = matfd->fctx; 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "MatFDColoringSetFunction" 281 /*@C 282 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 283 284 Collective on MatFDColoring 285 286 Input Parameters: 287 + coloring - the coloring context 288 . f - the function 289 - fctx - the optional user-defined function context 290 291 Level: intermediate 292 293 Notes: 294 In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 295 be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 296 with the TS solvers. 297 298 .keywords: Mat, Jacobian, finite differences, set, function 299 @*/ 300 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 301 { 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 304 matfd->f = f; 305 matfd->fctx = fctx; 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "MatFDColoringSetFromOptions" 311 /*@ 312 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 313 the options database. 314 315 Collective on MatFDColoring 316 317 The Jacobian, F'(u), is estimated with the differencing approximation 318 .vb 319 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 320 h = error_rel*u[i] if abs(u[i]) > umin 321 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 322 dx_{i} = (0, ... 1, .... 0) 323 .ve 324 325 Input Parameter: 326 . coloring - the coloring context 327 328 Options Database Keys: 329 + -mat_fd_coloring_err <err> - Sets <err> (square root 330 of relative error in the function) 331 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 332 . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 333 . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS) 334 . -mat_fd_coloring_view - Activates basic viewing 335 . -mat_fd_coloring_view_info - Activates viewing info 336 - -mat_fd_coloring_view_draw - Activates drawing 337 338 Level: intermediate 339 340 .keywords: Mat, finite differences, parameters 341 342 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 343 344 @*/ 345 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd) 346 { 347 PetscErrorCode ierr; 348 PetscTruth flg; 349 char value[3]; 350 351 PetscFunctionBegin; 352 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 353 354 ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 355 ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 356 ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 357 ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 358 ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&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_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 363 } 364 /* not used here; but so they are presented in the GUI */ 365 ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 366 ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 367 ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 368 PetscOptionsEnd();CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372 #undef __FUNCT__ 373 #define __FUNCT__ "MatFDColoringView_Private" 374 PetscErrorCode MatFDColoringView_Private(MatFDColoring fd) 375 { 376 PetscErrorCode ierr; 377 PetscTruth flg; 378 379 PetscFunctionBegin; 380 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 381 if (flg) { 382 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 383 } 384 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 385 if (flg) { 386 ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 387 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 388 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 389 } 390 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 391 if (flg) { 392 ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 393 ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 394 } 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "MatFDColoringCreate" 400 /*@ 401 MatFDColoringCreate - Creates a matrix coloring context for finite difference 402 computation of Jacobians. 403 404 Collective on Mat 405 406 Input Parameters: 407 + mat - the matrix containing the nonzero structure of the Jacobian 408 - iscoloring - the coloring of the matrix 409 410 Output Parameter: 411 . color - the new coloring context 412 413 Level: intermediate 414 415 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 416 MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 417 MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(), 418 MatFDColoringSetParameters() 419 @*/ 420 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 421 { 422 MatFDColoring c; 423 MPI_Comm comm; 424 PetscErrorCode ierr; 425 PetscInt M,N; 426 PetscMPIInt size; 427 428 PetscFunctionBegin; 429 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 430 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 431 if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 432 433 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 434 ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 435 436 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 437 c->ctype = iscoloring->ctype; 438 if (size == 1) c->ctype = iscoloring->ctype = IS_COLORING_LOCAL; 439 440 if (mat->ops->fdcoloringcreate) { 441 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 442 } else { 443 SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 444 } 445 446 ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr); 447 ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr); 448 ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 449 ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr); 450 451 c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 452 c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 453 c->freq = 1; 454 c->usersetsrecompute = PETSC_FALSE; 455 c->recompute = PETSC_FALSE; 456 c->currentcolor = -1; 457 c->htype = "wp"; 458 459 *color = c; 460 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 461 PetscFunctionReturn(0); 462 } 463 464 #undef __FUNCT__ 465 #define __FUNCT__ "MatFDColoringDestroy" 466 /*@ 467 MatFDColoringDestroy - Destroys a matrix coloring context that was created 468 via MatFDColoringCreate(). 469 470 Collective on MatFDColoring 471 472 Input Parameter: 473 . c - coloring context 474 475 Level: intermediate 476 477 .seealso: MatFDColoringCreate() 478 @*/ 479 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c) 480 { 481 PetscErrorCode ierr; 482 PetscInt i; 483 484 PetscFunctionBegin; 485 if (--c->refct > 0) PetscFunctionReturn(0); 486 487 for (i=0; i<c->ncolors; i++) { 488 ierr = PetscFree(c->columns[i]);CHKERRQ(ierr); 489 ierr = PetscFree(c->rows[i]);CHKERRQ(ierr); 490 ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr); 491 if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 492 } 493 ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 494 ierr = PetscFree(c->columns);CHKERRQ(ierr); 495 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 496 ierr = PetscFree(c->rows);CHKERRQ(ierr); 497 ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 498 ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr); 499 if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 500 if (c->w1) { 501 ierr = VecDestroy(c->w1);CHKERRQ(ierr); 502 ierr = VecDestroy(c->w2);CHKERRQ(ierr); 503 } 504 if (c->w3){ 505 ierr = VecDestroy(c->w3);CHKERRQ(ierr); 506 } 507 ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 508 PetscFunctionReturn(0); 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 513 /*@C 514 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 515 that are currently being perturbed. 516 517 Not Collective 518 519 Input Parameters: 520 . coloring - coloring context created with MatFDColoringCreate() 521 522 Output Parameters: 523 + n - the number of local columns being perturbed 524 - cols - the column indices, in global numbering 525 526 Level: intermediate 527 528 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 529 530 .keywords: coloring, Jacobian, finite differences 531 @*/ 532 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 533 { 534 PetscFunctionBegin; 535 if (coloring->currentcolor >= 0) { 536 *n = coloring->ncolumns[coloring->currentcolor]; 537 *cols = coloring->columns[coloring->currentcolor]; 538 } else { 539 *n = 0; 540 } 541 PetscFunctionReturn(0); 542 } 543 544 #include "petscda.h" /*I "petscda.h" I*/ 545 546 #undef __FUNCT__ 547 #define __FUNCT__ "MatFDColoringApply" 548 /*@ 549 MatFDColoringApply - Given a matrix for which a MatFDColoring context 550 has been created, computes the Jacobian for a function via finite differences. 551 552 Collective on MatFDColoring 553 554 Input Parameters: 555 + mat - location to store Jacobian 556 . coloring - coloring context created with MatFDColoringCreate() 557 . x1 - location at which Jacobian is to be computed 558 - sctx - optional context required by function (actually a SNES context) 559 560 Options Database Keys: 561 + -mat_fd_coloring_freq <freq> - Sets coloring frequency 562 . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS) 563 . -mat_fd_coloring_view - Activates basic viewing or coloring 564 . -mat_fd_coloring_view_draw - Activates drawing of coloring 565 - -mat_fd_coloring_view_info - Activates viewing of coloring info 566 567 Level: intermediate 568 569 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 570 571 .keywords: coloring, Jacobian, finite differences 572 @*/ 573 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 574 { 575 PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; 576 PetscErrorCode ierr; 577 PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 578 PetscScalar dx,*y,*xx,*w3_array; 579 PetscScalar *vscale_array; 580 PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 581 Vec w1=coloring->w1,w2=coloring->w2,w3; 582 void *fctx = coloring->fctx; 583 PetscTruth flg; 584 PetscInt ctype=coloring->ctype,N,col_start,col_end; 585 Vec x1_tmp; 586 /* remove ! */ 587 PetscMPIInt rank; 588 PetscInt prid=10; 589 /* ex5 590 PetscTruth fd_jacobian_ghost=PETSC_FALSE; 591 DA da; 592 */ 593 594 PetscFunctionBegin; 595 596 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 597 ierr = PetscOptionsGetInt(PETSC_NULL,"-prid",&prid,PETSC_NULL);CHKERRQ(ierr); 598 599 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 600 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 601 PetscValidHeaderSpecific(x1,VEC_COOKIE,3); 602 if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 603 604 if (coloring->usersetsrecompute) { 605 if (!coloring->recompute) { 606 *flag = SAME_PRECONDITIONER; 607 ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr); 608 PetscFunctionReturn(0); 609 } else { 610 coloring->recompute = PETSC_FALSE; 611 } 612 } 613 614 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 615 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 616 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 617 if (flg) { 618 ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 619 } else { 620 PetscTruth assembled; 621 ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 622 if (assembled) { 623 ierr = MatZeroEntries(J);CHKERRQ(ierr); 624 } 625 } 626 627 x1_tmp = x1; 628 629 if (ctype == IS_COLORING_GHOSTED && !coloring->vscale){ 630 ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 631 } 632 633 /* 634 This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets 635 coloring->F for the coarser grids from the finest 636 */ 637 if (coloring->F) { 638 ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 639 ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 640 if (m1 != m2) { 641 coloring->F = 0; 642 } 643 } 644 645 if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 646 ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 647 } 648 ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 649 650 /* Set w1 = F(x1) */ 651 if (coloring->F) { 652 w1 = coloring->F; /* use already computed value of function */ 653 coloring->F = 0; 654 } else { 655 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 656 ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 657 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 658 } 659 660 if (!coloring->w3) { 661 ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 662 ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 663 } 664 w3 = coloring->w3; 665 666 /* Compute all the local scale factors, including ghost points */ 667 ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 668 ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 669 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 670 if (ctype == IS_COLORING_GHOSTED){ 671 col_start = 0; col_end = N; 672 } else if (ctype == IS_COLORING_LOCAL){ 673 xx = xx - start; 674 vscale_array = vscale_array - start; 675 col_start = start; col_end = N + start; 676 } 677 for (col=col_start; col<col_end; col++){ 678 /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 679 if (coloring->htype[0] == 'w') { 680 dx = 1.0 + unorm; 681 } else { 682 dx = xx[col]; 683 } 684 if (dx == 0.0) dx = 1.0; 685 #if !defined(PETSC_USE_COMPLEX) 686 if (dx < umin && dx >= 0.0) dx = umin; 687 else if (dx < 0.0 && dx > -umin) dx = -umin; 688 #else 689 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 690 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 691 #endif 692 dx *= epsilon; 693 vscale_array[col] = 1.0/dx; 694 } 695 if (ctype == IS_COLORING_LOCAL) vscale_array = vscale_array + start; 696 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 697 if (ctype == IS_COLORING_LOCAL){ 698 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700 } 701 702 if (coloring->vscaleforrow) { 703 vscaleforrow = coloring->vscaleforrow; 704 } else { 705 SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 706 } 707 708 /* 709 Loop over each color 710 */ 711 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 712 for (k=0; k<coloring->ncolors; k++) { 713 coloring->currentcolor = k; 714 ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 715 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 716 if (ctype == IS_COLORING_LOCAL) w3_array = w3_array - start; 717 718 if (prid == rank) printf(" [%d] color %d \n -----------\n",rank,k); 719 /* 720 Loop over each column associated with color 721 adding the perturbation to the vector w3. 722 */ 723 for (l=0; l<coloring->ncolumns[k]; l++) { 724 col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 725 if (coloring->htype[0] == 'w') { 726 dx = 1.0 + unorm; 727 } else { 728 dx = xx[col]; 729 } 730 if (dx == 0.0) dx = 1.0; 731 #if !defined(PETSC_USE_COMPLEX) 732 if (dx < umin && dx >= 0.0) dx = umin; 733 else if (dx < 0.0 && dx > -umin) dx = -umin; 734 #else 735 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 736 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 737 #endif 738 dx *= epsilon; 739 if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 740 w3_array[col] += dx; 741 } 742 if (ctype == IS_COLORING_LOCAL) w3_array = w3_array + start; 743 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 744 745 /* 746 Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 747 w2 = F(x1 + dx) - F(x1) 748 */ 749 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 750 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 751 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 752 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 753 754 /* 755 Loop over rows of vector, putting results into Jacobian matrix 756 */ 757 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 758 for (l=0; l<coloring->nrows[k]; l++) { 759 row = coloring->rows[k][l]; /* local row index */ 760 col = coloring->columnsforrow[k][l]; /* global column index */ 761 y[row] *= vscale_array[vscaleforrow[k][l]]; 762 srow = row + start; 763 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 764 } 765 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 766 } /* endof for each color */ 767 if (ctype == IS_COLORING_LOCAL) xx = xx + start; 768 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 769 ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 770 771 coloring->currentcolor = -1; 772 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 773 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 774 /*ierr = MatView(J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 775 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 776 777 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 778 if (flg) { 779 ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 780 } 781 ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 782 783 PetscFunctionReturn(0); 784 } 785 786 #undef __FUNCT__ 787 #define __FUNCT__ "MatFDColoringApplyTS" 788 /*@ 789 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 790 has been created, computes the Jacobian for a function via finite differences. 791 792 Collective on Mat, MatFDColoring, and Vec 793 794 Input Parameters: 795 + mat - location to store Jacobian 796 . coloring - coloring context created with MatFDColoringCreate() 797 . x1 - location at which Jacobian is to be computed 798 - sctx - optional context required by function (actually a SNES context) 799 800 Options Database Keys: 801 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 802 803 Level: intermediate 804 805 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 806 807 .keywords: coloring, Jacobian, finite differences 808 @*/ 809 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 810 { 811 PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f; 812 PetscErrorCode ierr; 813 PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow; 814 PetscScalar dx,*y,*xx,*w3_array; 815 PetscScalar *vscale_array; 816 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 817 Vec w1=coloring->w1,w2=coloring->w2,w3; 818 void *fctx = coloring->fctx; 819 PetscTruth flg; 820 821 PetscFunctionBegin; 822 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 823 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 824 PetscValidHeaderSpecific(x1,VEC_COOKIE,4); 825 826 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 827 if (!coloring->w3) { 828 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 829 ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 830 } 831 w3 = coloring->w3; 832 833 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 834 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 835 if (flg) { 836 ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 837 } else { 838 PetscTruth assembled; 839 ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 840 if (assembled) { 841 ierr = MatZeroEntries(J);CHKERRQ(ierr); 842 } 843 } 844 845 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 846 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 847 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 848 ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 849 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 850 851 /* 852 Compute all the scale factors and share with other processors 853 */ 854 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 855 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 856 for (k=0; k<coloring->ncolors; k++) { 857 /* 858 Loop over each column associated with color adding the 859 perturbation to the vector w3. 860 */ 861 for (l=0; l<coloring->ncolumns[k]; l++) { 862 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 863 dx = xx[col]; 864 if (dx == 0.0) dx = 1.0; 865 #if !defined(PETSC_USE_COMPLEX) 866 if (dx < umin && dx >= 0.0) dx = umin; 867 else if (dx < 0.0 && dx > -umin) dx = -umin; 868 #else 869 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 870 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 871 #endif 872 dx *= epsilon; 873 vscale_array[col] = 1.0/dx; 874 } 875 } 876 vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 877 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 878 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 879 880 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 881 else vscaleforrow = coloring->columnsforrow; 882 883 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 884 /* 885 Loop over each color 886 */ 887 for (k=0; k<coloring->ncolors; k++) { 888 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 889 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 890 /* 891 Loop over each column associated with color adding the 892 perturbation to the vector w3. 893 */ 894 for (l=0; l<coloring->ncolumns[k]; l++) { 895 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 896 dx = xx[col]; 897 if (dx == 0.0) dx = 1.0; 898 #if !defined(PETSC_USE_COMPLEX) 899 if (dx < umin && dx >= 0.0) dx = umin; 900 else if (dx < 0.0 && dx > -umin) dx = -umin; 901 #else 902 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 903 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 904 #endif 905 dx *= epsilon; 906 w3_array[col] += dx; 907 } 908 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 909 910 /* 911 Evaluate function at x1 + dx (here dx is a vector of perturbations) 912 */ 913 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 914 ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 915 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 916 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 917 918 /* 919 Loop over rows of vector, putting results into Jacobian matrix 920 */ 921 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 922 for (l=0; l<coloring->nrows[k]; l++) { 923 row = coloring->rows[k][l]; 924 col = coloring->columnsforrow[k][l]; 925 y[row] *= vscale_array[vscaleforrow[k][l]]; 926 srow = row + start; 927 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 928 } 929 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 930 } 931 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 932 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 933 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 934 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 935 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 936 PetscFunctionReturn(0); 937 } 938 939 940 #undef __FUNCT__ 941 #define __FUNCT__ "MatFDColoringSetRecompute()" 942 /*@C 943 MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 944 is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 945 no additional Jacobian's will be computed (the same one will be used) until this is 946 called again. 947 948 Collective on MatFDColoring 949 950 Input Parameters: 951 . c - the coloring context 952 953 Level: intermediate 954 955 Notes: The MatFDColoringSetFrequency() is ignored once this is called 956 957 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 958 959 .keywords: Mat, finite differences, coloring 960 @*/ 961 PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c) 962 { 963 PetscFunctionBegin; 964 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1); 965 c->usersetsrecompute = PETSC_TRUE; 966 c->recompute = PETSC_TRUE; 967 PetscFunctionReturn(0); 968 } 969 970 971