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