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