1840b8ebdSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*f881d145SBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.41 1999/01/13 21:40:15 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 83b4fc646aSLois Curfman McInnes Notes: 84b4fc646aSLois Curfman McInnes The available visualization contexts include 85ef5ee4d1SLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 86ef5ee4d1SLois Curfman McInnes . VIEWER_STDOUT_WORLD - synchronized standard 87ef5ee4d1SLois Curfman McInnes output where only the first processor opens 88ef5ee4d1SLois Curfman McInnes the file. All other processors send their 89ef5ee4d1SLois Curfman McInnes data to the first processor to print. 90c655490fSBarry Smith - VIEWER_DRAW_WORLD - graphical display of nonzero structure 91bbf0e169SBarry Smith 92639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 93005c665bSBarry Smith 94b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 95bbf0e169SBarry Smith @*/ 96b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer) 97bbf0e169SBarry Smith { 98005c665bSBarry Smith ViewerType vtype; 99639f9d9dSBarry Smith int i,j,format,ierr; 100bbf0e169SBarry Smith 1013a40ed3dSBarry Smith PetscFunctionBegin; 102b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 103b4fc646aSLois Curfman McInnes if (viewer) {PetscValidHeader(viewer);} 104b4fc646aSLois Curfman McInnes else {viewer = VIEWER_STDOUT_SELF;} 105bbf0e169SBarry Smith 106005c665bSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 1073f1db9ecSBarry Smith if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 108b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr); 1093a40ed3dSBarry Smith PetscFunctionReturn(0); 1103f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 1110ef38995SBarry Smith ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n"); 1120ef38995SBarry Smith ViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel); 1130ef38995SBarry Smith ViewerASCIIPrintf(viewer," Umin=%g\n",c->umin); 1140ef38995SBarry Smith ViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors); 115ae09f205SBarry Smith 116ae09f205SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 117ae09f205SBarry Smith if (format != VIEWER_FORMAT_ASCII_INFO) { 118b4fc646aSLois Curfman McInnes for ( i=0; i<c->ncolors; i++ ) { 1190ef38995SBarry Smith ViewerASCIIPrintf(viewer," Information for color %d\n",i); 1200ef38995SBarry Smith ViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]); 121b4fc646aSLois Curfman McInnes for ( j=0; j<c->ncolumns[i]; j++ ) { 1220ef38995SBarry Smith ViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]); 123639f9d9dSBarry Smith } 1240ef38995SBarry Smith ViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]); 125b4fc646aSLois Curfman McInnes for ( j=0; j<c->nrows[i]; j++ ) { 1260ef38995SBarry Smith ViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]); 127b4fc646aSLois Curfman McInnes } 128bbf0e169SBarry Smith } 129bbf0e169SBarry Smith } 1305cd90555SBarry Smith } else { 1315cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 132bbf0e169SBarry Smith } 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 134639f9d9dSBarry Smith } 135639f9d9dSBarry Smith 1365615d1e5SSatish Balay #undef __FUNC__ 1375615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters" 138639f9d9dSBarry Smith /*@ 139b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 140b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 141639f9d9dSBarry Smith 142ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 143ef5ee4d1SLois Curfman McInnes 144ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 145ef5ee4d1SLois Curfman McInnes .vb 146ef5ee4d1SLois Curfman McInnes J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 147f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 148f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 149ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 150ef5ee4d1SLois Curfman McInnes .ve 151639f9d9dSBarry Smith 152639f9d9dSBarry Smith Input Parameters: 153ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 154639f9d9dSBarry Smith . error_rel - relative error 155f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 156fee21e36SBarry Smith 157b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 158b4fc646aSLois Curfman McInnes 159b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 160639f9d9dSBarry Smith @*/ 161639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 162639f9d9dSBarry Smith { 1633a40ed3dSBarry Smith PetscFunctionBegin; 164639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 165639f9d9dSBarry Smith 166639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 167639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1683a40ed3dSBarry Smith PetscFunctionReturn(0); 169639f9d9dSBarry Smith } 170639f9d9dSBarry Smith 1715615d1e5SSatish Balay #undef __FUNC__ 172005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 173005c665bSBarry Smith /*@ 174e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 175e0907662SLois Curfman McInnes matrices. 176005c665bSBarry Smith 177fee21e36SBarry Smith Collective on MatFDColoring 178fee21e36SBarry Smith 179ef5ee4d1SLois Curfman McInnes Input Parameters: 180ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 181ef5ee4d1SLois Curfman McInnes - freq - frequency (default is 1) 182ef5ee4d1SLois Curfman McInnes 183e0907662SLois Curfman McInnes Notes: 184e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 185e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 186e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 187e0907662SLois Curfman McInnes <freq> nonlinear iterations. 188e0907662SLois Curfman McInnes 189e0907662SLois Curfman McInnes Options Database Keys: 190ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 191005c665bSBarry Smith 192b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 193ef5ee4d1SLois Curfman McInnes 194ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency() 195005c665bSBarry Smith @*/ 196005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 197005c665bSBarry Smith { 1983a40ed3dSBarry Smith PetscFunctionBegin; 199005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 200005c665bSBarry Smith 201005c665bSBarry Smith matfd->freq = freq; 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 203005c665bSBarry Smith } 204005c665bSBarry Smith 205005c665bSBarry Smith #undef __FUNC__ 206ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency" 207ff0cfa39SBarry Smith /*@ 208ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 209ff0cfa39SBarry Smith matrices. 210ff0cfa39SBarry Smith 211ef5ee4d1SLois Curfman McInnes Not Collective 212ef5ee4d1SLois Curfman McInnes 213ff0cfa39SBarry Smith Input Parameters: 214ff0cfa39SBarry Smith . coloring - the coloring context 215ff0cfa39SBarry Smith 216ff0cfa39SBarry Smith Output Parameters: 217ff0cfa39SBarry Smith . freq - frequency (default is 1) 218ff0cfa39SBarry Smith 219ff0cfa39SBarry Smith Notes: 220ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 221ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 222ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 223ff0cfa39SBarry Smith <freq> nonlinear iterations. 224ff0cfa39SBarry Smith 225ff0cfa39SBarry Smith Options Database Keys: 226ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 227ff0cfa39SBarry Smith 228ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 229ef5ee4d1SLois Curfman McInnes 230ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency() 231ff0cfa39SBarry Smith @*/ 232ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 233ff0cfa39SBarry Smith { 2343a40ed3dSBarry Smith PetscFunctionBegin; 235ff0cfa39SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 236ff0cfa39SBarry Smith 237ff0cfa39SBarry Smith *freq = matfd->freq; 2383a40ed3dSBarry Smith PetscFunctionReturn(0); 239ff0cfa39SBarry Smith } 240ff0cfa39SBarry Smith 241ff0cfa39SBarry Smith #undef __FUNC__ 242005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 243d64ed03dSBarry Smith /*@C 244005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 245005c665bSBarry Smith 246fee21e36SBarry Smith Collective on MatFDColoring 247fee21e36SBarry Smith 248ef5ee4d1SLois Curfman McInnes Input Parameters: 249ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 250ef5ee4d1SLois Curfman McInnes . f - the function 251ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 252ef5ee4d1SLois Curfman McInnes 253*f881d145SBarry Smith Notes: 254*f881d145SBarry Smith In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 255*f881d145SBarry Smith be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 256*f881d145SBarry Smith with the TS solvers. 257*f881d145SBarry Smith 258b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 259005c665bSBarry Smith @*/ 260840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 261005c665bSBarry Smith { 2623a40ed3dSBarry Smith PetscFunctionBegin; 263005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 264005c665bSBarry Smith 265005c665bSBarry Smith matfd->f = f; 266005c665bSBarry Smith matfd->fctx = fctx; 267005c665bSBarry Smith 2683a40ed3dSBarry Smith PetscFunctionReturn(0); 269005c665bSBarry Smith } 270005c665bSBarry Smith 271005c665bSBarry Smith #undef __FUNC__ 272d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 273639f9d9dSBarry Smith /*@ 274b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 275639f9d9dSBarry Smith the options database. 276639f9d9dSBarry Smith 277fee21e36SBarry Smith Collective on MatFDColoring 278fee21e36SBarry Smith 279ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 280ef5ee4d1SLois Curfman McInnes .vb 281ef5ee4d1SLois Curfman McInnes J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 282f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 283f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 284ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 285ef5ee4d1SLois Curfman McInnes .ve 286ef5ee4d1SLois Curfman McInnes 287ef5ee4d1SLois Curfman McInnes Input Parameter: 288ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 289ef5ee4d1SLois Curfman McInnes 290b4fc646aSLois Curfman McInnes Options Database Keys: 291ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_error <err> - Sets <err> (square root 292ef5ee4d1SLois Curfman McInnes of relative error in the function) 293f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 294ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 295ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 296ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 297ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 298639f9d9dSBarry Smith 299b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 300639f9d9dSBarry Smith @*/ 301639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 302639f9d9dSBarry Smith { 303005c665bSBarry Smith int ierr,flag,freq = 1; 304639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 3053a40ed3dSBarry Smith 3063a40ed3dSBarry Smith PetscFunctionBegin; 307639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 308639f9d9dSBarry Smith 309639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 310639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 311639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 312005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 313005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 314005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 315639f9d9dSBarry Smith if (flag) { 316639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 317639f9d9dSBarry Smith } 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 319639f9d9dSBarry Smith } 320639f9d9dSBarry Smith 3215615d1e5SSatish Balay #undef __FUNC__ 322d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 323639f9d9dSBarry Smith /*@ 324639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 325639f9d9dSBarry Smith using coloring. 326639f9d9dSBarry Smith 327ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 328ef5ee4d1SLois Curfman McInnes 329639f9d9dSBarry Smith Input Parameter: 330639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 331639f9d9dSBarry Smith 332639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 333639f9d9dSBarry Smith @*/ 334639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 335639f9d9dSBarry Smith { 3363a40ed3dSBarry Smith PetscFunctionBegin; 337639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 338639f9d9dSBarry Smith 33976be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel); 34076be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin); 34176be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq); 34276be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n"); 34376be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n"); 34476be9ce4SBarry Smith (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n"); 3453a40ed3dSBarry Smith PetscFunctionReturn(0); 346005c665bSBarry Smith } 347005c665bSBarry Smith 348005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 349005c665bSBarry Smith { 350005c665bSBarry Smith int ierr,flg; 351005c665bSBarry Smith 3523a40ed3dSBarry Smith PetscFunctionBegin; 353005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 354005c665bSBarry Smith if (flg) { 355f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 356005c665bSBarry Smith } 357ae09f205SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 358ae09f205SBarry Smith if (flg) { 359f8590f6eSBarry Smith ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 360f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 361f8590f6eSBarry Smith ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 362ae09f205SBarry Smith } 363005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 364005c665bSBarry Smith if (flg) { 365c655490fSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm)); CHKERRQ(ierr); 366c655490fSBarry Smith ierr = ViewerFlush(VIEWER_DRAW_(fd->comm)); CHKERRQ(ierr); 367005c665bSBarry Smith } 3683a40ed3dSBarry Smith PetscFunctionReturn(0); 369bbf0e169SBarry Smith } 370bbf0e169SBarry Smith 3715615d1e5SSatish Balay #undef __FUNC__ 3725615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 373bbf0e169SBarry Smith /*@C 374639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 375639f9d9dSBarry Smith computation of Jacobians. 376bbf0e169SBarry Smith 377ef5ee4d1SLois Curfman McInnes Collective on Mat 378ef5ee4d1SLois Curfman McInnes 379639f9d9dSBarry Smith Input Parameters: 380ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 381ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 382bbf0e169SBarry Smith 383bbf0e169SBarry Smith Output Parameter: 384639f9d9dSBarry Smith . color - the new coloring context 385bbf0e169SBarry Smith 386b4fc646aSLois Curfman McInnes Options Database Keys: 387ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_view - Activates basic viewing or coloring 388ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 389ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 390639f9d9dSBarry Smith 391639f9d9dSBarry Smith .seealso: MatFDColoringDestroy() 392bbf0e169SBarry Smith @*/ 393639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 394bbf0e169SBarry Smith { 395639f9d9dSBarry Smith MatFDColoring c; 396639f9d9dSBarry Smith MPI_Comm comm; 397639f9d9dSBarry Smith int ierr,M,N; 398639f9d9dSBarry Smith 3993a40ed3dSBarry Smith PetscFunctionBegin; 400639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 401e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 402639f9d9dSBarry Smith 403*f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 4043f1db9ecSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm, 4053f1db9ecSBarry Smith MatFDColoringDestroy,MatFDColoringView); 406639f9d9dSBarry Smith PLogObjectCreate(c); 407639f9d9dSBarry Smith 408f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 409f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 410639f9d9dSBarry Smith } else { 411e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 412639f9d9dSBarry Smith } 413639f9d9dSBarry Smith 414639f9d9dSBarry Smith c->error_rel = 1.e-8; 415ae09f205SBarry Smith c->umin = 1.e-6; 416005c665bSBarry Smith c->freq = 1; 417005c665bSBarry Smith 418005c665bSBarry Smith ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 419639f9d9dSBarry Smith 420639f9d9dSBarry Smith *color = c; 421639f9d9dSBarry Smith 4223a40ed3dSBarry Smith PetscFunctionReturn(0); 423bbf0e169SBarry Smith } 424bbf0e169SBarry Smith 4255615d1e5SSatish Balay #undef __FUNC__ 426d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 427bbf0e169SBarry Smith /*@C 428639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 429639f9d9dSBarry Smith via MatFDColoringCreate(). 430bbf0e169SBarry Smith 431ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 432ef5ee4d1SLois Curfman McInnes 433b4fc646aSLois Curfman McInnes Input Parameter: 434639f9d9dSBarry Smith . c - coloring context 435bbf0e169SBarry Smith 436639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 437bbf0e169SBarry Smith @*/ 438639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 439bbf0e169SBarry Smith { 440263760aaSBarry Smith int i,ierr; 441bbf0e169SBarry Smith 4423a40ed3dSBarry Smith PetscFunctionBegin; 4433a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 444d4bb536fSBarry Smith 445639f9d9dSBarry Smith 446639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 447639f9d9dSBarry Smith if (c->columns[i]) PetscFree(c->columns[i]); 448639f9d9dSBarry Smith if (c->rows[i]) PetscFree(c->rows[i]); 449639f9d9dSBarry Smith if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 450bbf0e169SBarry Smith } 451639f9d9dSBarry Smith PetscFree(c->ncolumns); 452639f9d9dSBarry Smith PetscFree(c->columns); 453639f9d9dSBarry Smith PetscFree(c->nrows); 454639f9d9dSBarry Smith PetscFree(c->rows); 455639f9d9dSBarry Smith PetscFree(c->columnsforrow); 456639f9d9dSBarry Smith PetscFree(c->scale); 457005c665bSBarry Smith if (c->w1) { 458005c665bSBarry Smith ierr = VecDestroy(c->w1); CHKERRQ(ierr); 459005c665bSBarry Smith ierr = VecDestroy(c->w2); CHKERRQ(ierr); 460005c665bSBarry Smith ierr = VecDestroy(c->w3); CHKERRQ(ierr); 461005c665bSBarry Smith } 462639f9d9dSBarry Smith PLogObjectDestroy(c); 463639f9d9dSBarry Smith PetscHeaderDestroy(c); 4643a40ed3dSBarry Smith PetscFunctionReturn(0); 465bbf0e169SBarry Smith } 46643a90d84SBarry Smith 467005c665bSBarry Smith #include "snes.h" 468005c665bSBarry Smith 4695615d1e5SSatish Balay #undef __FUNC__ 4705615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 47143a90d84SBarry Smith /*@ 472e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 473e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 47443a90d84SBarry Smith 475fee21e36SBarry Smith Collective on MatFDColoring 476fee21e36SBarry Smith 477ef5ee4d1SLois Curfman McInnes Input Parameters: 478ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 479ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 480ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 481ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 482ef5ee4d1SLois Curfman McInnes 4838bba8e72SBarry Smith Options Database Keys: 484ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 4858bba8e72SBarry Smith 48643a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 48743a90d84SBarry Smith 48843a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 48943a90d84SBarry Smith @*/ 490005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 49143a90d84SBarry Smith { 492e0907662SLois Curfman McInnes int k,fg,ierr,N,start,end,l,row,col,srow; 49343a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 49443a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 49543a90d84SBarry Smith MPI_Comm comm = coloring->comm; 496005c665bSBarry Smith Vec w1,w2,w3; 497840b8ebdSBarry Smith int (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f; 498005c665bSBarry Smith void *fctx = coloring->fctx; 499005c665bSBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 501e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 502e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 503e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 504e0907662SLois Curfman McInnes 505005c665bSBarry Smith 506005c665bSBarry Smith if (!coloring->w1) { 507005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 508005c665bSBarry Smith PLogObjectParent(coloring,coloring->w1); 509005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 510005c665bSBarry Smith PLogObjectParent(coloring,coloring->w2); 511005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 512005c665bSBarry Smith PLogObjectParent(coloring,coloring->w3); 513005c665bSBarry Smith } 514005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 51543a90d84SBarry Smith 516e0907662SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 517e0907662SLois Curfman McInnes if (fg) { 518e0907662SLois Curfman McInnes PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 519e0907662SLois Curfman McInnes } else { 52043a90d84SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 521e0907662SLois Curfman McInnes } 52243a90d84SBarry Smith 52343a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 52443a90d84SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 52543a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 52643a90d84SBarry Smith 52743a90d84SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 52843a90d84SBarry Smith /* 52943a90d84SBarry Smith Loop over each color 53043a90d84SBarry Smith */ 53143a90d84SBarry Smith 5323b28642cSBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 53343a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 53443a90d84SBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 53543a90d84SBarry Smith /* 53643a90d84SBarry Smith Loop over each column associated with color adding the 53743a90d84SBarry Smith perturbation to the vector w3. 53843a90d84SBarry Smith */ 53943a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 54043a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 54143a90d84SBarry Smith dx = xx[col-start]; 542ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 5433a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 544ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 545ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 54643a90d84SBarry Smith #else 547e20fef11SSatish Balay if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0) dx = umin; 548e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 54943a90d84SBarry Smith #endif 55043a90d84SBarry Smith dx *= epsilon; 55143a90d84SBarry Smith wscale[col] = 1.0/dx; 5523b28642cSBarry Smith ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr); 55343a90d84SBarry Smith } 5543b28642cSBarry Smith 55543a90d84SBarry Smith /* 556e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 55743a90d84SBarry Smith */ 55843a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 55943a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 56043a90d84SBarry Smith /* Communicate scale to all processors */ 5613a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 562ca161407SBarry Smith ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 56343a90d84SBarry Smith #else 564ca161407SBarry Smith ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 56543a90d84SBarry Smith #endif 56643a90d84SBarry Smith /* 567e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 56843a90d84SBarry Smith */ 5693b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 57043a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 57143a90d84SBarry Smith row = coloring->rows[k][l]; 57243a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 57343a90d84SBarry Smith y[row] *= scale[col]; 57443a90d84SBarry Smith srow = row + start; 57543a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 57643a90d84SBarry Smith } 5773b28642cSBarry Smith ierr = VecRestoreArray(w2,&y); CHKERRQ(ierr); 57843a90d84SBarry Smith } 5793b28642cSBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 58043a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 58143a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5823a40ed3dSBarry Smith PetscFunctionReturn(0); 58343a90d84SBarry Smith } 584840b8ebdSBarry Smith 585840b8ebdSBarry Smith #include "ts.h" 586840b8ebdSBarry Smith 587840b8ebdSBarry Smith #undef __FUNC__ 588840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS" 589840b8ebdSBarry Smith /*@ 590840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 591840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 592840b8ebdSBarry Smith 593fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 594fee21e36SBarry Smith 595ef5ee4d1SLois Curfman McInnes Input Parameters: 5963b28642cSBarry Smith + mat - location to store Jacobian 597ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 598ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 599ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 600ef5ee4d1SLois Curfman McInnes 601840b8ebdSBarry Smith Options Database Keys: 602ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 603840b8ebdSBarry Smith 604840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 605840b8ebdSBarry Smith 606840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 607840b8ebdSBarry Smith @*/ 608840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx) 609840b8ebdSBarry Smith { 610840b8ebdSBarry Smith int k,fg,ierr,N,start,end,l,row,col,srow; 611840b8ebdSBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 612840b8ebdSBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 613840b8ebdSBarry Smith MPI_Comm comm = coloring->comm; 614840b8ebdSBarry Smith Vec w1,w2,w3; 615840b8ebdSBarry Smith int (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f; 616840b8ebdSBarry Smith void *fctx = coloring->fctx; 617840b8ebdSBarry Smith 6183a40ed3dSBarry Smith PetscFunctionBegin; 619840b8ebdSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE); 620840b8ebdSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 621840b8ebdSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE); 622840b8ebdSBarry Smith 623840b8ebdSBarry Smith if (!coloring->w1) { 624840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 625840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w1); 626840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 627840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w2); 628840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 629840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w3); 630840b8ebdSBarry Smith } 631840b8ebdSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 632840b8ebdSBarry Smith 633840b8ebdSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 634840b8ebdSBarry Smith if (fg) { 635840b8ebdSBarry Smith PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n"); 636840b8ebdSBarry Smith } else { 637840b8ebdSBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 638840b8ebdSBarry Smith } 639840b8ebdSBarry Smith 640840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 641840b8ebdSBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 642840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr); 643840b8ebdSBarry Smith 644840b8ebdSBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 645840b8ebdSBarry Smith /* 646840b8ebdSBarry Smith Loop over each color 647840b8ebdSBarry Smith */ 648840b8ebdSBarry Smith 6493b28642cSBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 650840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 651840b8ebdSBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 652840b8ebdSBarry Smith /* 653840b8ebdSBarry Smith Loop over each column associated with color adding the 654840b8ebdSBarry Smith perturbation to the vector w3. 655840b8ebdSBarry Smith */ 656840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 657840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 658840b8ebdSBarry Smith dx = xx[col-start]; 659840b8ebdSBarry Smith if (dx == 0.0) dx = 1.0; 6603a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 661840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 662840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 663840b8ebdSBarry Smith #else 664e20fef11SSatish Balay if (PetscAbsScalar(dx) < umin && PetscReal(dx) >= 0.0) dx = umin; 665e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 666840b8ebdSBarry Smith #endif 667840b8ebdSBarry Smith dx *= epsilon; 668840b8ebdSBarry Smith wscale[col] = 1.0/dx; 6693b28642cSBarry Smith ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES); CHKERRQ(ierr); 670840b8ebdSBarry Smith } 671840b8ebdSBarry Smith /* 672840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 673840b8ebdSBarry Smith */ 674840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr); 675840b8ebdSBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 676840b8ebdSBarry Smith /* Communicate scale to all processors */ 6773a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 678ca161407SBarry Smith ierr = MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 679840b8ebdSBarry Smith #else 680ca161407SBarry Smith ierr = MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 681840b8ebdSBarry Smith #endif 682840b8ebdSBarry Smith /* 683840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 684840b8ebdSBarry Smith */ 6853b28642cSBarry Smith ierr = VecGetArray(w2,&y); CHKERRQ(ierr); 686840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 687840b8ebdSBarry Smith row = coloring->rows[k][l]; 688840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 689840b8ebdSBarry Smith y[row] *= scale[col]; 690840b8ebdSBarry Smith srow = row + start; 691840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 692840b8ebdSBarry Smith } 6933b28642cSBarry Smith ierr = VecRestoreArray(w2,&y); CHKERRQ(ierr); 694840b8ebdSBarry Smith } 6953b28642cSBarry Smith ierr = VecRestoreArray(x1,&xx); CHKERRQ(ierr); 696840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 697840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 6983a40ed3dSBarry Smith PetscFunctionReturn(0); 699840b8ebdSBarry Smith } 7003b28642cSBarry Smith 7013b28642cSBarry Smith 7023b28642cSBarry Smith 703