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