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