1 /*$Id: fdmatrix.c,v 1.90 2001/08/06 21:16:12 bsmith Exp balay $*/ 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 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 314 315 @*/ 316 int MatFDColoringSetFromOptions(MatFDColoring matfd) 317 { 318 int ierr; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 322 323 ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 324 ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 325 ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 326 ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 327 /* not used here; but so they are presented in the GUI */ 328 ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 329 ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 330 ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 331 PetscOptionsEnd();CHKERRQ(ierr); 332 PetscFunctionReturn(0); 333 } 334 335 #undef __FUNCT__ 336 #define __FUNCT__ "MatFDColoringView_Private" 337 int MatFDColoringView_Private(MatFDColoring fd) 338 { 339 int ierr; 340 PetscTruth flg; 341 342 PetscFunctionBegin; 343 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 344 if (flg) { 345 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 346 } 347 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 348 if (flg) { 349 ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 350 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 351 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 352 } 353 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 354 if (flg) { 355 ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 356 ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 357 } 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "MatFDColoringCreate" 363 /*@C 364 MatFDColoringCreate - Creates a matrix coloring context for finite difference 365 computation of Jacobians. 366 367 Collective on Mat 368 369 Input Parameters: 370 + mat - the matrix containing the nonzero structure of the Jacobian 371 - iscoloring - the coloring of the matrix 372 373 Output Parameter: 374 . color - the new coloring context 375 376 Options Database Keys: 377 + -mat_fd_coloring_view - Activates basic viewing or coloring 378 . -mat_fd_coloring_view_draw - Activates drawing of coloring 379 - -mat_fd_coloring_view_info - Activates viewing of coloring info 380 381 Level: intermediate 382 383 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 384 MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 385 MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(), 386 MatFDColoringSetParameters() 387 @*/ 388 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 389 { 390 MatFDColoring c; 391 MPI_Comm comm; 392 int ierr,M,N; 393 394 PetscFunctionBegin; 395 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 396 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 397 if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 398 399 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 400 PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 401 PetscLogObjectCreate(c); 402 403 if (mat->ops->fdcoloringcreate) { 404 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 405 } else { 406 SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 407 } 408 409 c->error_rel = 1.e-8; 410 c->umin = 1.e-6; 411 c->freq = 1; 412 c->usersetsrecompute = PETSC_FALSE; 413 c->recompute = PETSC_FALSE; 414 415 ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 416 417 *color = c; 418 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421 422 #undef __FUNCT__ 423 #define __FUNCT__ "MatFDColoringDestroy" 424 /*@C 425 MatFDColoringDestroy - Destroys a matrix coloring context that was created 426 via MatFDColoringCreate(). 427 428 Collective on MatFDColoring 429 430 Input Parameter: 431 . c - coloring context 432 433 Level: intermediate 434 435 .seealso: MatFDColoringCreate() 436 @*/ 437 int MatFDColoringDestroy(MatFDColoring c) 438 { 439 int i,ierr; 440 441 PetscFunctionBegin; 442 if (--c->refct > 0) PetscFunctionReturn(0); 443 444 for (i=0; i<c->ncolors; i++) { 445 if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 446 if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 447 if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 448 if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 449 } 450 ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 451 ierr = PetscFree(c->columns);CHKERRQ(ierr); 452 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 453 ierr = PetscFree(c->rows);CHKERRQ(ierr); 454 ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 455 if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 456 if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 457 if (c->w1) { 458 ierr = VecDestroy(c->w1);CHKERRQ(ierr); 459 ierr = VecDestroy(c->w2);CHKERRQ(ierr); 460 ierr = VecDestroy(c->w3);CHKERRQ(ierr); 461 } 462 PetscLogObjectDestroy(c); 463 PetscHeaderDestroy(c); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "MatFDColoringApply" 469 /*@ 470 MatFDColoringApply - Given a matrix for which a MatFDColoring context 471 has been created, computes the Jacobian for a function via finite differences. 472 473 Collective on MatFDColoring 474 475 Input Parameters: 476 + mat - location to store Jacobian 477 . coloring - coloring context created with MatFDColoringCreate() 478 . x1 - location at which Jacobian is to be computed 479 - sctx - optional context required by function (actually a SNES context) 480 481 Options Database Keys: 482 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 483 484 Level: intermediate 485 486 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 487 488 .keywords: coloring, Jacobian, finite differences 489 @*/ 490 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 491 { 492 int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f; 493 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 494 PetscScalar dx,mone = -1.0,*y,*xx,*w3_array; 495 PetscScalar *vscale_array; 496 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 497 Vec w1,w2,w3; 498 void *fctx = coloring->fctx; 499 PetscTruth flg; 500 501 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(J,MAT_COOKIE); 504 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 505 PetscValidHeaderSpecific(x1,VEC_COOKIE); 506 507 if (coloring->usersetsrecompute) { 508 if (!coloring->recompute) { 509 *flag = SAME_PRECONDITIONER; 510 PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"); 511 PetscFunctionReturn(0); 512 } else { 513 coloring->recompute = PETSC_FALSE; 514 } 515 } 516 517 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 518 if (J->ops->fdcoloringapply) { 519 ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 520 } else { 521 522 if (!coloring->w1) { 523 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 524 PetscLogObjectParent(coloring,coloring->w1); 525 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 526 PetscLogObjectParent(coloring,coloring->w2); 527 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 528 PetscLogObjectParent(coloring,coloring->w3); 529 } 530 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 531 532 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 533 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 534 if (flg) { 535 PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 536 } else { 537 ierr = MatZeroEntries(J);CHKERRQ(ierr); 538 } 539 540 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 541 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 542 543 /* 544 This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets 545 coloring->F for the coarser grids from the finest 546 */ 547 if (coloring->F) { 548 ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 549 ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 550 if (m1 != m2) { 551 coloring->F = 0; 552 } 553 } 554 555 if (coloring->F) { 556 w1 = coloring->F; /* use already computed value of function */ 557 coloring->F = 0; 558 } else { 559 ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 560 } 561 562 /* 563 Compute all the scale factors and share with other processors 564 */ 565 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 566 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 567 for (k=0; k<coloring->ncolors; k++) { 568 /* 569 Loop over each column associated with color adding the 570 perturbation to the vector w3. 571 */ 572 for (l=0; l<coloring->ncolumns[k]; l++) { 573 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 574 dx = xx[col]; 575 if (dx == 0.0) dx = 1.0; 576 #if !defined(PETSC_USE_COMPLEX) 577 if (dx < umin && dx >= 0.0) dx = umin; 578 else if (dx < 0.0 && dx > -umin) dx = -umin; 579 #else 580 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 581 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 582 #endif 583 dx *= epsilon; 584 vscale_array[col] = 1.0/dx; 585 } 586 } 587 vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 588 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 589 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 590 591 /* ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); 592 ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/ 593 594 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 595 else vscaleforrow = coloring->columnsforrow; 596 597 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 598 /* 599 Loop over each color 600 */ 601 for (k=0; k<coloring->ncolors; k++) { 602 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 603 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 604 /* 605 Loop over each column associated with color adding the 606 perturbation to the vector w3. 607 */ 608 for (l=0; l<coloring->ncolumns[k]; l++) { 609 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 610 dx = xx[col]; 611 if (dx == 0.0) dx = 1.0; 612 #if !defined(PETSC_USE_COMPLEX) 613 if (dx < umin && dx >= 0.0) dx = umin; 614 else if (dx < 0.0 && dx > -umin) dx = -umin; 615 #else 616 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 617 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 618 #endif 619 dx *= epsilon; 620 if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter"); 621 w3_array[col] += dx; 622 } 623 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 624 625 /* 626 Evaluate function at x1 + dx (here dx is a vector of perturbations) 627 */ 628 629 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 630 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 631 632 /* 633 Loop over rows of vector, putting results into Jacobian matrix 634 */ 635 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 636 for (l=0; l<coloring->nrows[k]; l++) { 637 row = coloring->rows[k][l]; 638 col = coloring->columnsforrow[k][l]; 639 y[row] *= vscale_array[vscaleforrow[k][l]]; 640 srow = row + start; 641 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 642 } 643 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 644 } 645 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 646 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 647 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 648 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 649 } 650 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 651 652 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 653 if (flg) { 654 ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 655 } 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "MatFDColoringApplyTS" 661 /*@ 662 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 663 has been created, computes the Jacobian for a function via finite differences. 664 665 Collective on Mat, MatFDColoring, and Vec 666 667 Input Parameters: 668 + mat - location to store Jacobian 669 . coloring - coloring context created with MatFDColoringCreate() 670 . x1 - location at which Jacobian is to be computed 671 - sctx - optional context required by function (actually a SNES context) 672 673 Options Database Keys: 674 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 675 676 Level: intermediate 677 678 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 679 680 .keywords: coloring, Jacobian, finite differences 681 @*/ 682 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 683 { 684 int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 685 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 686 PetscScalar dx,mone = -1.0,*y,*xx,*w3_array; 687 PetscScalar *vscale_array; 688 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 689 Vec w1,w2,w3; 690 void *fctx = coloring->fctx; 691 PetscTruth flg; 692 693 PetscFunctionBegin; 694 PetscValidHeaderSpecific(J,MAT_COOKIE); 695 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 696 PetscValidHeaderSpecific(x1,VEC_COOKIE); 697 698 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 699 if (!coloring->w1) { 700 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 701 PetscLogObjectParent(coloring,coloring->w1); 702 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 703 PetscLogObjectParent(coloring,coloring->w2); 704 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 705 PetscLogObjectParent(coloring,coloring->w3); 706 } 707 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 708 709 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 710 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 711 if (flg) { 712 PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 713 } else { 714 ierr = MatZeroEntries(J);CHKERRQ(ierr); 715 } 716 717 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 718 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 719 ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 720 721 /* 722 Compute all the scale factors and share with other processors 723 */ 724 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 725 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 726 for (k=0; k<coloring->ncolors; k++) { 727 /* 728 Loop over each column associated with color adding the 729 perturbation to the vector w3. 730 */ 731 for (l=0; l<coloring->ncolumns[k]; l++) { 732 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 733 dx = xx[col]; 734 if (dx == 0.0) dx = 1.0; 735 #if !defined(PETSC_USE_COMPLEX) 736 if (dx < umin && dx >= 0.0) dx = umin; 737 else if (dx < 0.0 && dx > -umin) dx = -umin; 738 #else 739 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 740 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 741 #endif 742 dx *= epsilon; 743 vscale_array[col] = 1.0/dx; 744 } 745 } 746 vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 747 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 748 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 749 750 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 751 else vscaleforrow = coloring->columnsforrow; 752 753 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 754 /* 755 Loop over each color 756 */ 757 for (k=0; k<coloring->ncolors; k++) { 758 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 759 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 760 /* 761 Loop over each column associated with color adding the 762 perturbation to the vector w3. 763 */ 764 for (l=0; l<coloring->ncolumns[k]; l++) { 765 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 766 dx = xx[col]; 767 if (dx == 0.0) dx = 1.0; 768 #if !defined(PETSC_USE_COMPLEX) 769 if (dx < umin && dx >= 0.0) dx = umin; 770 else if (dx < 0.0 && dx > -umin) dx = -umin; 771 #else 772 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 773 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 774 #endif 775 dx *= epsilon; 776 w3_array[col] += dx; 777 } 778 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 779 780 /* 781 Evaluate function at x1 + dx (here dx is a vector of perturbations) 782 */ 783 ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 784 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 785 786 /* 787 Loop over rows of vector, putting results into Jacobian matrix 788 */ 789 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 790 for (l=0; l<coloring->nrows[k]; l++) { 791 row = coloring->rows[k][l]; 792 col = coloring->columnsforrow[k][l]; 793 y[row] *= vscale_array[vscaleforrow[k][l]]; 794 srow = row + start; 795 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 796 } 797 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 798 } 799 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 800 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 801 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 802 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 803 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 804 PetscFunctionReturn(0); 805 } 806 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "MatFDColoringSetRecompute()" 810 /*@C 811 MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 812 is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 813 no additional Jacobian's will be computed (the same one will be used) until this is 814 called again. 815 816 Collective on MatFDColoring 817 818 Input Parameters: 819 . c - the coloring context 820 821 Level: intermediate 822 823 Notes: The MatFDColoringSetFrequency() is ignored once this is called 824 825 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 826 827 .keywords: Mat, finite differences, coloring 828 @*/ 829 int MatFDColoringSetRecompute(MatFDColoring c) 830 { 831 PetscFunctionBegin; 832 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 833 c->usersetsrecompute = PETSC_TRUE; 834 c->recompute = PETSC_TRUE; 835 PetscFunctionReturn(0); 836 } 837 838 839