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