1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: fdmatrix.c,v 1.17 1997/09/27 23:29:32 curfman 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 = DrawSyncClear(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 = DrawSyncFlush(draw); CHKERRQ(ierr); 42 ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr); 43 if (pause >= 0) { PetscSleep(pause); return 0;} 44 ierr = DrawCheckResizedWindow(draw); 45 ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0); 46 while (button != BUTTON_RIGHT) { 47 ierr = DrawSyncClear(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 = DrawSyncGetMouseButton(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__ "MatFDColoringSetFunction" 194 /*@ 195 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 196 197 Input Parameters: 198 . coloring - the coloring context 199 . f - the function 200 . fctx - the function context 201 202 .keywords: Mat, Jacobian, finite differences, set, function 203 @*/ 204 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx) 205 { 206 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 207 208 matfd->f = f; 209 matfd->fctx = fctx; 210 211 return 0; 212 } 213 214 #undef __FUNC__ 215 #define __FUNC__ "MatFDColoringSetFromOptions" 216 /*@ 217 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 218 the options database. 219 220 The Jacobian is estimated with the differencing approximation 221 $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 222 $ h = error_rel*u[i] if u[i] > umin 223 $ = error_rel*umin else 224 $ 225 $ dx_{i} = (0, ... 1, .... 0) 226 227 Input Parameters: 228 . coloring - the coloring context 229 230 Options Database Keys: 231 $ -mat_fd_coloring_error <err>, where <err> is the square root 232 $ of relative error in the function 233 $ -mat_fd_coloring_umin <umin>, where umin is described above 234 $ -mat_fd_coloring_freq <freq> where <freq> is the frequency of 235 $ computing a new Jacobian 236 237 .keywords: Mat, finite differences, parameters 238 @*/ 239 int MatFDColoringSetFromOptions(MatFDColoring matfd) 240 { 241 int ierr,flag,freq = 1; 242 double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 243 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 244 245 ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 246 ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 247 ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 248 ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 249 ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 250 ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 251 if (flag) { 252 ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 253 } 254 return 0; 255 } 256 257 #undef __FUNC__ 258 #define __FUNC__ "MatFDColoringPrintHelp" 259 /*@ 260 MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 261 using coloring. 262 263 Input Parameter: 264 . fdcoloring - the MatFDColoring context 265 266 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 267 @*/ 268 int MatFDColoringPrintHelp(MatFDColoring fd) 269 { 270 PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 271 272 PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel); 273 PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin); 274 PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq); 275 PetscPrintf(fd->comm,"-mat_fd_coloring_view\n"); 276 PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n"); 277 PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n"); 278 return 0; 279 } 280 281 int MatFDColoringView_Private(MatFDColoring fd) 282 { 283 int ierr,flg; 284 285 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 286 if (flg) { 287 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 288 } 289 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 290 if (flg) { 291 ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 292 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 293 ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 294 } 295 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 296 if (flg) { 297 ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 298 ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 299 } 300 return 0; 301 } 302 303 #undef __FUNC__ 304 #define __FUNC__ "MatFDColoringCreate" 305 /*@C 306 MatFDColoringCreate - Creates a matrix coloring context for finite difference 307 computation of Jacobians. 308 309 Input Parameters: 310 . mat - the matrix containing the nonzero structure of the Jacobian 311 . iscoloring - the coloring of the matrix 312 313 Output Parameter: 314 . color - the new coloring context 315 316 Options Database Keys: 317 $ -mat_fd_coloring_view 318 $ -mat_fd_coloring_view_draw 319 $ -mat_fd_coloring_view_info 320 321 .seealso: MatFDColoringDestroy() 322 @*/ 323 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 324 { 325 MatFDColoring c; 326 MPI_Comm comm; 327 int ierr,M,N; 328 329 ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 330 if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 331 332 PetscObjectGetComm((PetscObject)mat,&comm); 333 PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 334 PLogObjectCreate(c); 335 336 if (mat->ops.fdcoloringcreate) { 337 ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 338 } else { 339 SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 340 } 341 342 c->error_rel = 1.e-8; 343 c->umin = 1.e-6; 344 c->freq = 1; 345 346 ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 347 348 *color = c; 349 350 return 0; 351 } 352 353 #undef __FUNC__ 354 #define __FUNC__ "MatFDColoringDestroy" 355 /*@C 356 MatFDColoringDestroy - Destroys a matrix coloring context that was created 357 via MatFDColoringCreate(). 358 359 Input Parameter: 360 . c - coloring context 361 362 .seealso: MatFDColoringCreate() 363 @*/ 364 int MatFDColoringDestroy(MatFDColoring c) 365 { 366 int i,ierr,flag; 367 368 if (--c->refct > 0) return 0; 369 370 ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag); 371 if (flag) { 372 ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr); 373 } 374 ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag); 375 if (flag) { 376 ierr = ViewerPushFormat(VIEWER_STDOUT_(c->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 377 ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr); 378 ierr = ViewerPopFormat(VIEWER_STDOUT_(c->comm)); 379 } 380 381 for ( i=0; i<c->ncolors; i++ ) { 382 if (c->columns[i]) PetscFree(c->columns[i]); 383 if (c->rows[i]) PetscFree(c->rows[i]); 384 if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 385 } 386 PetscFree(c->ncolumns); 387 PetscFree(c->columns); 388 PetscFree(c->nrows); 389 PetscFree(c->rows); 390 PetscFree(c->columnsforrow); 391 PetscFree(c->scale); 392 if (c->w1) { 393 ierr = VecDestroy(c->w1); CHKERRQ(ierr); 394 ierr = VecDestroy(c->w2); CHKERRQ(ierr); 395 ierr = VecDestroy(c->w3); CHKERRQ(ierr); 396 } 397 PLogObjectDestroy(c); 398 PetscHeaderDestroy(c); 399 return 0; 400 } 401 402 #include "snes.h" 403 404 #undef __FUNC__ 405 #define __FUNC__ "MatFDColoringApply" 406 /*@ 407 MatFDColoringApply - Given a matrix for which a MatFDColoring context 408 has been created, computes the Jacobian for a function via finite differences. 409 410 Input Parameters: 411 . mat - location to store Jacobian 412 . coloring - coloring context created with MatFDColoringCreate() 413 . x1 - location at which Jacobian is to be computed 414 . sctx - optional context required by function (actually a SNES context) 415 416 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 417 418 .keywords: coloring, Jacobian, finite differences 419 @*/ 420 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 421 { 422 int k,fg,ierr,N,start,end,l,row,col,srow; 423 Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 424 double epsilon = coloring->error_rel, umin = coloring->umin; 425 MPI_Comm comm = coloring->comm; 426 Vec w1,w2,w3; 427 int (*f)(void *,Vec,Vec,void *) = coloring->f; 428 void *fctx = coloring->fctx; 429 430 PetscValidHeaderSpecific(J,MAT_COOKIE); 431 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 432 PetscValidHeaderSpecific(x1,VEC_COOKIE); 433 434 /* 435 ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr); 436 if ((freq > 1) && ((it % freq) != 1)) { 437 PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq); 438 *flag = SAME_PRECONDITIONER; 439 return 0; 440 } else { 441 PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq); 442 *flag = SAME_NONZERO_PATTERN; 443 }*/ 444 445 if (!coloring->w1) { 446 ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 447 PLogObjectParent(coloring,coloring->w1); 448 ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 449 PLogObjectParent(coloring,coloring->w2); 450 ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 451 PLogObjectParent(coloring,coloring->w3); 452 } 453 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 454 455 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 456 if (fg) { 457 PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 458 } else { 459 ierr = MatZeroEntries(J); CHKERRQ(ierr); 460 } 461 462 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 463 ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 464 ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 465 ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 466 467 PetscMemzero(wscale,N*sizeof(Scalar)); 468 /* 469 Loop over each color 470 */ 471 472 for (k=0; k<coloring->ncolors; k++) { 473 ierr = VecCopy(x1,w3); CHKERRQ(ierr); 474 /* 475 Loop over each column associated with color adding the 476 perturbation to the vector w3. 477 */ 478 for (l=0; l<coloring->ncolumns[k]; l++) { 479 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 480 dx = xx[col-start]; 481 if (dx == 0.0) dx = 1.0; 482 #if !defined(PETSC_COMPLEX) 483 if (dx < umin && dx >= 0.0) dx = umin; 484 else if (dx < 0.0 && dx > -umin) dx = -umin; 485 #else 486 if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 487 else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 488 #endif 489 dx *= epsilon; 490 wscale[col] = 1.0/dx; 491 VecSetValues(w3,1,&col,&dx,ADD_VALUES); 492 } 493 VecRestoreArray(x1,&xx); 494 /* 495 Evaluate function at x1 + dx (here dx is a vector of perturbations) 496 */ 497 ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 498 ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 499 /* Communicate scale to all processors */ 500 #if !defined(PETSC_COMPLEX) 501 MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 502 #else 503 MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 504 #endif 505 /* 506 Loop over rows of vector, putting results into Jacobian matrix 507 */ 508 VecGetArray(w2,&y); 509 for (l=0; l<coloring->nrows[k]; l++) { 510 row = coloring->rows[k][l]; 511 col = coloring->columnsforrow[k][l]; 512 y[row] *= scale[col]; 513 srow = row + start; 514 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 515 } 516 VecRestoreArray(w2,&y); 517 } 518 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 519 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 520 return 0; 521 } 522