1 2 #ifdef PETSC_RCS_HEADER 3 static char vcid[] = "$Id: fdmatrix.c,v 1.27 1997/11/03 04:45:03 bsmith Exp bsmith $"; 4 #endif 5 6 /* 7 This is where the abstract matrix operations are defined that are 8 used for finite difference computations of Jacobians using coloring. 9 */ 10 11 #include "petsc.h" 12 #include "src/mat/matimpl.h" /*I "mat.h" I*/ 13 #include "src/vec/vecimpl.h" 14 #include "pinclude/pviewer.h" 15 16 #undef __FUNC__ 17 #define __FUNC__ "MatFDColoringView_Draw" 18 static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer) 19 { 20 int ierr,i,j,pause; 21 PetscTruth isnull; 22 Draw draw; 23 double xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0; 24 DrawButton button; 25 26 PetscFunctionBegin; 27 ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 28 ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 29 ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr); 30 31 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 32 xr += w; yr += h; xl = -w; yl = -h; 33 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 34 35 /* loop over colors */ 36 for (i=0; i<fd->ncolors; i++ ) { 37 for ( j=0; j<fd->nrows[i]; j++ ) { 38 y = fd->M - fd->rows[i][j] - fd->rstart; 39 x = fd->columnsforrow[i][j]; 40 DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 41 } 42 } 43 ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr); 44 ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr); 45 if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 46 ierr = DrawCheckResizedWindow(draw); 47 ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 48 while (button != BUTTON_RIGHT) { 49 ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr); 50 if (button == BUTTON_LEFT) scale = .5; 51 else if (button == BUTTON_CENTER) scale = 2.; 52 xl = scale*(xl + w - xc) + xc - w*scale; 53 xr = scale*(xr - w - xc) + xc + w*scale; 54 yl = scale*(yl + h - yc) + yc - h*scale; 55 yr = scale*(yr - h - yc) + yc + h*scale; 56 w *= scale; h *= scale; 57 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 58 /* loop over colors */ 59 for (i=0; i<fd->ncolors; i++ ) { 60 for ( j=0; j<fd->nrows[i]; j++ ) { 61 y = fd->M - fd->rows[i][j] - fd->rstart; 62 x = fd->columnsforrow[i][j]; 63 DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 64 } 65 } 66 ierr = DrawCheckResizedWindow(draw); 67 ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 68 } 69 70 PetscFunctionReturn(0); 71 } 72 73 #undef __FUNC__ 74 #define __FUNC__ "MatFDColoringView" 75 /*@C 76 MatFDColoringView - Views a finite difference coloring context. 77 78 Input Parameters: 79 . c - the coloring context 80 . viewer - visualization context 81 82 Notes: 83 The available visualization contexts include 84 $ VIEWER_STDOUT_SELF - standard output (default) 85 $ VIEWER_STDOUT_WORLD - synchronized standard 86 $ output where only the first processor opens 87 $ the file. All other processors send their 88 $ data to the first processor to print. 89 $ VIEWER_DRAWX_WORLD - graphical display of nonzero structure 90 91 .seealso: MatFDColoringCreate() 92 93 .keywords: Mat, finite differences, coloring, view 94 @*/ 95 int MatFDColoringView(MatFDColoring c,Viewer viewer) 96 { 97 ViewerType vtype; 98 int i,j,format,ierr; 99 FILE *fd; 100 101 PetscFunctionBegin; 102 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 103 if (viewer) {PetscValidHeader(viewer);} 104 else {viewer = VIEWER_STDOUT_SELF;} 105 106 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 107 if (vtype == DRAW_VIEWER) { 108 ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr); 109 PetscFunctionReturn(0); 110 } 111 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 112 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 113 PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n"); 114 PetscFPrintf(c->comm,fd," Error tolerance=%g\n",c->error_rel); 115 PetscFPrintf(c->comm,fd," Umin=%g\n",c->umin); 116 PetscFPrintf(c->comm,fd," Number of colors=%d\n",c->ncolors); 117 118 ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 119 if (format != VIEWER_FORMAT_ASCII_INFO) { 120 for ( i=0; i<c->ncolors; i++ ) { 121 PetscFPrintf(c->comm,fd," Information for color %d\n",i); 122 PetscFPrintf(c->comm,fd," Number of columns %d\n",c->ncolumns[i]); 123 for ( j=0; j<c->ncolumns[i]; j++ ) { 124 PetscFPrintf(c->comm,fd," %d\n",c->columns[i][j]); 125 } 126 PetscFPrintf(c->comm,fd," Number of rows %d\n",c->nrows[i]); 127 for ( j=0; j<c->nrows[i]; j++ ) { 128 PetscFPrintf(c->comm,fd," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]); 129 } 130 } 131 } 132 } 133 PetscFunctionReturn(0); 134 } 135 136 #undef __FUNC__ 137 #define __FUNC__ "MatFDColoringSetParameters" 138 /*@ 139 MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 140 a Jacobian matrix using finite differences. 141 142 $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 143 $ h = error_rel*u[i] if u[i] > umin 144 $ = error_rel*umin else 145 $ 146 $ dx_{i} = (0, ... 1, .... 0) 147 148 Input Parameters: 149 . coloring - the coloring context 150 . error_rel - relative error 151 . umin - minimum allowable u-value 152 153 .keywords: Mat, finite differences, coloring, set, parameters 154 155 .seealso: MatFDColoringCreate() 156 @*/ 157 int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 158 { 159 PetscFunctionBegin; 160 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 161 162 if (error != PETSC_DEFAULT) matfd->error_rel = error; 163 if (umin != PETSC_DEFAULT) matfd->umin = umin; 164 PetscFunctionReturn(0); 165 } 166 167 #undef __FUNC__ 168 #define __FUNC__ "MatFDColoringSetFrequency" 169 /*@ 170 MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 171 matrices. 172 173 Input Parameters: 174 . coloring - the coloring context 175 . freq - frequency (default is 1) 176 177 Notes: 178 Using a modified Newton strategy, where the Jacobian remains fixed for several 179 iterations, can be cost effective in terms of overall nonlinear solution 180 efficiency. This parameter indicates that a new Jacobian will be computed every 181 <freq> nonlinear iterations. 182 183 Options Database Keys: 184 $ -mat_fd_coloring_freq <freq> 185 186 .keywords: Mat, finite differences, coloring, set, frequency 187 @*/ 188 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 189 { 190 PetscFunctionBegin; 191 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 192 193 matfd->freq = freq; 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNC__ 198 #define __FUNC__ "MatFDColoringGetFrequency" 199 /*@ 200 MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 201 matrices. 202 203 Input Parameters: 204 . coloring - the coloring context 205 206 Output Parameters: 207 . freq - frequency (default is 1) 208 209 Notes: 210 Using a modified Newton strategy, where the Jacobian remains fixed for several 211 iterations, can be cost effective in terms of overall nonlinear solution 212 efficiency. This parameter indicates that a new Jacobian will be computed every 213 <freq> nonlinear iterations. 214 215 Options Database Keys: 216 $ -mat_fd_coloring_freq <freq> 217 218 .keywords: Mat, finite differences, coloring, get, frequency 219 @*/ 220 int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 221 { 222 PetscFunctionBegin; 223 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 224 225 *freq = matfd->freq; 226 PetscFunctionReturn(0); 227 } 228 229 #undef __FUNC__ 230 #define __FUNC__ "MatFDColoringSetFunction" 231 /*@C 232 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 233 234 Input Parameters: 235 . coloring - the coloring context 236 . f - the function 237 . fctx - the function context 238 239 .keywords: Mat, Jacobian, finite differences, set, function 240 @*/ 241 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 242 { 243 PetscFunctionBegin; 244 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 245 246 matfd->f = f; 247 matfd->fctx = fctx; 248 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNC__ 253 #define __FUNC__ "MatFDColoringSetFromOptions" 254 /*@ 255 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 256 the options database. 257 258 The Jacobian is estimated with the differencing approximation 259 $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 260 $ h = error_rel*u[i] if u[i] > umin 261 $ = error_rel*umin else 262 $ 263 $ dx_{i} = (0, ... 1, .... 0) 264 265 Input Parameters: 266 . coloring - the coloring context 267 268 Options Database Keys: 269 $ -mat_fd_coloring_error <err>, where <err> is the square root 270 $ of relative error in the function 271 $ -mat_fd_coloring_umin <umin>, where umin is described above 272 $ -mat_fd_coloring_freq <freq> where <freq> is the frequency of 273 $ computing a new Jacobian 274 $ -mat_fd_coloring_view 275 $ -mat_fd_coloring_view_info 276 $ -mat_fd_coloring_view_draw 277 278 .keywords: Mat, finite differences, parameters 279 @*/ 280 int MatFDColoringSetFromOptions(MatFDColoring matfd) 281 { 282 int ierr,flag,freq = 1; 283 double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 284 285 PetscFunctionBegin; 286 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 287 288 ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 289 ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 290 ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 291 ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 292 ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 293 ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 294 if (flag) { 295 ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 296 } 297 PetscFunctionReturn(0); 298 } 299 300 #undef __FUNC__ 301 #define __FUNC__ "MatFDColoringPrintHelp" 302 /*@ 303 MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 304 using coloring. 305 306 Input Parameter: 307 . fdcoloring - the MatFDColoring context 308 309 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 310 @*/ 311 int MatFDColoringPrintHelp(MatFDColoring fd) 312 { 313 PetscFunctionBegin; 314 PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 315 316 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel); 317 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin); 318 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq); 319 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n"); 320 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n"); 321 (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n"); 322 PetscFunctionReturn(0); 323 } 324 325 int MatFDColoringView_Private(MatFDColoring fd) 326 { 327 int ierr,flg; 328 329 PetscFunctionBegin; 330 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 331 if (flg) { 332 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 333 } 334 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 335 if (flg) { 336 ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 337 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 338 ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339 } 340 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 341 if (flg) { 342 ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 343 ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 344 } 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNC__ 349 #define __FUNC__ "MatFDColoringCreate" 350 /*@C 351 MatFDColoringCreate - Creates a matrix coloring context for finite difference 352 computation of Jacobians. 353 354 Input Parameters: 355 . mat - the matrix containing the nonzero structure of the Jacobian 356 . iscoloring - the coloring of the matrix 357 358 Output Parameter: 359 . color - the new coloring context 360 361 Options Database Keys: 362 $ -mat_fd_coloring_view 363 $ -mat_fd_coloring_view_draw 364 $ -mat_fd_coloring_view_info 365 366 .seealso: MatFDColoringDestroy() 367 @*/ 368 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 369 { 370 MatFDColoring c; 371 MPI_Comm comm; 372 int ierr,M,N; 373 374 PetscFunctionBegin; 375 ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 376 if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 377 378 PetscObjectGetComm((PetscObject)mat,&comm); 379 PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 380 PLogObjectCreate(c); 381 382 if (mat->ops.fdcoloringcreate) { 383 ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 384 } else { 385 SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 386 } 387 388 c->error_rel = 1.e-8; 389 c->umin = 1.e-6; 390 c->freq = 1; 391 392 ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 393 394 *color = c; 395 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNC__ 400 #define __FUNC__ "MatFDColoringDestroy" 401 /*@C 402 MatFDColoringDestroy - Destroys a matrix coloring context that was created 403 via MatFDColoringCreate(). 404 405 Input Parameter: 406 . c - coloring context 407 408 .seealso: MatFDColoringCreate() 409 @*/ 410 int MatFDColoringDestroy(MatFDColoring c) 411 { 412 int i,ierr; 413 414 PetscFunctionBegin; 415 if (--c->refct > 0) PetscFunctionReturn(0); 416 417 418 for ( i=0; i<c->ncolors; i++ ) { 419 if (c->columns[i]) PetscFree(c->columns[i]); 420 if (c->rows[i]) PetscFree(c->rows[i]); 421 if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 422 } 423 PetscFree(c->ncolumns); 424 PetscFree(c->columns); 425 PetscFree(c->nrows); 426 PetscFree(c->rows); 427 PetscFree(c->columnsforrow); 428 PetscFree(c->scale); 429 if (c->w1) { 430 ierr = VecDestroy(c->w1); CHKERRQ(ierr); 431 ierr = VecDestroy(c->w2); CHKERRQ(ierr); 432 ierr = VecDestroy(c->w3); CHKERRQ(ierr); 433 } 434 PLogObjectDestroy(c); 435 PetscHeaderDestroy(c); 436 PetscFunctionReturn(0); 437 } 438 439 #include "snes.h" 440 441 #undef __FUNC__ 442 #define __FUNC__ "MatFDColoringApply" 443 /*@ 444 MatFDColoringApply - Given a matrix for which a MatFDColoring context 445 has been created, computes the Jacobian for a function via finite differences. 446 447 Input Parameters: 448 . mat - location to store Jacobian 449 . coloring - coloring context created with MatFDColoringCreate() 450 . x1 - location at which Jacobian is to be computed 451 . sctx - optional context required by function (actually a SNES context) 452 453 Options Database Keys: 454 $ -mat_fd_coloring_freq <freq> 455 456 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 457 458 .keywords: coloring, Jacobian, finite differences 459 @*/ 460 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 461 { 462 int k,fg,ierr,N,start,end,l,row,col,srow; 463 Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 464 double epsilon = coloring->error_rel, umin = coloring->umin; 465 MPI_Comm comm = coloring->comm; 466 Vec w1,w2,w3; 467 int (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f; 468 void *fctx = coloring->fctx; 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(J,MAT_COOKIE); 472 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 473 PetscValidHeaderSpecific(x1,VEC_COOKIE); 474 475 476 if (!coloring->w1) { 477 ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 478 PLogObjectParent(coloring,coloring->w1); 479 ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 480 PLogObjectParent(coloring,coloring->w2); 481 ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 482 PLogObjectParent(coloring,coloring->w3); 483 } 484 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 485 486 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 487 if (fg) { 488 PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 489 } else { 490 ierr = MatZeroEntries(J); CHKERRQ(ierr); 491 } 492 493 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 494 ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 495 ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 496 ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 497 498 PetscMemzero(wscale,N*sizeof(Scalar)); 499 /* 500 Loop over each color 501 */ 502 503 for (k=0; k<coloring->ncolors; k++) { 504 ierr = VecCopy(x1,w3); CHKERRQ(ierr); 505 /* 506 Loop over each column associated with color adding the 507 perturbation to the vector w3. 508 */ 509 for (l=0; l<coloring->ncolumns[k]; l++) { 510 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 511 dx = xx[col-start]; 512 if (dx == 0.0) dx = 1.0; 513 #if !defined(USE_PETSC_COMPLEX) 514 if (dx < umin && dx >= 0.0) dx = umin; 515 else if (dx < 0.0 && dx > -umin) dx = -umin; 516 #else 517 if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 518 else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 519 #endif 520 dx *= epsilon; 521 wscale[col] = 1.0/dx; 522 VecSetValues(w3,1,&col,&dx,ADD_VALUES); 523 } 524 VecRestoreArray(x1,&xx); 525 /* 526 Evaluate function at x1 + dx (here dx is a vector of perturbations) 527 */ 528 ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 529 ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 530 /* Communicate scale to all processors */ 531 #if !defined(USE_PETSC_COMPLEX) 532 ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 533 #else 534 ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 535 #endif 536 /* 537 Loop over rows of vector, putting results into Jacobian matrix 538 */ 539 VecGetArray(w2,&y); 540 for (l=0; l<coloring->nrows[k]; l++) { 541 row = coloring->rows[k][l]; 542 col = coloring->columnsforrow[k][l]; 543 y[row] *= scale[col]; 544 srow = row + start; 545 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 546 } 547 VecRestoreArray(w2,&y); 548 } 549 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 550 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 551 PetscFunctionReturn(0); 552 } 553 554 #include "ts.h" 555 556 #undef __FUNC__ 557 #define __FUNC__ "MatFDColoringApplyTS" 558 /*@ 559 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 560 has been created, computes the Jacobian for a function via finite differences. 561 562 Input Parameters: 563 . mat - location to store Jacobian 564 . coloring - coloring context created with MatFDColoringCreate() 565 . x1 - location at which Jacobian is to be computed 566 . sctx - optional context required by function (actually a SNES context) 567 568 Options Database Keys: 569 $ -mat_fd_coloring_freq <freq> 570 571 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 572 573 .keywords: coloring, Jacobian, finite differences 574 @*/ 575 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx) 576 { 577 int k,fg,ierr,N,start,end,l,row,col,srow; 578 Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 579 double epsilon = coloring->error_rel, umin = coloring->umin; 580 MPI_Comm comm = coloring->comm; 581 Vec w1,w2,w3; 582 int (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f; 583 void *fctx = coloring->fctx; 584 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(J,MAT_COOKIE); 587 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 588 PetscValidHeaderSpecific(x1,VEC_COOKIE); 589 590 591 if (!coloring->w1) { 592 ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 593 PLogObjectParent(coloring,coloring->w1); 594 ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 595 PLogObjectParent(coloring,coloring->w2); 596 ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 597 PLogObjectParent(coloring,coloring->w3); 598 } 599 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 600 601 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 602 if (fg) { 603 PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n"); 604 } else { 605 ierr = MatZeroEntries(J); CHKERRQ(ierr); 606 } 607 608 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 609 ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 610 ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 611 ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr); 612 613 PetscMemzero(wscale,N*sizeof(Scalar)); 614 /* 615 Loop over each color 616 */ 617 618 for (k=0; k<coloring->ncolors; k++) { 619 ierr = VecCopy(x1,w3); CHKERRQ(ierr); 620 /* 621 Loop over each column associated with color adding the 622 perturbation to the vector w3. 623 */ 624 for (l=0; l<coloring->ncolumns[k]; l++) { 625 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 626 dx = xx[col-start]; 627 if (dx == 0.0) dx = 1.0; 628 #if !defined(USE_PETSC_COMPLEX) 629 if (dx < umin && dx >= 0.0) dx = umin; 630 else if (dx < 0.0 && dx > -umin) dx = -umin; 631 #else 632 if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 633 else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 634 #endif 635 dx *= epsilon; 636 wscale[col] = 1.0/dx; 637 VecSetValues(w3,1,&col,&dx,ADD_VALUES); 638 } 639 VecRestoreArray(x1,&xx); 640 /* 641 Evaluate function at x1 + dx (here dx is a vector of perturbations) 642 */ 643 ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr); 644 ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 645 /* Communicate scale to all processors */ 646 #if !defined(USE_PETSC_COMPLEX) 647 ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 648 #else 649 ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 650 #endif 651 /* 652 Loop over rows of vector, putting results into Jacobian matrix 653 */ 654 VecGetArray(w2,&y); 655 for (l=0; l<coloring->nrows[k]; l++) { 656 row = coloring->rows[k][l]; 657 col = coloring->columnsforrow[k][l]; 658 y[row] *= scale[col]; 659 srow = row + start; 660 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 661 } 662 VecRestoreArray(w2,&y); 663 } 664 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 665 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 666 PetscFunctionReturn(0); 667 } 668