1840b8ebdSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*6831982aSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.52 1999/10/13 20:37:47 bsmith Exp bsmith $"; 4bbf0e169SBarry Smith #endif 5bbf0e169SBarry Smith 6bbf0e169SBarry Smith /* 7639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 8639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 9bbf0e169SBarry Smith */ 10bbf0e169SBarry Smith 11bbf0e169SBarry Smith #include "petsc.h" 12bbf0e169SBarry Smith #include "src/mat/matimpl.h" /*I "mat.h" I*/ 13bbf0e169SBarry Smith #include "src/vec/vecimpl.h" 14bbf0e169SBarry Smith 155615d1e5SSatish Balay #undef __FUNC__ 16d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw" 17005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer) 18005c665bSBarry Smith { 19005c665bSBarry Smith int ierr,i,j,pause; 20005c665bSBarry Smith PetscTruth isnull; 21005c665bSBarry Smith Draw draw; 22d4bb536fSBarry Smith double xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0; 23005c665bSBarry Smith DrawButton button; 24005c665bSBarry Smith 253a40ed3dSBarry Smith PetscFunctionBegin; 2677ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 273a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 282bdab257SBarry Smith ierr = DrawSynchronizedClear(draw);CHKERRQ(ierr); 29005c665bSBarry Smith 30005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 31005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 32005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 33005c665bSBarry Smith 34005c665bSBarry Smith /* loop over colors */ 35005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 36005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 37005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 38005c665bSBarry Smith x = fd->columnsforrow[i][j]; 395cd90555SBarry Smith ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 40005c665bSBarry Smith } 41005c665bSBarry Smith } 422bdab257SBarry Smith ierr = DrawSynchronizedFlush(draw);CHKERRQ(ierr); 43005c665bSBarry Smith ierr = DrawGetPause(draw,&pause);CHKERRQ(ierr); 443a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 455cd90555SBarry Smith ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr); 465cd90555SBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr); 47005c665bSBarry Smith while (button != BUTTON_RIGHT) { 482bdab257SBarry Smith ierr = DrawSynchronizedClear(draw);CHKERRQ(ierr); 49005c665bSBarry Smith if (button == BUTTON_LEFT) scale = .5; 50005c665bSBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 51005c665bSBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 52005c665bSBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 53005c665bSBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 54005c665bSBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 55005c665bSBarry Smith w *= scale; h *= scale; 56005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 57005c665bSBarry Smith /* loop over colors */ 58005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 59005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 60005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 61005c665bSBarry Smith x = fd->columnsforrow[i][j]; 625cd90555SBarry Smith ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 63005c665bSBarry Smith } 64005c665bSBarry Smith } 655cd90555SBarry Smith ierr = DrawCheckResizedWindow(draw);CHKERRQ(ierr); 665cd90555SBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);CHKERRQ(ierr); 67005c665bSBarry Smith } 68005c665bSBarry Smith 693a40ed3dSBarry Smith PetscFunctionReturn(0); 70005c665bSBarry Smith } 71005c665bSBarry Smith 72005c665bSBarry Smith #undef __FUNC__ 73d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView" 74bbf0e169SBarry Smith /*@C 75639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 76bbf0e169SBarry Smith 77fee21e36SBarry Smith Collective on MatFDColoring unless Viewer is VIEWER_STDOUT_SELF 78fee21e36SBarry Smith 79ef5ee4d1SLois Curfman McInnes Input Parameters: 80ef5ee4d1SLois Curfman McInnes + c - the coloring context 81ef5ee4d1SLois Curfman McInnes - viewer - visualization context 82ef5ee4d1SLois Curfman McInnes 8315091d37SBarry Smith Level: intermediate 8415091d37SBarry Smith 85b4fc646aSLois Curfman McInnes Notes: 86b4fc646aSLois Curfman McInnes The available visualization contexts include 87ef5ee4d1SLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 88ef5ee4d1SLois Curfman McInnes . VIEWER_STDOUT_WORLD - synchronized standard 89ef5ee4d1SLois Curfman McInnes output where only the first processor opens 90ef5ee4d1SLois Curfman McInnes the file. All other processors send their 91ef5ee4d1SLois Curfman McInnes data to the first processor to print. 92c655490fSBarry Smith - VIEWER_DRAW_WORLD - graphical display of nonzero structure 93bbf0e169SBarry Smith 94639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 95005c665bSBarry Smith 96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 97bbf0e169SBarry Smith @*/ 98b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer) 99bbf0e169SBarry Smith { 100639f9d9dSBarry Smith int i,j,format,ierr; 101*6831982aSBarry Smith PetscTruth isdraw,isascii; 102bbf0e169SBarry Smith 1033a40ed3dSBarry Smith PetscFunctionBegin; 104b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 1050f5bd95cSBarry Smith if (!viewer) viewer = VIEWER_STDOUT_SELF; 1060f5bd95cSBarry Smith PetscValidHeaderSpecific(viewer,VIEWER_COOKIE); 107*6831982aSBarry Smith PetscCheckSameComm(c,viewer); 108bbf0e169SBarry Smith 109*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 110*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 1110f5bd95cSBarry Smith if (isdraw) { 112b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 1130f5bd95cSBarry Smith } else if (isascii) { 114d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 115d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr); 116d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Umin=%g\n",c->umin);CHKERRQ(ierr); 117d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors);CHKERRQ(ierr); 118ae09f205SBarry Smith 119ae09f205SBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 120ae09f205SBarry Smith if (format != VIEWER_FORMAT_ASCII_INFO) { 121b4fc646aSLois Curfman McInnes for ( i=0; i<c->ncolors; i++ ) { 122d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Information for color %d\n",i);CHKERRQ(ierr); 123d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr); 124b4fc646aSLois Curfman McInnes for ( j=0; j<c->ncolumns[i]; j++ ) { 125d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]);CHKERRQ(ierr); 126639f9d9dSBarry Smith } 127d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr); 128b4fc646aSLois Curfman McInnes for ( j=0; j<c->nrows[i]; j++ ) { 129d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 130b4fc646aSLois Curfman McInnes } 131bbf0e169SBarry Smith } 132bbf0e169SBarry Smith } 133*6831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 1345cd90555SBarry Smith } else { 1350f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 136bbf0e169SBarry Smith } 1373a40ed3dSBarry Smith PetscFunctionReturn(0); 138639f9d9dSBarry Smith } 139639f9d9dSBarry Smith 1405615d1e5SSatish Balay #undef __FUNC__ 1415615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters" 142639f9d9dSBarry Smith /*@ 143b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 144b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 145639f9d9dSBarry Smith 146ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 147ef5ee4d1SLois Curfman McInnes 148ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 149ef5ee4d1SLois Curfman McInnes .vb 15065f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 151f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 152f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 153ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 154ef5ee4d1SLois Curfman McInnes .ve 155639f9d9dSBarry Smith 156639f9d9dSBarry Smith Input Parameters: 157ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 158639f9d9dSBarry Smith . error_rel - relative error 159f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 160fee21e36SBarry Smith 16115091d37SBarry Smith Level: advanced 16215091d37SBarry Smith 163b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 164b4fc646aSLois Curfman McInnes 165b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 166639f9d9dSBarry Smith @*/ 167639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 168639f9d9dSBarry Smith { 1693a40ed3dSBarry Smith PetscFunctionBegin; 170639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 171639f9d9dSBarry Smith 172639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 173639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1743a40ed3dSBarry Smith PetscFunctionReturn(0); 175639f9d9dSBarry Smith } 176639f9d9dSBarry Smith 1775615d1e5SSatish Balay #undef __FUNC__ 178005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 179005c665bSBarry Smith /*@ 180e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 181e0907662SLois Curfman McInnes matrices. 182005c665bSBarry Smith 183fee21e36SBarry Smith Collective on MatFDColoring 184fee21e36SBarry Smith 185ef5ee4d1SLois Curfman McInnes Input Parameters: 186ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 187ef5ee4d1SLois Curfman McInnes - freq - frequency (default is 1) 188ef5ee4d1SLois Curfman McInnes 18915091d37SBarry Smith Options Database Keys: 19015091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 19115091d37SBarry Smith 19215091d37SBarry Smith Level: advanced 19315091d37SBarry Smith 194e0907662SLois Curfman McInnes Notes: 195e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 196e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 197e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 198e0907662SLois Curfman McInnes <freq> nonlinear iterations. 199e0907662SLois Curfman McInnes 200b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 201ef5ee4d1SLois Curfman McInnes 202ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency() 203005c665bSBarry Smith @*/ 204005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 205005c665bSBarry Smith { 2063a40ed3dSBarry Smith PetscFunctionBegin; 207005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 208005c665bSBarry Smith 209005c665bSBarry Smith matfd->freq = freq; 2103a40ed3dSBarry Smith PetscFunctionReturn(0); 211005c665bSBarry Smith } 212005c665bSBarry Smith 213005c665bSBarry Smith #undef __FUNC__ 214ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency" 215ff0cfa39SBarry Smith /*@ 216ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 217ff0cfa39SBarry Smith matrices. 218ff0cfa39SBarry Smith 219ef5ee4d1SLois Curfman McInnes Not Collective 220ef5ee4d1SLois Curfman McInnes 221ff0cfa39SBarry Smith Input Parameters: 222ff0cfa39SBarry Smith . coloring - the coloring context 223ff0cfa39SBarry Smith 224ff0cfa39SBarry Smith Output Parameters: 225ff0cfa39SBarry Smith . freq - frequency (default is 1) 226ff0cfa39SBarry Smith 22715091d37SBarry Smith Options Database Keys: 22815091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 22915091d37SBarry Smith 23015091d37SBarry Smith Level: advanced 23115091d37SBarry Smith 232ff0cfa39SBarry Smith Notes: 233ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 234ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 235ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 236ff0cfa39SBarry Smith <freq> nonlinear iterations. 237ff0cfa39SBarry Smith 238ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 239ef5ee4d1SLois Curfman McInnes 240ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency() 241ff0cfa39SBarry Smith @*/ 242ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 243ff0cfa39SBarry Smith { 2443a40ed3dSBarry Smith PetscFunctionBegin; 245ff0cfa39SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 246ff0cfa39SBarry Smith 247ff0cfa39SBarry Smith *freq = matfd->freq; 2483a40ed3dSBarry Smith PetscFunctionReturn(0); 249ff0cfa39SBarry Smith } 250ff0cfa39SBarry Smith 251ff0cfa39SBarry Smith #undef __FUNC__ 252005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 253d64ed03dSBarry Smith /*@C 254005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 255005c665bSBarry Smith 256fee21e36SBarry Smith Collective on MatFDColoring 257fee21e36SBarry Smith 258ef5ee4d1SLois Curfman McInnes Input Parameters: 259ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 260ef5ee4d1SLois Curfman McInnes . f - the function 261ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 262ef5ee4d1SLois Curfman McInnes 26315091d37SBarry Smith Level: intermediate 26415091d37SBarry Smith 265f881d145SBarry Smith Notes: 266f881d145SBarry Smith In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 267f881d145SBarry Smith be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 268f881d145SBarry Smith with the TS solvers. 269f881d145SBarry Smith 270b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 271005c665bSBarry Smith @*/ 272840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 273005c665bSBarry Smith { 2743a40ed3dSBarry Smith PetscFunctionBegin; 275005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 276005c665bSBarry Smith 277005c665bSBarry Smith matfd->f = f; 278005c665bSBarry Smith matfd->fctx = fctx; 279005c665bSBarry Smith 2803a40ed3dSBarry Smith PetscFunctionReturn(0); 281005c665bSBarry Smith } 282005c665bSBarry Smith 283005c665bSBarry Smith #undef __FUNC__ 284d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 285639f9d9dSBarry Smith /*@ 286b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 287639f9d9dSBarry Smith the options database. 288639f9d9dSBarry Smith 289fee21e36SBarry Smith Collective on MatFDColoring 290fee21e36SBarry Smith 29165f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 292ef5ee4d1SLois Curfman McInnes .vb 29365f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 294f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 295f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 296ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 297ef5ee4d1SLois Curfman McInnes .ve 298ef5ee4d1SLois Curfman McInnes 299ef5ee4d1SLois Curfman McInnes Input Parameter: 300ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 301ef5ee4d1SLois Curfman McInnes 302b4fc646aSLois Curfman McInnes Options Database Keys: 303ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_error <err> - Sets <err> (square root 304ef5ee4d1SLois Curfman McInnes of relative error in the function) 305f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 306ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 307ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 308ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 309ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 310639f9d9dSBarry Smith 31115091d37SBarry Smith Level: intermediate 31215091d37SBarry Smith 313b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 314639f9d9dSBarry Smith @*/ 315639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 316639f9d9dSBarry Smith { 317005c665bSBarry Smith int ierr,flag,freq = 1; 318639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 3193a40ed3dSBarry Smith 3203a40ed3dSBarry Smith PetscFunctionBegin; 321639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 322639f9d9dSBarry Smith 323639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 324639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 325639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin);CHKERRQ(ierr); 326005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 327005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 328005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag);CHKERRQ(ierr); 329639f9d9dSBarry Smith if (flag) { 330639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd);CHKERRQ(ierr); 331639f9d9dSBarry Smith } 3323a40ed3dSBarry Smith PetscFunctionReturn(0); 333639f9d9dSBarry Smith } 334639f9d9dSBarry Smith 3355615d1e5SSatish Balay #undef __FUNC__ 336d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 337639f9d9dSBarry Smith /*@ 338639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 339639f9d9dSBarry Smith using coloring. 340639f9d9dSBarry Smith 341ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 342ef5ee4d1SLois Curfman McInnes 343639f9d9dSBarry Smith Input Parameter: 344639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 345639f9d9dSBarry Smith 34615091d37SBarry Smith Level: intermediate 34715091d37SBarry Smith 348639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 349639f9d9dSBarry Smith @*/ 350639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 351639f9d9dSBarry Smith { 352d132466eSBarry Smith int ierr; 353d132466eSBarry Smith 3543a40ed3dSBarry Smith PetscFunctionBegin; 355639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 356639f9d9dSBarry Smith 357d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);CHKERRQ(ierr); 358d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);CHKERRQ(ierr); 359d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);CHKERRQ(ierr); 360d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");CHKERRQ(ierr); 361d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");CHKERRQ(ierr); 362d132466eSBarry Smith ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");CHKERRQ(ierr); 3633a40ed3dSBarry Smith PetscFunctionReturn(0); 364005c665bSBarry Smith } 365005c665bSBarry Smith 366433994e6SBarry Smith #undef __FUNC__ 367433994e6SBarry Smith #define __FUNC__ "MatFDColoringView_Private" 368005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 369005c665bSBarry Smith { 370005c665bSBarry Smith int ierr,flg; 371005c665bSBarry Smith 3723a40ed3dSBarry Smith PetscFunctionBegin; 373005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 374005c665bSBarry Smith if (flg) { 375f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 376005c665bSBarry Smith } 377ae09f205SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 378ae09f205SBarry Smith if (flg) { 379f8590f6eSBarry Smith ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 380f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 381f8590f6eSBarry Smith ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 382ae09f205SBarry Smith } 383005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 384005c665bSBarry Smith if (flg) { 385c655490fSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 386c655490fSBarry Smith ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 387005c665bSBarry Smith } 3883a40ed3dSBarry Smith PetscFunctionReturn(0); 389bbf0e169SBarry Smith } 390bbf0e169SBarry Smith 3915615d1e5SSatish Balay #undef __FUNC__ 3925615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 393bbf0e169SBarry Smith /*@C 394639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 395639f9d9dSBarry Smith computation of Jacobians. 396bbf0e169SBarry Smith 397ef5ee4d1SLois Curfman McInnes Collective on Mat 398ef5ee4d1SLois Curfman McInnes 399639f9d9dSBarry Smith Input Parameters: 400ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 401ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 402bbf0e169SBarry Smith 403bbf0e169SBarry Smith Output Parameter: 404639f9d9dSBarry Smith . color - the new coloring context 405bbf0e169SBarry Smith 406b4fc646aSLois Curfman McInnes Options Database Keys: 407ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_view - Activates basic viewing or coloring 408ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 409ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 410639f9d9dSBarry Smith 41115091d37SBarry Smith Level: intermediate 41215091d37SBarry Smith 413*6831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 414*6831982aSBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 415bbf0e169SBarry Smith @*/ 416639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 417bbf0e169SBarry Smith { 418639f9d9dSBarry Smith MatFDColoring c; 419639f9d9dSBarry Smith MPI_Comm comm; 420639f9d9dSBarry Smith int ierr,M,N; 421639f9d9dSBarry Smith 4223a40ed3dSBarry Smith PetscFunctionBegin; 423639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 424e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 425639f9d9dSBarry Smith 426f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 4273f1db9ecSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm, 4283f1db9ecSBarry Smith MatFDColoringDestroy,MatFDColoringView); 429639f9d9dSBarry Smith PLogObjectCreate(c); 430639f9d9dSBarry Smith 431f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 432f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 433639f9d9dSBarry Smith } else { 434e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 435639f9d9dSBarry Smith } 436639f9d9dSBarry Smith 437639f9d9dSBarry Smith c->error_rel = 1.e-8; 438ae09f205SBarry Smith c->umin = 1.e-6; 439005c665bSBarry Smith c->freq = 1; 440005c665bSBarry Smith 441005c665bSBarry Smith ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 442639f9d9dSBarry Smith 443639f9d9dSBarry Smith *color = c; 444639f9d9dSBarry Smith 4453a40ed3dSBarry Smith PetscFunctionReturn(0); 446bbf0e169SBarry Smith } 447bbf0e169SBarry Smith 4485615d1e5SSatish Balay #undef __FUNC__ 449d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 450bbf0e169SBarry Smith /*@C 451639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 452639f9d9dSBarry Smith via MatFDColoringCreate(). 453bbf0e169SBarry Smith 454ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 455ef5ee4d1SLois Curfman McInnes 456b4fc646aSLois Curfman McInnes Input Parameter: 457639f9d9dSBarry Smith . c - coloring context 458bbf0e169SBarry Smith 45915091d37SBarry Smith Level: intermediate 46015091d37SBarry Smith 461639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 462bbf0e169SBarry Smith @*/ 463639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 464bbf0e169SBarry Smith { 465263760aaSBarry Smith int i,ierr; 466bbf0e169SBarry Smith 4673a40ed3dSBarry Smith PetscFunctionBegin; 4683a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 469d4bb536fSBarry Smith 470639f9d9dSBarry Smith 471639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 472606d414cSSatish Balay if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 473606d414cSSatish Balay if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 474606d414cSSatish Balay if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 475bbf0e169SBarry Smith } 476606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 477606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 478606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 479606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 480606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 481606d414cSSatish Balay ierr = PetscFree(c->scale);CHKERRQ(ierr); 482005c665bSBarry Smith if (c->w1) { 483005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 484005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 485005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 486005c665bSBarry Smith } 487639f9d9dSBarry Smith PLogObjectDestroy(c); 488639f9d9dSBarry Smith PetscHeaderDestroy(c); 4893a40ed3dSBarry Smith PetscFunctionReturn(0); 490bbf0e169SBarry Smith } 49143a90d84SBarry Smith 492005c665bSBarry Smith 4935615d1e5SSatish Balay #undef __FUNC__ 4945615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 49543a90d84SBarry Smith /*@ 496e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 497e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 49843a90d84SBarry Smith 499fee21e36SBarry Smith Collective on MatFDColoring 500fee21e36SBarry Smith 501ef5ee4d1SLois Curfman McInnes Input Parameters: 502ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 503ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 504ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 505ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 506ef5ee4d1SLois Curfman McInnes 5078bba8e72SBarry Smith Options Database Keys: 508ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 5098bba8e72SBarry Smith 51015091d37SBarry Smith Level: intermediate 51115091d37SBarry Smith 51243a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 51343a90d84SBarry Smith 51443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 51543a90d84SBarry Smith @*/ 516005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 51743a90d84SBarry Smith { 518e0907662SLois Curfman McInnes int k,fg,ierr,N,start,end,l,row,col,srow; 51943a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 52043a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 52143a90d84SBarry Smith MPI_Comm comm = coloring->comm; 522005c665bSBarry Smith Vec w1,w2,w3; 523840b8ebdSBarry Smith int (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f; 524005c665bSBarry Smith void *fctx = coloring->fctx; 525005c665bSBarry Smith 5263a40ed3dSBarry Smith PetscFunctionBegin; 527e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 528e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 529e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 530e0907662SLois Curfman McInnes 531005c665bSBarry Smith 532005c665bSBarry Smith if (!coloring->w1) { 533005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 534005c665bSBarry Smith PLogObjectParent(coloring,coloring->w1); 535005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 536005c665bSBarry Smith PLogObjectParent(coloring,coloring->w2); 537005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 538005c665bSBarry Smith PLogObjectParent(coloring,coloring->w3); 539005c665bSBarry Smith } 540005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 54143a90d84SBarry Smith 542e0907662SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg);CHKERRQ(ierr); 543e0907662SLois Curfman McInnes if (fg) { 544e0907662SLois Curfman McInnes PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 545e0907662SLois Curfman McInnes } else { 54643a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 547e0907662SLois Curfman McInnes } 54843a90d84SBarry Smith 54943a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 55043a90d84SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 55143a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 55243a90d84SBarry Smith 553549d3d68SSatish Balay ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr); 55443a90d84SBarry Smith /* 55543a90d84SBarry Smith Loop over each color 55643a90d84SBarry Smith */ 55743a90d84SBarry Smith 5583b28642cSBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 55943a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 56043a90d84SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 56143a90d84SBarry Smith /* 56243a90d84SBarry Smith Loop over each column associated with color adding the 56343a90d84SBarry Smith perturbation to the vector w3. 56443a90d84SBarry Smith */ 56543a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 56643a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 56743a90d84SBarry Smith dx = xx[col-start]; 568ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 569aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 570ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 571ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 57243a90d84SBarry Smith #else 573e20fef11SSatish Balay if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0) dx = umin; 574e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 57543a90d84SBarry Smith #endif 57643a90d84SBarry Smith dx *= epsilon; 57743a90d84SBarry Smith wscale[col] = 1.0/dx; 5783b28642cSBarry Smith ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr); 57943a90d84SBarry Smith } 5803b28642cSBarry Smith 58143a90d84SBarry Smith /* 582e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 58343a90d84SBarry Smith */ 58443a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 58543a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 58643a90d84SBarry Smith /* Communicate scale to all processors */ 587*6831982aSBarry Smith ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 58843a90d84SBarry Smith /* 589e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 59043a90d84SBarry Smith */ 5913b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 59243a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 59343a90d84SBarry Smith row = coloring->rows[k][l]; 59443a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 59543a90d84SBarry Smith y[row] *= scale[col]; 59643a90d84SBarry Smith srow = row + start; 59743a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 59843a90d84SBarry Smith } 5993b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 60043a90d84SBarry Smith } 6013b28642cSBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 60243a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60343a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6043a40ed3dSBarry Smith PetscFunctionReturn(0); 60543a90d84SBarry Smith } 606840b8ebdSBarry Smith 607840b8ebdSBarry Smith #undef __FUNC__ 608840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS" 609840b8ebdSBarry Smith /*@ 610840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 611840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 612840b8ebdSBarry Smith 613fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 614fee21e36SBarry Smith 615ef5ee4d1SLois Curfman McInnes Input Parameters: 6163b28642cSBarry Smith + mat - location to store Jacobian 617ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 618ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 619ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 620ef5ee4d1SLois Curfman McInnes 621840b8ebdSBarry Smith Options Database Keys: 622ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 623840b8ebdSBarry Smith 62415091d37SBarry Smith Level: intermediate 62515091d37SBarry Smith 626840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 627840b8ebdSBarry Smith 628840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 629840b8ebdSBarry Smith @*/ 630840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx) 631840b8ebdSBarry Smith { 632840b8ebdSBarry Smith int k,fg,ierr,N,start,end,l,row,col,srow; 633840b8ebdSBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 634840b8ebdSBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 635840b8ebdSBarry Smith MPI_Comm comm = coloring->comm; 636840b8ebdSBarry Smith Vec w1,w2,w3; 637840b8ebdSBarry Smith int (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f; 638840b8ebdSBarry Smith void *fctx = coloring->fctx; 639840b8ebdSBarry Smith 6403a40ed3dSBarry Smith PetscFunctionBegin; 641840b8ebdSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE); 642840b8ebdSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 643840b8ebdSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE); 644840b8ebdSBarry Smith 645840b8ebdSBarry Smith if (!coloring->w1) { 646840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 647840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w1); 648840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 649840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w2); 650840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 651840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w3); 652840b8ebdSBarry Smith } 653840b8ebdSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 654840b8ebdSBarry Smith 655840b8ebdSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg);CHKERRQ(ierr); 656840b8ebdSBarry Smith if (fg) { 657840b8ebdSBarry Smith PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n"); 658840b8ebdSBarry Smith } else { 659840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 660840b8ebdSBarry Smith } 661840b8ebdSBarry Smith 662840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 663840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 664840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 665840b8ebdSBarry Smith 666549d3d68SSatish Balay ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr); 667840b8ebdSBarry Smith /* 668840b8ebdSBarry Smith Loop over each color 669840b8ebdSBarry Smith */ 670840b8ebdSBarry Smith 6713b28642cSBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 672840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 673840b8ebdSBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 674840b8ebdSBarry Smith /* 675840b8ebdSBarry Smith Loop over each column associated with color adding the 676840b8ebdSBarry Smith perturbation to the vector w3. 677840b8ebdSBarry Smith */ 678840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 679840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 680840b8ebdSBarry Smith dx = xx[col-start]; 681840b8ebdSBarry Smith if (dx == 0.0) dx = 1.0; 682aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 683840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 684840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 685840b8ebdSBarry Smith #else 686e20fef11SSatish Balay if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0) dx = umin; 687e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 688840b8ebdSBarry Smith #endif 689840b8ebdSBarry Smith dx *= epsilon; 690840b8ebdSBarry Smith wscale[col] = 1.0/dx; 6913b28642cSBarry Smith ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr); 692840b8ebdSBarry Smith } 693840b8ebdSBarry Smith /* 694840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 695840b8ebdSBarry Smith */ 696840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 697840b8ebdSBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 698840b8ebdSBarry Smith /* Communicate scale to all processors */ 699*6831982aSBarry Smith ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 700840b8ebdSBarry Smith /* 701840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 702840b8ebdSBarry Smith */ 7033b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 704840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 705840b8ebdSBarry Smith row = coloring->rows[k][l]; 706840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 707840b8ebdSBarry Smith y[row] *= scale[col]; 708840b8ebdSBarry Smith srow = row + start; 709840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 710840b8ebdSBarry Smith } 7113b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 712840b8ebdSBarry Smith } 7133b28642cSBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 714840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 715840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7163a40ed3dSBarry Smith PetscFunctionReturn(0); 717840b8ebdSBarry Smith } 7183b28642cSBarry Smith 7193b28642cSBarry Smith 7203b28642cSBarry Smith 721