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