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