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