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