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