1 /*$Id: fdmatrix.c,v 1.92 2001/08/21 21:03:06 bsmith Exp $*/ 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 int MAT_FDCOLORING_COOKIE; 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatFDColoringSetF" 15 int 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 int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 25 { 26 MatFDColoring fd = (MatFDColoring)Aa; 27 int ierr,i,j; 28 PetscReal x,y; 29 30 PetscFunctionBegin; 31 32 /* loop over colors */ 33 for (i=0; i<fd->ncolors; i++) { 34 for (j=0; j<fd->nrows[i]; j++) { 35 y = fd->M - fd->rows[i][j] - fd->rstart; 36 x = fd->columnsforrow[i][j]; 37 ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 38 } 39 } 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "MatFDColoringView_Draw" 45 static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 46 { 47 int ierr; 48 PetscTruth isnull; 49 PetscDraw draw; 50 PetscReal xr,yr,xl,yl,h,w; 51 52 PetscFunctionBegin; 53 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 54 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 55 56 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 57 58 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 59 xr += w; yr += h; xl = -w; yl = -h; 60 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 61 ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 62 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "MatFDColoringView" 68 /*@C 69 MatFDColoringView - Views a finite difference coloring context. 70 71 Collective on MatFDColoring 72 73 Input Parameters: 74 + c - the coloring context 75 - viewer - visualization context 76 77 Level: intermediate 78 79 Notes: 80 The available visualization contexts include 81 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 82 . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 83 output where only the first processor opens 84 the file. All other processors send their 85 data to the first processor to print. 86 - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 87 88 Notes: 89 Since PETSc uses only a small number of basic colors (currently 33), if the coloring 90 involves more than 33 then some seemingly identical colors are displayed making it look 91 like an illegal coloring. This is just a graphical artifact. 92 93 .seealso: MatFDColoringCreate() 94 95 .keywords: Mat, finite differences, coloring, view 96 @*/ 97 int MatFDColoringView(MatFDColoring c,PetscViewer viewer) 98 { 99 int i,j,ierr; 100 PetscTruth isdraw,isascii; 101 PetscViewerFormat format; 102 103 PetscFunctionBegin; 104 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 105 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 106 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 107 PetscCheckSameComm(c,viewer); 108 109 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 110 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 111 if (isdraw) { 112 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 113 } else if (isascii) { 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 int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 168 { 169 PetscFunctionBegin; 170 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 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 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 205 { 206 PetscFunctionBegin; 207 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 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 int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 243 { 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 246 247 *freq = matfd->freq; 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "MatFDColoringSetFunction" 253 /*@C 254 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 255 256 Collective on MatFDColoring 257 258 Input Parameters: 259 + coloring - the coloring context 260 . f - the function 261 - fctx - the optional user-defined function context 262 263 Level: intermediate 264 265 Notes: 266 In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 267 be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 268 with the TS solvers. 269 270 .keywords: Mat, Jacobian, finite differences, set, function 271 @*/ 272 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 273 { 274 PetscFunctionBegin; 275 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 276 277 matfd->f = f; 278 matfd->fctx = fctx; 279 280 PetscFunctionReturn(0); 281 } 282 283 #undef __FUNCT__ 284 #define __FUNCT__ "MatFDColoringSetFromOptions" 285 /*@ 286 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 287 the options database. 288 289 Collective on MatFDColoring 290 291 The Jacobian, F'(u), is estimated with the differencing approximation 292 .vb 293 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 294 h = error_rel*u[i] if abs(u[i]) > umin 295 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 296 dx_{i} = (0, ... 1, .... 0) 297 .ve 298 299 Input Parameter: 300 . coloring - the coloring context 301 302 Options Database Keys: 303 + -mat_fd_coloring_err <err> - Sets <err> (square root 304 of relative error in the function) 305 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 306 . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 307 . -mat_fd_coloring_view - Activates basic viewing 308 . -mat_fd_coloring_view_info - Activates viewing info 309 - -mat_fd_coloring_view_draw - Activates drawing 310 311 Level: intermediate 312 313 .keywords: Mat, finite differences, parameters 314 315 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 316 317 @*/ 318 int MatFDColoringSetFromOptions(MatFDColoring matfd) 319 { 320 int ierr; 321 322 PetscFunctionBegin; 323 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 324 325 ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 326 ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 327 ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 328 ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 329 /* not used here; but so they are presented in the GUI */ 330 ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 331 ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 332 ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 333 PetscOptionsEnd();CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "MatFDColoringView_Private" 339 int MatFDColoringView_Private(MatFDColoring fd) 340 { 341 int ierr; 342 PetscTruth flg; 343 344 PetscFunctionBegin; 345 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 346 if (flg) { 347 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 348 } 349 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 350 if (flg) { 351 ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 352 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 353 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 354 } 355 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 356 if (flg) { 357 ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 358 ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 359 } 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "MatFDColoringCreate" 365 /*@C 366 MatFDColoringCreate - Creates a matrix coloring context for finite difference 367 computation of Jacobians. 368 369 Collective on Mat 370 371 Input Parameters: 372 + mat - the matrix containing the nonzero structure of the Jacobian 373 - iscoloring - the coloring of the matrix 374 375 Output Parameter: 376 . color - the new coloring context 377 378 Level: intermediate 379 380 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 381 MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 382 MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(), 383 MatFDColoringSetParameters() 384 @*/ 385 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 386 { 387 MatFDColoring c; 388 MPI_Comm comm; 389 int ierr,M,N; 390 391 PetscFunctionBegin; 392 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 393 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 394 if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 395 396 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 397 PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 398 PetscLogObjectCreate(c); 399 400 if (mat->ops->fdcoloringcreate) { 401 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 402 } else { 403 SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 404 } 405 406 c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 407 c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 408 c->freq = 1; 409 c->usersetsrecompute = PETSC_FALSE; 410 c->recompute = PETSC_FALSE; 411 c->currentcolor = -1; 412 413 *color = c; 414 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "MatFDColoringDestroy" 420 /*@C 421 MatFDColoringDestroy - Destroys a matrix coloring context that was created 422 via MatFDColoringCreate(). 423 424 Collective on MatFDColoring 425 426 Input Parameter: 427 . c - coloring context 428 429 Level: intermediate 430 431 .seealso: MatFDColoringCreate() 432 @*/ 433 int MatFDColoringDestroy(MatFDColoring c) 434 { 435 int i,ierr; 436 437 PetscFunctionBegin; 438 if (--c->refct > 0) PetscFunctionReturn(0); 439 440 for (i=0; i<c->ncolors; i++) { 441 if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 442 if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 443 if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 444 if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 445 } 446 ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 447 ierr = PetscFree(c->columns);CHKERRQ(ierr); 448 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 449 ierr = PetscFree(c->rows);CHKERRQ(ierr); 450 ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 451 if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 452 if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 453 if (c->w1) { 454 ierr = VecDestroy(c->w1);CHKERRQ(ierr); 455 ierr = VecDestroy(c->w2);CHKERRQ(ierr); 456 ierr = VecDestroy(c->w3);CHKERRQ(ierr); 457 } 458 PetscLogObjectDestroy(c); 459 PetscHeaderDestroy(c); 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 465 /*@C 466 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 467 that are currently being perturbed. 468 469 Not Collective 470 471 Input Parameters: 472 . coloring - coloring context created with MatFDColoringCreate() 473 474 Output Parameters: 475 + n - the number of local columns being perturbed 476 - cols - the column indices, in global numbering 477 478 Level: intermediate 479 480 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 481 482 .keywords: coloring, Jacobian, finite differences 483 @*/ 484 EXTERN int MatFDColoringGetPerturbedColumns(MatFDColoring coloring,int *n,int **cols) 485 { 486 PetscFunctionBegin; 487 if (coloring->currentcolor >= 0) { 488 *n = coloring->ncolumns[coloring->currentcolor]; 489 *cols = coloring->columns[coloring->currentcolor]; 490 } else { 491 *n = 0; 492 } 493 PetscFunctionReturn(0); 494 } 495 496 #undef __FUNCT__ 497 #define __FUNCT__ "MatFDColoringApply" 498 /*@ 499 MatFDColoringApply - Given a matrix for which a MatFDColoring context 500 has been created, computes the Jacobian for a function via finite differences. 501 502 Collective on MatFDColoring 503 504 Input Parameters: 505 + mat - location to store Jacobian 506 . coloring - coloring context created with MatFDColoringCreate() 507 . x1 - location at which Jacobian is to be computed 508 - sctx - optional context required by function (actually a SNES context) 509 510 Options Database Keys: 511 + -mat_fd_coloring_freq <freq> - Sets coloring frequency 512 . -mat_fd_coloring_view - Activates basic viewing or coloring 513 . -mat_fd_coloring_view_draw - Activates drawing of coloring 514 - -mat_fd_coloring_view_info - Activates viewing of coloring info 515 516 Level: intermediate 517 518 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 519 520 .keywords: coloring, Jacobian, finite differences 521 @*/ 522 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 523 { 524 int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f; 525 int k,ierr,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); 536 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 537 PetscValidHeaderSpecific(x1,VEC_COOKIE); 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 = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 592 } 593 594 /* 595 Compute all the scale factors and share with other processors 596 */ 597 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 598 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 599 for (k=0; k<coloring->ncolors; k++) { 600 /* 601 Loop over each column associated with color adding the 602 perturbation to the vector w3. 603 */ 604 for (l=0; l<coloring->ncolumns[k]; l++) { 605 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 606 dx = xx[col]; 607 if (dx == 0.0) dx = 1.0; 608 #if !defined(PETSC_USE_COMPLEX) 609 if (dx < umin && dx >= 0.0) dx = umin; 610 else if (dx < 0.0 && dx > -umin) dx = -umin; 611 #else 612 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 613 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 614 #endif 615 dx *= epsilon; 616 vscale_array[col] = 1.0/dx; 617 } 618 } 619 vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 620 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 621 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 622 623 /* ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); 624 ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/ 625 626 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 627 else vscaleforrow = coloring->columnsforrow; 628 629 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 630 /* 631 Loop over each color 632 */ 633 for (k=0; k<coloring->ncolors; k++) { 634 coloring->currentcolor = k; 635 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 636 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 637 /* 638 Loop over each column associated with color adding the 639 perturbation to the vector w3. 640 */ 641 for (l=0; l<coloring->ncolumns[k]; l++) { 642 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 643 dx = xx[col]; 644 if (dx == 0.0) dx = 1.0; 645 #if !defined(PETSC_USE_COMPLEX) 646 if (dx < umin && dx >= 0.0) dx = umin; 647 else if (dx < 0.0 && dx > -umin) dx = -umin; 648 #else 649 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 650 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 651 #endif 652 dx *= epsilon; 653 if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter"); 654 w3_array[col] += dx; 655 } 656 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 657 658 /* 659 Evaluate function at x1 + dx (here dx is a vector of perturbations) 660 */ 661 662 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 663 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 664 665 /* 666 Loop over rows of vector, putting results into Jacobian matrix 667 */ 668 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 669 for (l=0; l<coloring->nrows[k]; l++) { 670 row = coloring->rows[k][l]; 671 col = coloring->columnsforrow[k][l]; 672 y[row] *= vscale_array[vscaleforrow[k][l]]; 673 srow = row + start; 674 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 675 } 676 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 677 } 678 coloring->currentcolor = -1; 679 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 680 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 681 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 682 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 683 } 684 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 685 686 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 687 if (flg) { 688 ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 689 } 690 ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 691 692 PetscFunctionReturn(0); 693 } 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "MatFDColoringApplyTS" 697 /*@ 698 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 699 has been created, computes the Jacobian for a function via finite differences. 700 701 Collective on Mat, MatFDColoring, and Vec 702 703 Input Parameters: 704 + mat - location to store Jacobian 705 . coloring - coloring context created with MatFDColoringCreate() 706 . x1 - location at which Jacobian is to be computed 707 - sctx - optional context required by function (actually a SNES context) 708 709 Options Database Keys: 710 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 711 712 Level: intermediate 713 714 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 715 716 .keywords: coloring, Jacobian, finite differences 717 @*/ 718 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 719 { 720 int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 721 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 722 PetscScalar dx,mone = -1.0,*y,*xx,*w3_array; 723 PetscScalar *vscale_array; 724 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 725 Vec w1,w2,w3; 726 void *fctx = coloring->fctx; 727 PetscTruth flg; 728 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(J,MAT_COOKIE); 731 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 732 PetscValidHeaderSpecific(x1,VEC_COOKIE); 733 734 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 735 if (!coloring->w1) { 736 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 737 PetscLogObjectParent(coloring,coloring->w1); 738 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 739 PetscLogObjectParent(coloring,coloring->w2); 740 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 741 PetscLogObjectParent(coloring,coloring->w3); 742 } 743 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 744 745 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 746 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 747 if (flg) { 748 PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 749 } else { 750 ierr = MatZeroEntries(J);CHKERRQ(ierr); 751 } 752 753 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 754 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 755 ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 756 757 /* 758 Compute all the scale factors and share with other processors 759 */ 760 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 761 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 762 for (k=0; k<coloring->ncolors; k++) { 763 /* 764 Loop over each column associated with color adding the 765 perturbation to the vector w3. 766 */ 767 for (l=0; l<coloring->ncolumns[k]; l++) { 768 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 769 dx = xx[col]; 770 if (dx == 0.0) dx = 1.0; 771 #if !defined(PETSC_USE_COMPLEX) 772 if (dx < umin && dx >= 0.0) dx = umin; 773 else if (dx < 0.0 && dx > -umin) dx = -umin; 774 #else 775 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 776 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 777 #endif 778 dx *= epsilon; 779 vscale_array[col] = 1.0/dx; 780 } 781 } 782 vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 783 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 784 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 785 786 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 787 else vscaleforrow = coloring->columnsforrow; 788 789 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 790 /* 791 Loop over each color 792 */ 793 for (k=0; k<coloring->ncolors; k++) { 794 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 795 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 796 /* 797 Loop over each column associated with color adding the 798 perturbation to the vector w3. 799 */ 800 for (l=0; l<coloring->ncolumns[k]; l++) { 801 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 802 dx = xx[col]; 803 if (dx == 0.0) dx = 1.0; 804 #if !defined(PETSC_USE_COMPLEX) 805 if (dx < umin && dx >= 0.0) dx = umin; 806 else if (dx < 0.0 && dx > -umin) dx = -umin; 807 #else 808 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 809 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 810 #endif 811 dx *= epsilon; 812 w3_array[col] += dx; 813 } 814 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 815 816 /* 817 Evaluate function at x1 + dx (here dx is a vector of perturbations) 818 */ 819 ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 820 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 821 822 /* 823 Loop over rows of vector, putting results into Jacobian matrix 824 */ 825 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 826 for (l=0; l<coloring->nrows[k]; l++) { 827 row = coloring->rows[k][l]; 828 col = coloring->columnsforrow[k][l]; 829 y[row] *= vscale_array[vscaleforrow[k][l]]; 830 srow = row + start; 831 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 832 } 833 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 834 } 835 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 836 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 837 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 838 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 839 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 843 844 #undef __FUNCT__ 845 #define __FUNCT__ "MatFDColoringSetRecompute()" 846 /*@C 847 MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 848 is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 849 no additional Jacobian's will be computed (the same one will be used) until this is 850 called again. 851 852 Collective on MatFDColoring 853 854 Input Parameters: 855 . c - the coloring context 856 857 Level: intermediate 858 859 Notes: The MatFDColoringSetFrequency() is ignored once this is called 860 861 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 862 863 .keywords: Mat, finite differences, coloring 864 @*/ 865 int MatFDColoringSetRecompute(MatFDColoring c) 866 { 867 PetscFunctionBegin; 868 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 869 c->usersetsrecompute = PETSC_TRUE; 870 c->recompute = PETSC_TRUE; 871 PetscFunctionReturn(0); 872 } 873 874 875