1 /*$Id: fdmatrix.c,v 1.85 2001/04/10 19:35:59 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__ "MatFDColoringView_Draw_Zoom" 13 static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 14 { 15 MatFDColoring fd = (MatFDColoring)Aa; 16 int ierr,i,j; 17 PetscReal x,y; 18 19 PetscFunctionBegin; 20 21 /* loop over colors */ 22 for (i=0; i<fd->ncolors; i++) { 23 for (j=0; j<fd->nrows[i]; j++) { 24 y = fd->M - fd->rows[i][j] - fd->rstart; 25 x = fd->columnsforrow[i][j]; 26 ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 27 } 28 } 29 PetscFunctionReturn(0); 30 } 31 32 #undef __FUNCT__ 33 #define __FUNCT__ "MatFDColoringView_Draw" 34 static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 35 { 36 int ierr; 37 PetscTruth isnull; 38 PetscDraw draw; 39 PetscReal xr,yr,xl,yl,h,w; 40 41 PetscFunctionBegin; 42 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 43 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 44 45 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 46 47 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 48 xr += w; yr += h; xl = -w; yl = -h; 49 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 50 ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 51 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 55 #undef __FUNCT__ 56 #define __FUNCT__ "MatFDColoringView" 57 /*@C 58 MatFDColoringView - Views a finite difference coloring context. 59 60 Collective on MatFDColoring 61 62 Input Parameters: 63 + c - the coloring context 64 - viewer - visualization context 65 66 Level: intermediate 67 68 Notes: 69 The available visualization contexts include 70 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 71 . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 72 output where only the first processor opens 73 the file. All other processors send their 74 data to the first processor to print. 75 - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 76 77 Notes: 78 Since PETSc uses only a small number of basic colors (currently 33), if the coloring 79 involves more than 33 then some seemingly identical colors are displayed making it look 80 like an illegal coloring. This is just a graphical artifact. 81 82 .seealso: MatFDColoringCreate() 83 84 .keywords: Mat, finite differences, coloring, view 85 @*/ 86 int MatFDColoringView(MatFDColoring c,PetscViewer viewer) 87 { 88 int i,j,ierr; 89 PetscTruth isdraw,isascii; 90 PetscViewerFormat format; 91 92 PetscFunctionBegin; 93 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 94 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 95 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 96 PetscCheckSameComm(c,viewer); 97 98 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 99 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 100 if (isdraw) { 101 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 102 } else if (isascii) { 103 ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 104 ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr); 105 ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",c->umin);CHKERRQ(ierr); 106 ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors);CHKERRQ(ierr); 107 108 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 109 if (format != PETSC_VIEWER_ASCII_INFO) { 110 for (i=0; i<c->ncolors; i++) { 111 ierr = PetscViewerASCIIPrintf(viewer," Information for color %d\n",i);CHKERRQ(ierr); 112 ierr = PetscViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr); 113 for (j=0; j<c->ncolumns[i]; j++) { 114 ierr = PetscViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]);CHKERRQ(ierr); 115 } 116 ierr = PetscViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr); 117 for (j=0; j<c->nrows[i]; j++) { 118 ierr = PetscViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 119 } 120 } 121 } 122 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 123 } else { 124 SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 125 } 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "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 __FUNCT__ 167 #define __FUNCT__ "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(), MatFDColoringSetRecompute() 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 __FUNCT__ 203 #define __FUNCT__ "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 __FUNCT__ 241 #define __FUNCT__ "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 __FUNCT__ 273 #define __FUNCT__ "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_err <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; 307 308 PetscFunctionBegin; 309 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 310 311 ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 312 ierr = PetscOptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 313 ierr = PetscOptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 314 ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 315 /* not used here; but so they are presented in the GUI */ 316 ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 317 ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 318 ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 319 PetscOptionsEnd();CHKERRQ(ierr); 320 PetscFunctionReturn(0); 321 } 322 323 #undef __FUNCT__ 324 #define __FUNCT__ "MatFDColoringView_Private" 325 int MatFDColoringView_Private(MatFDColoring fd) 326 { 327 int ierr; 328 PetscTruth flg; 329 330 PetscFunctionBegin; 331 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 332 if (flg) { 333 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 334 } 335 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 336 if (flg) { 337 ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 338 ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 340 } 341 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 342 if (flg) { 343 ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 344 ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 345 } 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "MatFDColoringCreate" 351 /*@C 352 MatFDColoringCreate - Creates a matrix coloring context for finite difference 353 computation of Jacobians. 354 355 Collective on Mat 356 357 Input Parameters: 358 + mat - the matrix containing the nonzero structure of the Jacobian 359 - iscoloring - the coloring of the matrix 360 361 Output Parameter: 362 . color - the new coloring context 363 364 Options Database Keys: 365 + -mat_fd_coloring_view - Activates basic viewing or coloring 366 . -mat_fd_coloring_view_draw - Activates drawing of coloring 367 - -mat_fd_coloring_view_info - Activates viewing of coloring info 368 369 Level: intermediate 370 371 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 372 MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 373 @*/ 374 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 375 { 376 MatFDColoring c; 377 MPI_Comm comm; 378 int ierr,M,N; 379 380 PetscFunctionBegin; 381 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 382 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 383 if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 384 385 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 386 PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 387 PetscLogObjectCreate(c); 388 389 if (mat->ops->fdcoloringcreate) { 390 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 391 } else { 392 SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 393 } 394 395 c->error_rel = 1.e-8; 396 c->umin = 1.e-6; 397 c->freq = 1; 398 c->usersetsrecompute = PETSC_FALSE; 399 c->recompute = PETSC_FALSE; 400 401 ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 402 403 *color = c; 404 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 405 PetscFunctionReturn(0); 406 } 407 408 #undef __FUNCT__ 409 #define __FUNCT__ "MatFDColoringDestroy" 410 /*@C 411 MatFDColoringDestroy - Destroys a matrix coloring context that was created 412 via MatFDColoringCreate(). 413 414 Collective on MatFDColoring 415 416 Input Parameter: 417 . c - coloring context 418 419 Level: intermediate 420 421 .seealso: MatFDColoringCreate() 422 @*/ 423 int MatFDColoringDestroy(MatFDColoring c) 424 { 425 int i,ierr; 426 427 PetscFunctionBegin; 428 if (--c->refct > 0) PetscFunctionReturn(0); 429 430 for (i=0; i<c->ncolors; i++) { 431 if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 432 if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 433 if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 434 if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 435 } 436 ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 437 ierr = PetscFree(c->columns);CHKERRQ(ierr); 438 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 439 ierr = PetscFree(c->rows);CHKERRQ(ierr); 440 ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 441 if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 442 if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 443 if (c->w1) { 444 ierr = VecDestroy(c->w1);CHKERRQ(ierr); 445 ierr = VecDestroy(c->w2);CHKERRQ(ierr); 446 ierr = VecDestroy(c->w3);CHKERRQ(ierr); 447 } 448 PetscLogObjectDestroy(c); 449 PetscHeaderDestroy(c); 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "MatFDColoringApply" 455 /*@ 456 MatFDColoringApply - Given a matrix for which a MatFDColoring context 457 has been created, computes the Jacobian for a function via finite differences. 458 459 Collective on MatFDColoring 460 461 Input Parameters: 462 + mat - location to store Jacobian 463 . coloring - coloring context created with MatFDColoringCreate() 464 . x1 - location at which Jacobian is to be computed 465 - sctx - optional context required by function (actually a SNES context) 466 467 Options Database Keys: 468 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 469 470 Level: intermediate 471 472 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 473 474 .keywords: coloring, Jacobian, finite differences 475 @*/ 476 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 477 { 478 int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f; 479 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 480 Scalar dx,mone = -1.0,*y,*xx,*w3_array; 481 Scalar *vscale_array; 482 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 483 Vec w1,w2,w3; 484 void *fctx = coloring->fctx; 485 PetscTruth flg; 486 487 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(J,MAT_COOKIE); 490 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 491 PetscValidHeaderSpecific(x1,VEC_COOKIE); 492 493 if (coloring->usersetsrecompute) { 494 if (!coloring->recompute) { 495 *flag = SAME_PRECONDITIONER; 496 PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"); 497 PetscFunctionReturn(0); 498 } else { 499 coloring->recompute = PETSC_FALSE; 500 } 501 } 502 503 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 504 if (J->ops->fdcoloringapply) { 505 ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 506 } else { 507 508 if (!coloring->w1) { 509 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 510 PetscLogObjectParent(coloring,coloring->w1); 511 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 512 PetscLogObjectParent(coloring,coloring->w2); 513 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 514 PetscLogObjectParent(coloring,coloring->w3); 515 } 516 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 517 518 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 519 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 520 if (flg) { 521 PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 522 } else { 523 ierr = MatZeroEntries(J);CHKERRQ(ierr); 524 } 525 526 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 527 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 528 ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 529 530 /* 531 Compute all the scale factors and share with other processors 532 */ 533 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 534 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 535 for (k=0; k<coloring->ncolors; k++) { 536 /* 537 Loop over each column associated with color adding the 538 perturbation to the vector w3. 539 */ 540 for (l=0; l<coloring->ncolumns[k]; l++) { 541 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 542 dx = xx[col]; 543 if (dx == 0.0) dx = 1.0; 544 #if !defined(PETSC_USE_COMPLEX) 545 if (dx < umin && dx >= 0.0) dx = umin; 546 else if (dx < 0.0 && dx > -umin) dx = -umin; 547 #else 548 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 549 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 550 #endif 551 dx *= epsilon; 552 vscale_array[col] = 1.0/dx; 553 } 554 } 555 vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 556 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 557 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 558 559 /* ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); 560 ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/ 561 562 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 563 else vscaleforrow = coloring->columnsforrow; 564 565 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 566 /* 567 Loop over each color 568 */ 569 for (k=0; k<coloring->ncolors; k++) { 570 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 571 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 572 /* 573 Loop over each column associated with color adding the 574 perturbation to the vector w3. 575 */ 576 for (l=0; l<coloring->ncolumns[k]; l++) { 577 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 578 dx = xx[col]; 579 if (dx == 0.0) dx = 1.0; 580 #if !defined(PETSC_USE_COMPLEX) 581 if (dx < umin && dx >= 0.0) dx = umin; 582 else if (dx < 0.0 && dx > -umin) dx = -umin; 583 #else 584 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 585 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 586 #endif 587 dx *= epsilon; 588 if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter"); 589 w3_array[col] += dx; 590 } 591 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 592 593 /* 594 Evaluate function at x1 + dx (here dx is a vector of perturbations) 595 */ 596 597 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 598 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 599 600 /* 601 Loop over rows of vector, putting results into Jacobian matrix 602 */ 603 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 604 for (l=0; l<coloring->nrows[k]; l++) { 605 row = coloring->rows[k][l]; 606 col = coloring->columnsforrow[k][l]; 607 y[row] *= vscale_array[vscaleforrow[k][l]]; 608 srow = row + start; 609 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 610 } 611 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 612 } 613 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 614 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 615 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 616 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 617 } 618 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 619 620 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 621 if (flg) { 622 ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 623 } 624 PetscFunctionReturn(0); 625 } 626 627 #undef __FUNCT__ 628 #define __FUNCT__ "MatFDColoringApplyTS" 629 /*@ 630 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 631 has been created, computes the Jacobian for a function via finite differences. 632 633 Collective on Mat, MatFDColoring, and Vec 634 635 Input Parameters: 636 + mat - location to store Jacobian 637 . coloring - coloring context created with MatFDColoringCreate() 638 . x1 - location at which Jacobian is to be computed 639 - sctx - optional context required by function (actually a SNES context) 640 641 Options Database Keys: 642 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 643 644 Level: intermediate 645 646 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 647 648 .keywords: coloring, Jacobian, finite differences 649 @*/ 650 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 651 { 652 int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 653 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 654 Scalar dx,mone = -1.0,*y,*xx,*w3_array; 655 Scalar *vscale_array; 656 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 657 Vec w1,w2,w3; 658 void *fctx = coloring->fctx; 659 PetscTruth flg; 660 661 PetscFunctionBegin; 662 PetscValidHeaderSpecific(J,MAT_COOKIE); 663 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 664 PetscValidHeaderSpecific(x1,VEC_COOKIE); 665 666 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 667 if (!coloring->w1) { 668 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 669 PetscLogObjectParent(coloring,coloring->w1); 670 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 671 PetscLogObjectParent(coloring,coloring->w2); 672 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 673 PetscLogObjectParent(coloring,coloring->w3); 674 } 675 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 676 677 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 678 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 679 if (flg) { 680 PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 681 } else { 682 ierr = MatZeroEntries(J);CHKERRQ(ierr); 683 } 684 685 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 686 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 687 ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 688 689 /* 690 Compute all the scale factors and share with other processors 691 */ 692 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 693 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 694 for (k=0; k<coloring->ncolors; k++) { 695 /* 696 Loop over each column associated with color adding the 697 perturbation to the vector w3. 698 */ 699 for (l=0; l<coloring->ncolumns[k]; l++) { 700 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 701 dx = xx[col]; 702 if (dx == 0.0) dx = 1.0; 703 #if !defined(PETSC_USE_COMPLEX) 704 if (dx < umin && dx >= 0.0) dx = umin; 705 else if (dx < 0.0 && dx > -umin) dx = -umin; 706 #else 707 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 708 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 709 #endif 710 dx *= epsilon; 711 vscale_array[col] = 1.0/dx; 712 } 713 } 714 vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 715 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 716 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 717 718 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 719 else vscaleforrow = coloring->columnsforrow; 720 721 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 722 /* 723 Loop over each color 724 */ 725 for (k=0; k<coloring->ncolors; k++) { 726 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 727 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 728 /* 729 Loop over each column associated with color adding the 730 perturbation to the vector w3. 731 */ 732 for (l=0; l<coloring->ncolumns[k]; l++) { 733 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 734 dx = xx[col]; 735 if (dx == 0.0) dx = 1.0; 736 #if !defined(PETSC_USE_COMPLEX) 737 if (dx < umin && dx >= 0.0) dx = umin; 738 else if (dx < 0.0 && dx > -umin) dx = -umin; 739 #else 740 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 741 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 742 #endif 743 dx *= epsilon; 744 w3_array[col] += dx; 745 } 746 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 747 748 /* 749 Evaluate function at x1 + dx (here dx is a vector of perturbations) 750 */ 751 ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 752 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 753 754 /* 755 Loop over rows of vector, putting results into Jacobian matrix 756 */ 757 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 758 for (l=0; l<coloring->nrows[k]; l++) { 759 row = coloring->rows[k][l]; 760 col = coloring->columnsforrow[k][l]; 761 y[row] *= vscale_array[vscaleforrow[k][l]]; 762 srow = row + start; 763 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 764 } 765 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 766 } 767 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 768 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 769 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 770 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 771 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 772 PetscFunctionReturn(0); 773 } 774 775 776 #undef __FUNCT__ 777 #define __FUNCT__ "MatFDColoringSetRecompute()" 778 /*@C 779 MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 780 is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 781 no additional Jacobian's will be computed (the same one will be used) until this is 782 called again. 783 784 Collective on MatFDColoring 785 786 Input Parameters: 787 . c - the coloring context 788 789 Level: intermediate 790 791 Notes: The MatFDColoringSetFrequency() is ignored once this is called 792 793 .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 794 795 .keywords: Mat, finite differences, coloring 796 @*/ 797 int MatFDColoringSetRecompute(MatFDColoring c) 798 { 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 801 c->usersetsrecompute = PETSC_TRUE; 802 c->recompute = PETSC_TRUE; 803 PetscFunctionReturn(0); 804 } 805 806 807