1 /*$Id: fdmatrix.c,v 1.75 2000/08/25 16:45:21 balay Exp balay $*/ 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 "petsc.h" 9 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 10 11 #undef __FUNC__ 12 #define __FUNC__ /*<a name="MatFDColoringView_Draw_Zoom"></a>*/"MatFDColoringView_Draw_Zoom" 13 static int MatFDColoringView_Draw_Zoom(Draw draw,void *Aa) 14 { 15 MatFDColoring fd = (MatFDColoring)Aa; 16 int ierr,i,j; 17 PetscReal x,y; 18 19 PetscFunctionBegin; 20 21 /* loop over colors */ 22 for (i=0; i<fd->ncolors; i++) { 23 for (j=0; j<fd->nrows[i]; j++) { 24 y = fd->M - fd->rows[i][j] - fd->rstart; 25 x = fd->columnsforrow[i][j]; 26 ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 27 } 28 } 29 PetscFunctionReturn(0); 30 } 31 32 #undef __FUNC__ 33 #define __FUNC__ /*<a name="MatFDColoringView_Draw"></a>*/"MatFDColoringView_Draw" 34 static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer) 35 { 36 int ierr; 37 PetscTruth isnull; 38 Draw draw; 39 PetscReal xr,yr,xl,yl,h,w; 40 41 PetscFunctionBegin; 42 ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 43 ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 44 45 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 46 47 xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 48 xr += w; yr += h; xl = -w; yl = -h; 49 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 50 ierr = DrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 51 ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 55 #undef __FUNC__ 56 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView" 57 /*@C 58 MatFDColoringView - Views a finite difference coloring context. 59 60 Collective on MatFDColoring 61 62 Input Parameters: 63 + c - the coloring context 64 - viewer - visualization context 65 66 Level: intermediate 67 68 Notes: 69 The available visualization contexts include 70 + VIEWER_STDOUT_SELF - standard output (default) 71 . VIEWER_STDOUT_WORLD - synchronized standard 72 output where only the first processor opens 73 the file. All other processors send their 74 data to the first processor to print. 75 - VIEWER_DRAW_WORLD - graphical display of nonzero structure 76 77 Notes: 78 Since PETSc uses only a small number of basic colors (currently 33), if the coloring 79 involves moe than 33 then some seemingly identical colors are displayed making it look 80 like an illegal coloring. This is just a graphical artifact. 81 82 .seealso: MatFDColoringCreate() 83 84 .keywords: Mat, finite differences, coloring, view 85 @*/ 86 int MatFDColoringView(MatFDColoring c,Viewer viewer) 87 { 88 int i,j,format,ierr; 89 PetscTruth isdraw,isascii; 90 91 PetscFunctionBegin; 92 PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 93 if (!viewer) viewer = VIEWER_STDOUT_(c->comm); 94 PetscValidHeaderSpecific(viewer,VIEWER_COOKIE); 95 PetscCheckSameComm(c,viewer); 96 97 ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 98 ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 99 if (isdraw) { 100 ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 101 } else if (isascii) { 102 ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 103 ierr = ViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr); 104 ierr = ViewerASCIIPrintf(viewer," Umin=%g\n",c->umin);CHKERRQ(ierr); 105 ierr = ViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors);CHKERRQ(ierr); 106 107 ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 108 if (format != VIEWER_FORMAT_ASCII_INFO) { 109 for (i=0; i<c->ncolors; i++) { 110 ierr = ViewerASCIIPrintf(viewer," Information for color %d\n",i);CHKERRQ(ierr); 111 ierr = ViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr); 112 for (j=0; j<c->ncolumns[i]; j++) { 113 ierr = ViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]);CHKERRQ(ierr); 114 } 115 ierr = ViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr); 116 for (j=0; j<c->nrows[i]; j++) { 117 ierr = ViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 118 } 119 } 120 } 121 ierr = ViewerFlush(viewer);CHKERRQ(ierr); 122 } else { 123 SETERRQ1(1,1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 124 } 125 PetscFunctionReturn(0); 126 } 127 128 #undef __FUNC__ 129 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetParameters" 130 /*@ 131 MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 132 a Jacobian matrix using finite differences. 133 134 Collective on MatFDColoring 135 136 The Jacobian is estimated with the differencing approximation 137 .vb 138 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 139 h = error_rel*u[i] if abs(u[i]) > umin 140 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 141 dx_{i} = (0, ... 1, .... 0) 142 .ve 143 144 Input Parameters: 145 + coloring - the coloring context 146 . error_rel - relative error 147 - umin - minimum allowable u-value magnitude 148 149 Level: advanced 150 151 .keywords: Mat, finite differences, coloring, set, parameters 152 153 .seealso: MatFDColoringCreate() 154 @*/ 155 int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 156 { 157 PetscFunctionBegin; 158 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 159 160 if (error != PETSC_DEFAULT) matfd->error_rel = error; 161 if (umin != PETSC_DEFAULT) matfd->umin = umin; 162 PetscFunctionReturn(0); 163 } 164 165 #undef __FUNC__ 166 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFrequency" 167 /*@ 168 MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 169 matrices. 170 171 Collective on MatFDColoring 172 173 Input Parameters: 174 + coloring - the coloring context 175 - freq - frequency (default is 1) 176 177 Options Database Keys: 178 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 179 180 Level: advanced 181 182 Notes: 183 Using a modified Newton strategy, where the Jacobian remains fixed for several 184 iterations, can be cost effective in terms of overall nonlinear solution 185 efficiency. This parameter indicates that a new Jacobian will be computed every 186 <freq> nonlinear iterations. 187 188 .keywords: Mat, finite differences, coloring, set, frequency 189 190 .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency() 191 @*/ 192 int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 193 { 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 196 197 matfd->freq = freq; 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNC__ 202 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringGetFrequency" 203 /*@ 204 MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 205 matrices. 206 207 Not Collective 208 209 Input Parameters: 210 . coloring - the coloring context 211 212 Output Parameters: 213 . freq - frequency (default is 1) 214 215 Options Database Keys: 216 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 217 218 Level: advanced 219 220 Notes: 221 Using a modified Newton strategy, where the Jacobian remains fixed for several 222 iterations, can be cost effective in terms of overall nonlinear solution 223 efficiency. This parameter indicates that a new Jacobian will be computed every 224 <freq> nonlinear iterations. 225 226 .keywords: Mat, finite differences, coloring, get, frequency 227 228 .seealso: MatFDColoringSetFrequency() 229 @*/ 230 int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 231 { 232 PetscFunctionBegin; 233 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 234 235 *freq = matfd->freq; 236 PetscFunctionReturn(0); 237 } 238 239 #undef __FUNC__ 240 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFunction" 241 /*@C 242 MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 243 244 Collective on MatFDColoring 245 246 Input Parameters: 247 + coloring - the coloring context 248 . f - the function 249 - fctx - the optional user-defined function context 250 251 Level: intermediate 252 253 Notes: 254 In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 255 be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 256 with the TS solvers. 257 258 .keywords: Mat, Jacobian, finite differences, set, function 259 @*/ 260 int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 261 { 262 PetscFunctionBegin; 263 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 264 265 matfd->f = f; 266 matfd->fctx = fctx; 267 268 PetscFunctionReturn(0); 269 } 270 271 #undef __FUNC__ 272 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFromOptions" 273 /*@ 274 MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 275 the options database. 276 277 Collective on MatFDColoring 278 279 The Jacobian, F'(u), is estimated with the differencing approximation 280 .vb 281 F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 282 h = error_rel*u[i] if abs(u[i]) > umin 283 = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 284 dx_{i} = (0, ... 1, .... 0) 285 .ve 286 287 Input Parameter: 288 . coloring - the coloring context 289 290 Options Database Keys: 291 + -mat_fd_coloring_err <err> - Sets <err> (square root 292 of relative error in the function) 293 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 294 . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 295 . -mat_fd_coloring_view - Activates basic viewing 296 . -mat_fd_coloring_view_info - Activates viewing info 297 - -mat_fd_coloring_view_draw - Activates drawing 298 299 Level: intermediate 300 301 .keywords: Mat, finite differences, parameters 302 @*/ 303 int MatFDColoringSetFromOptions(MatFDColoring matfd) 304 { 305 int ierr; 306 307 PetscFunctionBegin; 308 PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 309 310 ierr = OptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option");CHKERRQ(ierr); 311 ierr = OptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 312 ierr = OptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 313 ierr = OptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 314 /* not used here; but so they are presented in the GUI */ 315 ierr = OptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 316 ierr = OptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 317 ierr = OptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 318 OptionsEnd();CHKERRQ(ierr); 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNC__ 323 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView_Private" 324 int MatFDColoringView_Private(MatFDColoring fd) 325 { 326 int ierr; 327 PetscTruth flg; 328 329 PetscFunctionBegin; 330 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 331 if (flg) { 332 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 333 } 334 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 335 if (flg) { 336 ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 337 ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 338 ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339 } 340 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 341 if (flg) { 342 ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 343 ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 344 } 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNC__ 349 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringCreate" 350 /*@C 351 MatFDColoringCreate - Creates a matrix coloring context for finite difference 352 computation of Jacobians. 353 354 Collective on Mat 355 356 Input Parameters: 357 + mat - the matrix containing the nonzero structure of the Jacobian 358 - iscoloring - the coloring of the matrix 359 360 Output Parameter: 361 . color - the new coloring context 362 363 Options Database Keys: 364 + -mat_fd_coloring_view - Activates basic viewing or coloring 365 . -mat_fd_coloring_view_draw - Activates drawing of coloring 366 - -mat_fd_coloring_view_info - Activates viewing of coloring info 367 368 Level: intermediate 369 370 .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 371 MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 372 @*/ 373 int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 374 { 375 MatFDColoring c; 376 MPI_Comm comm; 377 int ierr,M,N; 378 379 PetscFunctionBegin; 380 ierr = PLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 381 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 382 if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 383 384 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 385 PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 386 PLogObjectCreate(c); 387 388 if (mat->ops->fdcoloringcreate) { 389 ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 390 } else { 391 SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 392 } 393 394 c->error_rel = 1.e-8; 395 c->umin = 1.e-6; 396 c->freq = 1; 397 398 ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 399 400 *color = c; 401 ierr = PLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 402 PetscFunctionReturn(0); 403 } 404 405 #undef __FUNC__ 406 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringDestroy" 407 /*@C 408 MatFDColoringDestroy - Destroys a matrix coloring context that was created 409 via MatFDColoringCreate(). 410 411 Collective on MatFDColoring 412 413 Input Parameter: 414 . c - coloring context 415 416 Level: intermediate 417 418 .seealso: MatFDColoringCreate() 419 @*/ 420 int MatFDColoringDestroy(MatFDColoring c) 421 { 422 int i,ierr; 423 424 PetscFunctionBegin; 425 if (--c->refct > 0) PetscFunctionReturn(0); 426 427 for (i=0; i<c->ncolors; i++) { 428 if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 429 if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 430 if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 431 if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 432 } 433 ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 434 ierr = PetscFree(c->columns);CHKERRQ(ierr); 435 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 436 ierr = PetscFree(c->rows);CHKERRQ(ierr); 437 ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 438 if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 439 if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 440 if (c->w1) { 441 ierr = VecDestroy(c->w1);CHKERRQ(ierr); 442 ierr = VecDestroy(c->w2);CHKERRQ(ierr); 443 ierr = VecDestroy(c->w3);CHKERRQ(ierr); 444 } 445 PLogObjectDestroy(c); 446 PetscHeaderDestroy(c); 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNC__ 451 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApply" 452 /*@ 453 MatFDColoringApply - Given a matrix for which a MatFDColoring context 454 has been created, computes the Jacobian for a function via finite differences. 455 456 Collective on MatFDColoring 457 458 Input Parameters: 459 + mat - location to store Jacobian 460 . coloring - coloring context created with MatFDColoringCreate() 461 . x1 - location at which Jacobian is to be computed 462 - sctx - optional context required by function (actually a SNES context) 463 464 Options Database Keys: 465 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 466 467 Level: intermediate 468 469 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 470 471 .keywords: coloring, Jacobian, finite differences 472 @*/ 473 int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 474 { 475 int (*f)(void *,Vec,Vec,void*)= (int (*)(void *,Vec,Vec,void *))coloring->f; 476 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 477 Scalar dx,mone = -1.0,*y,*xx,*w3_array; 478 Scalar *vscale_array; 479 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 480 Vec w1,w2,w3; 481 void *fctx = coloring->fctx; 482 PetscTruth flg; 483 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(J,MAT_COOKIE); 487 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 488 PetscValidHeaderSpecific(x1,VEC_COOKIE); 489 490 ierr = PLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 491 if (!coloring->w1) { 492 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 493 PLogObjectParent(coloring,coloring->w1); 494 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 495 PLogObjectParent(coloring,coloring->w2); 496 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 497 PLogObjectParent(coloring,coloring->w3); 498 } 499 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 500 501 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 502 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 503 if (flg) { 504 PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 505 } else { 506 ierr = MatZeroEntries(J);CHKERRQ(ierr); 507 } 508 509 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 510 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 511 ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 512 513 /* 514 Compute all the scale factors and share with other processors 515 */ 516 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 517 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 518 for (k=0; k<coloring->ncolors; k++) { 519 /* 520 Loop over each column associated with color adding the 521 perturbation to the vector w3. 522 */ 523 for (l=0; l<coloring->ncolumns[k]; l++) { 524 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 525 dx = xx[col]; 526 if (dx == 0.0) dx = 1.0; 527 #if !defined(PETSC_USE_COMPLEX) 528 if (dx < umin && dx >= 0.0) dx = umin; 529 else if (dx < 0.0 && dx > -umin) dx = -umin; 530 #else 531 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 532 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 533 #endif 534 dx *= epsilon; 535 vscale_array[col] = 1.0/dx; 536 } 537 } 538 vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 539 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 540 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 541 542 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 543 else vscaleforrow = coloring->columnsforrow; 544 545 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 546 /* 547 Loop over each color 548 */ 549 for (k=0; k<coloring->ncolors; k++) { 550 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 551 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 552 /* 553 Loop over each column associated with color adding the 554 perturbation to the vector w3. 555 */ 556 for (l=0; l<coloring->ncolumns[k]; l++) { 557 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 558 dx = xx[col]; 559 if (dx == 0.0) dx = 1.0; 560 #if !defined(PETSC_USE_COMPLEX) 561 if (dx < umin && dx >= 0.0) dx = umin; 562 else if (dx < 0.0 && dx > -umin) dx = -umin; 563 #else 564 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 565 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 566 #endif 567 dx *= epsilon; 568 if (!PetscAbsScalar(dx)) SETERRQ(1,1,"Computed 0 differencing parameter"); 569 w3_array[col] += dx; 570 } 571 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 572 573 /* 574 Evaluate function at x1 + dx (here dx is a vector of perturbations) 575 */ 576 577 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 578 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 579 580 /* 581 Loop over rows of vector, putting results into Jacobian matrix 582 */ 583 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 584 for (l=0; l<coloring->nrows[k]; l++) { 585 row = coloring->rows[k][l]; 586 col = coloring->columnsforrow[k][l]; 587 y[row] *= vscale_array[vscaleforrow[k][l]]; 588 srow = row + start; 589 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 590 } 591 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 592 } 593 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 594 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 595 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 597 ierr = PLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 598 599 ierr = OptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 600 if (flg) { 601 ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 602 } 603 PetscFunctionReturn(0); 604 } 605 606 #undef __FUNC__ 607 #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApplyTS" 608 /*@ 609 MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 610 has been created, computes the Jacobian for a function via finite differences. 611 612 Collective on Mat, MatFDColoring, and Vec 613 614 Input Parameters: 615 + mat - location to store Jacobian 616 . coloring - coloring context created with MatFDColoringCreate() 617 . x1 - location at which Jacobian is to be computed 618 - sctx - optional context required by function (actually a SNES context) 619 620 Options Database Keys: 621 . -mat_fd_coloring_freq <freq> - Sets coloring frequency 622 623 Level: intermediate 624 625 .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 626 627 .keywords: coloring, Jacobian, finite differences 628 @*/ 629 int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 630 { 631 int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 632 int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 633 Scalar dx,mone = -1.0,*y,*xx,*w3_array; 634 Scalar *vscale_array; 635 PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 636 Vec w1,w2,w3; 637 void *fctx = coloring->fctx; 638 PetscTruth flg; 639 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(J,MAT_COOKIE); 642 PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 643 PetscValidHeaderSpecific(x1,VEC_COOKIE); 644 645 ierr = PLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 646 if (!coloring->w1) { 647 ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 648 PLogObjectParent(coloring,coloring->w1); 649 ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 650 PLogObjectParent(coloring,coloring->w2); 651 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 652 PLogObjectParent(coloring,coloring->w3); 653 } 654 w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 655 656 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 657 ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 658 if (flg) { 659 PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 660 } else { 661 ierr = MatZeroEntries(J);CHKERRQ(ierr); 662 } 663 664 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 665 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 666 ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 667 668 /* 669 Compute all the scale factors and share with other processors 670 */ 671 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 672 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 673 for (k=0; k<coloring->ncolors; k++) { 674 /* 675 Loop over each column associated with color adding the 676 perturbation to the vector w3. 677 */ 678 for (l=0; l<coloring->ncolumns[k]; l++) { 679 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 680 dx = xx[col]; 681 if (dx == 0.0) dx = 1.0; 682 #if !defined(PETSC_USE_COMPLEX) 683 if (dx < umin && dx >= 0.0) dx = umin; 684 else if (dx < 0.0 && dx > -umin) dx = -umin; 685 #else 686 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 687 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 688 #endif 689 dx *= epsilon; 690 vscale_array[col] = 1.0/dx; 691 } 692 } 693 vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 694 ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695 ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 696 697 if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 698 else vscaleforrow = coloring->columnsforrow; 699 700 ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 701 /* 702 Loop over each color 703 */ 704 for (k=0; k<coloring->ncolors; k++) { 705 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 706 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 707 /* 708 Loop over each column associated with color adding the 709 perturbation to the vector w3. 710 */ 711 for (l=0; l<coloring->ncolumns[k]; l++) { 712 col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 713 dx = xx[col]; 714 if (dx == 0.0) dx = 1.0; 715 #if !defined(PETSC_USE_COMPLEX) 716 if (dx < umin && dx >= 0.0) dx = umin; 717 else if (dx < 0.0 && dx > -umin) dx = -umin; 718 #else 719 if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 720 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 721 #endif 722 dx *= epsilon; 723 w3_array[col] += dx; 724 } 725 w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 726 727 /* 728 Evaluate function at x1 + dx (here dx is a vector of perturbations) 729 */ 730 ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 731 ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 732 733 /* 734 Loop over rows of vector, putting results into Jacobian matrix 735 */ 736 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 737 for (l=0; l<coloring->nrows[k]; l++) { 738 row = coloring->rows[k][l]; 739 col = coloring->columnsforrow[k][l]; 740 y[row] *= vscale_array[vscaleforrow[k][l]]; 741 srow = row + start; 742 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 743 } 744 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 745 } 746 ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 747 xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 748 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 749 ierr = PLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 750 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 751 PetscFunctionReturn(0); 752 } 753 754 755 756