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