1 2 /* 3 This is where the abstract matrix operations are defined that are 4 used for finite difference computations of Jacobians using coloring. 5 */ 6 7 #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatFDColoringSetF" 11 PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F) 12 { 13 PetscFunctionBegin; 14 fd->F = F; 15 PetscFunctionReturn(0); 16 } 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 20 static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 21 { 22 MatFDColoring fd = (MatFDColoring)Aa; 23 PetscErrorCode ierr; 24 PetscInt i,j; 25 PetscReal x,y; 26 27 PetscFunctionBegin; 28 29 /* loop over colors */ 30 for (i=0; i<fd->ncolors; i++) { 31 for (j=0; j<fd->nrows[i]; j++) { 32 y = fd->M - fd->rows[i][j] - fd->rstart; 33 x = fd->columnsforrow[i][j]; 34 ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 35 } 36 } 37 PetscFunctionReturn(0); 38 } 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "MatFDColoringView_Draw" 42 static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 43 { 44 PetscErrorCode ierr; 45 PetscBool isnull; 46 PetscDraw draw; 47 PetscReal xr,yr,xl,yl,h,w; 48 49 PetscFunctionBegin; 50 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 51 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 52 53 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 54 55 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 56 xr += w; yr += h; xl = -w; yl = -h; 57 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58 ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 59 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "MatFDColoringView" 65 /*@C 66 MatFDColoringView - Views a finite difference coloring context. 67 68 Collective on MatFDColoring 69 70 Input Parameters: 71 + c - the coloring context 72 - viewer - visualization context 73 74 Level: intermediate 75 76 Notes: 77 The available visualization contexts include 78 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 79 . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 80 output where only the first processor opens 81 the file. All other processors send their 82 data to the first processor to print. 83 - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 84 85 Notes: 86 Since PETSc uses only a small number of basic colors (currently 33), if the coloring 87 involves more than 33 then some seemingly identical colors are displayed making it look 88 like an illegal coloring. This is just a graphical artifact. 89 90 .seealso: MatFDColoringCreate() 91 92 .keywords: Mat, finite differences, coloring, view 93 @*/ 94 PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 95 { 96 PetscErrorCode ierr; 97 PetscInt i,j; 98 PetscBool isdraw,iascii; 99 PetscViewerFormat format; 100 101 PetscFunctionBegin; 102 PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 103 if (!viewer) { 104 ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr); 105 } 106 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 107 PetscCheckSameComm(c,1,viewer,2); 108 109 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 110 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 111 if (isdraw) { 112 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 113 } else if (iascii) { 114 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer,"MatFDColoring Object");CHKERRQ(ierr); 115 ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 116 ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 117 ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 118 119 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 120 if (format != PETSC_VIEWER_ASCII_INFO) { 121 for (i=0; i<c->ncolors; i++) { 122 ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 124 for (j=0; j<c->ncolumns[i]; j++) { 125 ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 126 } 127 ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 128 for (j=0; j<c->nrows[i]; j++) { 129 ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 130 } 131 } 132 } 133 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 134 } else { 135 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 136 } 137 PetscFunctionReturn(0); 138 } 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "MatFDColoringSetParameters" 142 /*@ 143 MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 144 a Jacobian matrix using finite differences. 145 146 Logically Collective on MatFDColoring 147 148 The Jacobian is estimated with the differencing approximation 149 .vb 150 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 151 h = error_rel*u[i] if abs(u[i]) > umin 152 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 153 dx_{i} = (0, ... 1, .... 0) 154 .ve 155 156 Input Parameters: 157 + coloring - the coloring context 158 . error_rel - relative error 159 - umin - minimum allowable u-value magnitude 160 161 Level: advanced 162 163 .keywords: Mat, finite differences, coloring, set, parameters 164 165 .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 166 167 @*/ 168 PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 169 { 170 PetscFunctionBegin; 171 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 172 PetscValidLogicalCollectiveReal(matfd,error,2); 173 PetscValidLogicalCollectiveReal(matfd,umin,3); 174 if (error != PETSC_DEFAULT) matfd->error_rel = error; 175 if (umin != PETSC_DEFAULT) matfd->umin = umin; 176 PetscFunctionReturn(0); 177 } 178 179 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "MatFDColoringGetFunction" 183 /*@C 184 MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 185 186 Not Collective 187 188 Input Parameters: 189 . coloring - the coloring context 190 191 Output Parameters: 192 + f - the function 193 - fctx - the optional user-defined function context 194 195 Level: intermediate 196 197 .keywords: Mat, Jacobian, finite differences, set, function 198 199 .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 200 201 @*/ 202 PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 203 { 204 PetscFunctionBegin; 205 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 206 if (f) *f = matfd->f; 207 if (fctx) *fctx = matfd->fctx; 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNCT__ 212 #define __FUNCT__ "MatFDColoringSetFunction" 213 /*@C 214 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 215 216 Logically Collective on MatFDColoring 217 218 Input Parameters: 219 + coloring - the coloring context 220 . f - the function 221 - fctx - the optional user-defined function context 222 223 Calling sequence of (*f) function: 224 For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 225 If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 226 227 Level: advanced 228 229 Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 230 SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by 231 calling MatFDColoringApply() 232 233 Fortran Notes: 234 In Fortran you must call MatFDColoringSetFunction() for a coloring object to 235 be used without SNES or within the SNES solvers. 236 237 .keywords: Mat, Jacobian, finite differences, set, function 238 239 .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 240 241 @*/ 242 PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 243 { 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 246 matfd->f = f; 247 matfd->fctx = fctx; 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "MatFDColoringSetFromOptions" 253 /*@ 254 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 255 the options database. 256 257 Collective on MatFDColoring 258 259 The Jacobian, F'(u), is estimated with the differencing approximation 260 .vb 261 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 262 h = error_rel*u[i] if abs(u[i]) > umin 263 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 264 dx_{i} = (0, ... 1, .... 0) 265 .ve 266 267 Input Parameter: 268 . coloring - the coloring context 269 270 Options Database Keys: 271 + -mat_fd_coloring_err <err> - Sets <err> (square root 272 of relative error in the function) 273 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 274 . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 275 . -mat_fd_coloring_view - Activates basic viewing 276 . -mat_fd_coloring_view_info - Activates viewing info 277 - -mat_fd_coloring_view_draw - Activates drawing 278 279 Level: intermediate 280 281 .keywords: Mat, finite differences, parameters 282 283 .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 284 285 @*/ 286 PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 287 { 288 PetscErrorCode ierr; 289 PetscBool flg; 290 char value[3]; 291 292 PetscFunctionBegin; 293 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 294 295 ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr); 296 ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 297 ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 298 ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr); 299 if (flg) { 300 if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 301 else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 302 else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 303 } 304 /* not used here; but so they are presented in the GUI */ 305 ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 306 ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 307 ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 308 309 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 310 ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr); 311 PetscOptionsEnd();CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "MatFDColoringView_Private" 317 PetscErrorCode MatFDColoringView_Private(MatFDColoring fd) 318 { 319 PetscErrorCode ierr; 320 PetscBool flg = PETSC_FALSE; 321 PetscViewer viewer; 322 323 PetscFunctionBegin; 324 ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr); 325 ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr); 326 if (flg) { 327 ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 328 } 329 flg = PETSC_FALSE; 330 ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr); 331 if (flg) { 332 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 333 ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 334 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 335 } 336 flg = PETSC_FALSE; 337 ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr); 338 if (flg) { 339 ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 340 ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 341 } 342 PetscFunctionReturn(0); 343 } 344 345 #undef __FUNCT__ 346 #define __FUNCT__ "MatFDColoringCreate" 347 /*@ 348 MatFDColoringCreate - Creates a matrix coloring context for finite difference 349 computation of Jacobians. 350 351 Collective on Mat 352 353 Input Parameters: 354 + mat - the matrix containing the nonzero structure of the Jacobian 355 - iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 356 357 Output Parameter: 358 . color - the new coloring context 359 360 Level: intermediate 361 362 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 363 MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 364 MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring() 365 @*/ 366 PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 367 { 368 MatFDColoring c; 369 MPI_Comm comm; 370 PetscErrorCode ierr; 371 PetscInt M,N; 372 373 PetscFunctionBegin; 374 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 375 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 376 if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices"); 377 378 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 379 ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 380 381 c->ctype = iscoloring->ctype; 382 383 if (mat->ops->fdcoloringcreate) { 384 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 385 } else SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 386 387 ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr); 388 ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr); 389 ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 390 ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr); 391 392 c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 393 c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 394 c->currentcolor = -1; 395 c->htype = "wp"; 396 397 *color = c; 398 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 399 PetscFunctionReturn(0); 400 } 401 402 #undef __FUNCT__ 403 #define __FUNCT__ "MatFDColoringDestroy" 404 /*@ 405 MatFDColoringDestroy - Destroys a matrix coloring context that was created 406 via MatFDColoringCreate(). 407 408 Collective on MatFDColoring 409 410 Input Parameter: 411 . c - coloring context 412 413 Level: intermediate 414 415 .seealso: MatFDColoringCreate() 416 @*/ 417 PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 418 { 419 PetscErrorCode ierr; 420 PetscInt i; 421 422 PetscFunctionBegin; 423 if (!*c) PetscFunctionReturn(0); 424 if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);} 425 426 for (i=0; i<(*c)->ncolors; i++) { 427 ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr); 428 ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr); 429 ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr); 430 if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);} 431 } 432 ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr); 433 ierr = PetscFree((*c)->columns);CHKERRQ(ierr); 434 ierr = PetscFree((*c)->nrows);CHKERRQ(ierr); 435 ierr = PetscFree((*c)->rows);CHKERRQ(ierr); 436 ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr); 437 ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr); 438 ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr); 439 ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr); 440 ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr); 441 ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr); 442 ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 446 #undef __FUNCT__ 447 #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 448 /*@C 449 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 450 that are currently being perturbed. 451 452 Not Collective 453 454 Input Parameters: 455 . coloring - coloring context created with MatFDColoringCreate() 456 457 Output Parameters: 458 + n - the number of local columns being perturbed 459 - cols - the column indices, in global numbering 460 461 Level: intermediate 462 463 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 464 465 .keywords: coloring, Jacobian, finite differences 466 @*/ 467 PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 468 { 469 PetscFunctionBegin; 470 if (coloring->currentcolor >= 0) { 471 *n = coloring->ncolumns[coloring->currentcolor]; 472 *cols = coloring->columns[coloring->currentcolor]; 473 } else { 474 *n = 0; 475 } 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "MatFDColoringApply" 481 /*@ 482 MatFDColoringApply - Given a matrix for which a MatFDColoring context 483 has been created, computes the Jacobian for a function via finite differences. 484 485 Collective on MatFDColoring 486 487 Input Parameters: 488 + mat - location to store Jacobian 489 . coloring - coloring context created with MatFDColoringCreate() 490 . x1 - location at which Jacobian is to be computed 491 - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 492 493 Options Database Keys: 494 + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 495 . -mat_fd_coloring_view - Activates basic viewing or coloring 496 . -mat_fd_coloring_view_draw - Activates drawing of coloring 497 - -mat_fd_coloring_view_info - Activates viewing of coloring info 498 499 Level: intermediate 500 501 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 502 503 .keywords: coloring, Jacobian, finite differences 504 @*/ 505 PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 506 { 507 PetscErrorCode ierr; 508 509 PetscFunctionBegin; 510 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 511 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 512 PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 513 if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 514 if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 515 ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 516 PetscFunctionReturn(0); 517 } 518 519 #undef __FUNCT__ 520 #define __FUNCT__ "MatFDColoringApply_AIJ" 521 PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 522 { 523 PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; 524 PetscErrorCode ierr; 525 PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 526 PetscScalar dx,*y,*xx,*w3_array; 527 PetscScalar *vscale_array; 528 PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 529 Vec w1=coloring->w1,w2=coloring->w2,w3; 530 void *fctx = coloring->fctx; 531 PetscBool flg = PETSC_FALSE; 532 PetscInt ctype=coloring->ctype,N,col_start=0,col_end=0; 533 Vec x1_tmp; 534 535 PetscFunctionBegin; 536 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 537 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 538 PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 539 if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 540 541 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 542 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 543 ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr); 544 if (flg) { 545 ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 546 } else { 547 PetscBool assembled; 548 ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 549 if (assembled) { 550 ierr = MatZeroEntries(J);CHKERRQ(ierr); 551 } 552 } 553 554 x1_tmp = x1; 555 if (!coloring->vscale){ 556 ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 557 } 558 559 /* 560 This is a horrible, horrible, hack. 561 */ 562 if (coloring->F) { 563 ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 564 ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 565 if (m1 != m2) { 566 coloring->F = 0; 567 } 568 } 569 570 if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 571 ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 572 } 573 ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 574 575 /* Set w1 = F(x1) */ 576 if (coloring->F) { 577 w1 = coloring->F; /* use already computed value of function */ 578 coloring->F = 0; 579 } else { 580 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 581 ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 582 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 583 } 584 585 if (!coloring->w3) { 586 ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 587 ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 588 } 589 w3 = coloring->w3; 590 591 /* Compute all the local scale factors, including ghost points */ 592 ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 593 ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 594 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 595 if (ctype == IS_COLORING_GHOSTED){ 596 col_start = 0; col_end = N; 597 } else if (ctype == IS_COLORING_GLOBAL){ 598 xx = xx - start; 599 vscale_array = vscale_array - start; 600 col_start = start; col_end = N + start; 601 } 602 for (col=col_start; col<col_end; col++){ 603 /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 604 if (coloring->htype[0] == 'w') { 605 dx = 1.0 + unorm; 606 } else { 607 dx = xx[col]; 608 } 609 if (dx == (PetscScalar)0.0) dx = 1.0; 610 #if !defined(PETSC_USE_COMPLEX) 611 if (dx < umin && dx >= 0.0) dx = umin; 612 else if (dx < 0.0 && dx > -umin) dx = -umin; 613 #else 614 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 615 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 616 #endif 617 dx *= epsilon; 618 vscale_array[col] = (PetscScalar)1.0/dx; 619 } 620 if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start; 621 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 622 if (ctype == IS_COLORING_GLOBAL){ 623 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 624 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 625 } 626 627 if (coloring->vscaleforrow) { 628 vscaleforrow = coloring->vscaleforrow; 629 } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 630 631 /* 632 Loop over each color 633 */ 634 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 635 for (k=0; k<coloring->ncolors; k++) { 636 coloring->currentcolor = k; 637 ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 638 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 639 if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 640 /* 641 Loop over each column associated with color 642 adding the perturbation to the vector w3. 643 */ 644 for (l=0; l<coloring->ncolumns[k]; l++) { 645 col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 646 if (coloring->htype[0] == 'w') { 647 dx = 1.0 + unorm; 648 } else { 649 dx = xx[col]; 650 } 651 if (dx == (PetscScalar)0.0) dx = 1.0; 652 #if !defined(PETSC_USE_COMPLEX) 653 if (dx < umin && dx >= 0.0) dx = umin; 654 else if (dx < 0.0 && dx > -umin) dx = -umin; 655 #else 656 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 657 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 658 #endif 659 dx *= epsilon; 660 if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 661 w3_array[col] += dx; 662 } 663 if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 664 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 665 666 /* 667 Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 668 w2 = F(x1 + dx) - F(x1) 669 */ 670 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 671 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 672 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 673 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 674 675 /* 676 Loop over rows of vector, putting results into Jacobian matrix 677 */ 678 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 679 for (l=0; l<coloring->nrows[k]; l++) { 680 row = coloring->rows[k][l]; /* local row index */ 681 col = coloring->columnsforrow[k][l]; /* global column index */ 682 y[row] *= vscale_array[vscaleforrow[k][l]]; 683 srow = row + start; 684 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 685 } 686 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 687 } /* endof for each color */ 688 if (ctype == IS_COLORING_GLOBAL) xx = xx + start; 689 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 690 ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 691 692 coloring->currentcolor = -1; 693 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 694 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 695 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 696 697 flg = PETSC_FALSE; 698 ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr); 699 if (flg) { 700 ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr); 701 } 702 ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 703 PetscFunctionReturn(0); 704 } 705