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