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);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 SNESComputeJacobianDefaultColor()) 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 prefix[],const char optionname[]) 318 { 319 PetscErrorCode ierr; 320 PetscBool flg; 321 PetscViewer viewer; 322 PetscViewerFormat format; 323 324 PetscFunctionBegin; 325 if (prefix) { 326 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 327 } else { 328 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 329 } 330 if (flg) { 331 ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 332 ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 333 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 334 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 335 } 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "MatFDColoringCreate" 341 /*@ 342 MatFDColoringCreate - Creates a matrix coloring context for finite difference 343 computation of Jacobians. 344 345 Collective on Mat 346 347 Input Parameters: 348 + mat - the matrix containing the nonzero structure of the Jacobian 349 - iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 350 351 Output Parameter: 352 . color - the new coloring context 353 354 Level: intermediate 355 356 .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(), 357 MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 358 MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring() 359 @*/ 360 PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 361 { 362 MatFDColoring c; 363 MPI_Comm comm; 364 PetscErrorCode ierr; 365 PetscInt M,N; 366 367 PetscFunctionBegin; 368 ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 369 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 370 if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices"); 371 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 372 ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 373 374 c->ctype = iscoloring->ctype; 375 376 if (mat->ops->fdcoloringcreate) { 377 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 378 } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 379 380 ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr); 381 ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr); 382 ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 383 ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr); 384 385 c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 386 c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 387 c->currentcolor = -1; 388 c->htype = "wp"; 389 c->fset = PETSC_FALSE; 390 391 *color = c; 392 ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr); 393 ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 394 PetscFunctionReturn(0); 395 } 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "MatFDColoringDestroy" 399 /*@ 400 MatFDColoringDestroy - Destroys a matrix coloring context that was created 401 via MatFDColoringCreate(). 402 403 Collective on MatFDColoring 404 405 Input Parameter: 406 . c - coloring context 407 408 Level: intermediate 409 410 .seealso: MatFDColoringCreate() 411 @*/ 412 PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 413 { 414 PetscErrorCode ierr; 415 PetscInt i; 416 417 PetscFunctionBegin; 418 if (!*c) PetscFunctionReturn(0); 419 if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);} 420 421 for (i=0; i<(*c)->ncolors; i++) { 422 ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr); 423 ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr); 424 ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr); 425 if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);} 426 } 427 ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr); 428 ierr = PetscFree((*c)->columns);CHKERRQ(ierr); 429 ierr = PetscFree((*c)->nrows);CHKERRQ(ierr); 430 ierr = PetscFree((*c)->rows);CHKERRQ(ierr); 431 ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr); 432 ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr); 433 ierr = PetscFree((*c)->den2sp);CHKERRQ(ierr); 434 ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr); 435 ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr); 436 ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr); 437 ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr); 438 ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 444 /*@C 445 MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 446 that are currently being perturbed. 447 448 Not Collective 449 450 Input Parameters: 451 . coloring - coloring context created with MatFDColoringCreate() 452 453 Output Parameters: 454 + n - the number of local columns being perturbed 455 - cols - the column indices, in global numbering 456 457 Level: intermediate 458 459 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 460 461 .keywords: coloring, Jacobian, finite differences 462 @*/ 463 PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 464 { 465 PetscFunctionBegin; 466 if (coloring->currentcolor >= 0) { 467 *n = coloring->ncolumns[coloring->currentcolor]; 468 *cols = coloring->columns[coloring->currentcolor]; 469 } else { 470 *n = 0; 471 } 472 PetscFunctionReturn(0); 473 } 474 475 #undef __FUNCT__ 476 #define __FUNCT__ "MatFDColoringApply" 477 /*@ 478 MatFDColoringApply - Given a matrix for which a MatFDColoring context 479 has been created, computes the Jacobian for a function via finite differences. 480 481 Collective on MatFDColoring 482 483 Input Parameters: 484 + mat - location to store Jacobian 485 . coloring - coloring context created with MatFDColoringCreate() 486 . x1 - location at which Jacobian is to be computed 487 - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 488 489 Options Database Keys: 490 + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 491 . -mat_fd_coloring_view - Activates basic viewing or coloring 492 . -mat_fd_coloring_view draw - Activates drawing of coloring 493 - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 494 495 Level: intermediate 496 497 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 498 499 .keywords: coloring, Jacobian, finite differences 500 @*/ 501 PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 502 { 503 PetscErrorCode ierr; 504 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 507 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 508 PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 509 if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 510 if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 511 ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 512 PetscFunctionReturn(0); 513 } 514 515 /* #define JACOBIANCOLOROPT */ 516 #if defined(JACOBIANCOLOROPT) 517 #include <petsctime.h> 518 #endif 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; 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 #if defined(JACOBIANCOLOROPT) 535 PetscLogDouble t0,t1,time_setvalues=0.0; 536 #endif 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 540 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 541 PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 542 if (!f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 543 544 ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 545 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 546 ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 547 if (flg) { 548 ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 549 } else { 550 PetscBool assembled; 551 ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 552 if (assembled) { 553 ierr = MatZeroEntries(J);CHKERRQ(ierr); 554 } 555 } 556 557 x1_tmp = x1; 558 if (!coloring->vscale) { 559 ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 560 } 561 562 if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 563 ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 564 } 565 ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 566 567 /* Set w1 = F(x1) */ 568 if (!coloring->fset) { 569 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 570 ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 571 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 572 } else { 573 coloring->fset = PETSC_FALSE; 574 } 575 576 if (!coloring->w3) { 577 ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 578 ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 579 } 580 w3 = coloring->w3; 581 582 /* Compute all the local scale factors, including ghost points */ 583 ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 584 ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 585 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 586 if (ctype == IS_COLORING_GHOSTED) { 587 col_start = 0; col_end = N; 588 } else if (ctype == IS_COLORING_GLOBAL) { 589 xx = xx - start; 590 vscale_array = vscale_array - start; 591 col_start = start; col_end = N + start; 592 } 593 for (col=col_start; col<col_end; col++) { 594 /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 595 if (coloring->htype[0] == 'w') { 596 dx = 1.0 + unorm; 597 } else { 598 dx = xx[col]; 599 } 600 if (dx == (PetscScalar)0.0) dx = 1.0; 601 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 602 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 603 dx *= epsilon; 604 vscale_array[col] = (PetscScalar)1.0/dx; 605 } 606 if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start; 607 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 608 if (ctype == IS_COLORING_GLOBAL) { 609 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 610 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 611 } 612 613 if (coloring->vscaleforrow) { 614 vscaleforrow = coloring->vscaleforrow; 615 } else SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 616 617 /* 618 Loop over each color 619 */ 620 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 621 for (k=0; k<coloring->ncolors; k++) { 622 coloring->currentcolor = k; 623 624 ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 625 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 626 if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 627 /* 628 Loop over each column associated with color 629 adding the perturbation to the vector w3. 630 */ 631 for (l=0; l<coloring->ncolumns[k]; l++) { 632 col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 633 if (coloring->htype[0] == 'w') { 634 dx = 1.0 + unorm; 635 } else { 636 dx = xx[col]; 637 } 638 if (dx == (PetscScalar)0.0) dx = 1.0; 639 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 640 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 641 dx *= epsilon; 642 if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 643 w3_array[col] += dx; 644 } 645 if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 646 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 647 648 /* 649 Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 650 w2 = F(x1 + dx) - F(x1) 651 */ 652 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 653 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 654 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 655 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 656 657 /* 658 Loop over rows of vector, putting results into Jacobian matrix 659 */ 660 #if defined(JACOBIANCOLOROPT) 661 ierr = PetscTime(&t0);CHKERRQ(ierr); 662 #endif 663 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 664 for (l=0; l<coloring->nrows[k]; l++) { 665 row = coloring->rows[k][l]; /* local row index */ 666 col = coloring->columnsforrow[k][l]; /* global column index */ 667 y[row] *= vscale_array[vscaleforrow[k][l]]; 668 srow = row + start; 669 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 670 } 671 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 672 #if defined(JACOBIANCOLOROPT) 673 ierr = PetscTime(&t1);CHKERRQ(ierr); 674 time_setvalues += t1-t0; 675 #endif 676 } /* endof for each color */ 677 #if defined(JACOBIANCOLOROPT) 678 printf(" MatFDColoringApply_AIJ: time_setvalues %g\n",time_setvalues); 679 #endif 680 if (ctype == IS_COLORING_GLOBAL) xx = xx + start; 681 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 682 ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 683 684 coloring->currentcolor = -1; 685 686 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 687 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 688 ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 689 690 ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr); 691 PetscFunctionReturn(0); 692 } 693