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