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