1840b8ebdSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*3a40ed3dSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.24 1997/10/14 20:15:42 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 #include "pinclude/pviewer.h" 15bbf0e169SBarry Smith 165615d1e5SSatish Balay #undef __FUNC__ 17d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView_Draw" 18005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer) 19005c665bSBarry Smith { 20005c665bSBarry Smith int ierr,i,j,pause; 21005c665bSBarry Smith PetscTruth isnull; 22005c665bSBarry Smith Draw draw; 23d4bb536fSBarry Smith double xr,yr,xl,yl,h,w,x,y,xc,yc,scale = 0.0; 24005c665bSBarry Smith DrawButton button; 25005c665bSBarry Smith 26*3a40ed3dSBarry Smith PetscFunctionBegin; 27005c665bSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 28*3a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 292bdab257SBarry Smith ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr); 30005c665bSBarry Smith 31005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 32005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 33005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 34005c665bSBarry Smith 35005c665bSBarry Smith /* loop over colors */ 36005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 37005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 38005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 39005c665bSBarry Smith x = fd->columnsforrow[i][j]; 40005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 41005c665bSBarry Smith } 42005c665bSBarry Smith } 432bdab257SBarry Smith ierr = DrawSynchronizedFlush(draw); CHKERRQ(ierr); 44005c665bSBarry Smith ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr); 45*3a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 46005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 472bdab257SBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 48005c665bSBarry Smith while (button != BUTTON_RIGHT) { 492bdab257SBarry Smith ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr); 50005c665bSBarry Smith if (button == BUTTON_LEFT) scale = .5; 51005c665bSBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 52005c665bSBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 53005c665bSBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 54005c665bSBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 55005c665bSBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 56005c665bSBarry Smith w *= scale; h *= scale; 57005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 58005c665bSBarry Smith /* loop over colors */ 59005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 60005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 61005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 62005c665bSBarry Smith x = fd->columnsforrow[i][j]; 63005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 64005c665bSBarry Smith } 65005c665bSBarry Smith } 66005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 672bdab257SBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 68005c665bSBarry Smith } 69005c665bSBarry Smith 70*3a40ed3dSBarry Smith PetscFunctionReturn(0); 71005c665bSBarry Smith } 72005c665bSBarry Smith 73005c665bSBarry Smith #undef __FUNC__ 74d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView" 75bbf0e169SBarry Smith /*@C 76639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 77bbf0e169SBarry Smith 78b4fc646aSLois Curfman McInnes Input Parameters: 79b4fc646aSLois Curfman McInnes . c - the coloring context 80b4fc646aSLois Curfman McInnes . viewer - visualization context 81b4fc646aSLois Curfman McInnes 82b4fc646aSLois Curfman McInnes Notes: 83b4fc646aSLois Curfman McInnes The available visualization contexts include 84b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_SELF - standard output (default) 85b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_WORLD - synchronized standard 86b4fc646aSLois Curfman McInnes $ output where only the first processor opens 87b4fc646aSLois Curfman McInnes $ the file. All other processors send their 88b4fc646aSLois Curfman McInnes $ data to the first processor to print. 89b4fc646aSLois Curfman McInnes $ VIEWER_DRAWX_WORLD - graphical display of nonzero structure 90bbf0e169SBarry Smith 91639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 92005c665bSBarry Smith 93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 94bbf0e169SBarry Smith @*/ 95b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer) 96bbf0e169SBarry Smith { 97005c665bSBarry Smith ViewerType vtype; 98639f9d9dSBarry Smith int i,j,format,ierr; 99b4fc646aSLois Curfman McInnes FILE *fd; 100bbf0e169SBarry Smith 101*3a40ed3dSBarry 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); 107005c665bSBarry Smith if (vtype == DRAW_VIEWER) { 108b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr); 109*3a40ed3dSBarry Smith PetscFunctionReturn(0); 110005c665bSBarry Smith } 111b4fc646aSLois Curfman McInnes else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 112b4fc646aSLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 113b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n"); 114b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Error tolerance=%g\n",c->error_rel); 115b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Umin=%g\n",c->umin); 116b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of colors=%d\n",c->ncolors); 117ae09f205SBarry Smith 118ae09f205SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 119ae09f205SBarry Smith if (format != VIEWER_FORMAT_ASCII_INFO) { 120b4fc646aSLois Curfman McInnes for ( i=0; i<c->ncolors; i++ ) { 121b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Information for color %d\n",i); 122b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of columns %d\n",c->ncolumns[i]); 123b4fc646aSLois Curfman McInnes for ( j=0; j<c->ncolumns[i]; j++ ) { 124b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d\n",c->columns[i][j]); 125639f9d9dSBarry Smith } 126b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of rows %d\n",c->nrows[i]); 127b4fc646aSLois Curfman McInnes for ( j=0; j<c->nrows[i]; j++ ) { 128b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]); 129b4fc646aSLois Curfman McInnes } 130bbf0e169SBarry Smith } 131bbf0e169SBarry Smith } 132bbf0e169SBarry Smith } 133*3a40ed3dSBarry 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 142639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 143639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 144639f9d9dSBarry Smith $ = error_rel*umin else 145639f9d9dSBarry Smith $ 146639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 147639f9d9dSBarry Smith 148639f9d9dSBarry Smith Input Parameters: 149b4fc646aSLois Curfman McInnes . coloring - the coloring context 150639f9d9dSBarry Smith . error_rel - relative error 151639f9d9dSBarry Smith . umin - minimum allowable u-value 152639f9d9dSBarry Smith 153b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 154b4fc646aSLois Curfman McInnes 155b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 156639f9d9dSBarry Smith @*/ 157639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 158639f9d9dSBarry Smith { 159*3a40ed3dSBarry Smith PetscFunctionBegin; 160639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 161639f9d9dSBarry Smith 162639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 163639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 164*3a40ed3dSBarry Smith PetscFunctionReturn(0); 165639f9d9dSBarry Smith } 166639f9d9dSBarry Smith 1675615d1e5SSatish Balay #undef __FUNC__ 168005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 169005c665bSBarry Smith /*@ 170e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 171e0907662SLois Curfman McInnes matrices. 172005c665bSBarry Smith 173005c665bSBarry Smith Input Parameters: 174b4fc646aSLois Curfman McInnes . coloring - the coloring context 175e0907662SLois Curfman McInnes . freq - frequency (default is 1) 176e0907662SLois Curfman McInnes 177e0907662SLois Curfman McInnes Notes: 178e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 179e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 180e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 181e0907662SLois Curfman McInnes <freq> nonlinear iterations. 182e0907662SLois Curfman McInnes 183e0907662SLois Curfman McInnes Options Database Keys: 184e0907662SLois Curfman McInnes $ -mat_fd_coloring_freq <freq> 185005c665bSBarry Smith 186b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 187005c665bSBarry Smith @*/ 188005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 189005c665bSBarry Smith { 190*3a40ed3dSBarry Smith PetscFunctionBegin; 191005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 192005c665bSBarry Smith 193005c665bSBarry Smith matfd->freq = freq; 194*3a40ed3dSBarry Smith PetscFunctionReturn(0); 195005c665bSBarry Smith } 196005c665bSBarry Smith 197005c665bSBarry Smith #undef __FUNC__ 198ff0cfa39SBarry Smith #define __FUNC__ "MatFDColoringGetFrequency" 199ff0cfa39SBarry Smith /*@ 200ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 201ff0cfa39SBarry Smith matrices. 202ff0cfa39SBarry Smith 203ff0cfa39SBarry Smith Input Parameters: 204ff0cfa39SBarry Smith . coloring - the coloring context 205ff0cfa39SBarry Smith 206ff0cfa39SBarry Smith Output Parameters: 207ff0cfa39SBarry Smith . freq - frequency (default is 1) 208ff0cfa39SBarry Smith 209ff0cfa39SBarry Smith Notes: 210ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 211ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 212ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 213ff0cfa39SBarry Smith <freq> nonlinear iterations. 214ff0cfa39SBarry Smith 215ff0cfa39SBarry Smith Options Database Keys: 216ff0cfa39SBarry Smith $ -mat_fd_coloring_freq <freq> 217ff0cfa39SBarry Smith 218ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 219ff0cfa39SBarry Smith @*/ 220ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 221ff0cfa39SBarry Smith { 222*3a40ed3dSBarry Smith PetscFunctionBegin; 223ff0cfa39SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 224ff0cfa39SBarry Smith 225ff0cfa39SBarry Smith *freq = matfd->freq; 226*3a40ed3dSBarry Smith PetscFunctionReturn(0); 227ff0cfa39SBarry Smith } 228ff0cfa39SBarry Smith 229ff0cfa39SBarry Smith #undef __FUNC__ 230005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 231005c665bSBarry Smith /*@ 232005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 233005c665bSBarry Smith 234005c665bSBarry Smith Input Parameters: 235b4fc646aSLois Curfman McInnes . coloring - the coloring context 236005c665bSBarry Smith . f - the function 237005c665bSBarry Smith . fctx - the function context 238005c665bSBarry Smith 239b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 240005c665bSBarry Smith @*/ 241840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 242005c665bSBarry Smith { 243*3a40ed3dSBarry Smith PetscFunctionBegin; 244005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 245005c665bSBarry Smith 246005c665bSBarry Smith matfd->f = f; 247005c665bSBarry Smith matfd->fctx = fctx; 248005c665bSBarry Smith 249*3a40ed3dSBarry Smith PetscFunctionReturn(0); 250005c665bSBarry Smith } 251005c665bSBarry Smith 252005c665bSBarry Smith #undef __FUNC__ 253d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 254639f9d9dSBarry Smith /*@ 255b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 256639f9d9dSBarry Smith the options database. 257639f9d9dSBarry Smith 258b4fc646aSLois Curfman McInnes The Jacobian is estimated with the differencing approximation 259639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 260639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 261ae09f205SBarry Smith $ = error_rel*umin else 262639f9d9dSBarry Smith $ 263639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 264639f9d9dSBarry Smith 265639f9d9dSBarry Smith Input Parameters: 266b4fc646aSLois Curfman McInnes . coloring - the coloring context 267639f9d9dSBarry Smith 268b4fc646aSLois Curfman McInnes Options Database Keys: 269b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_error <err>, where <err> is the square root 270b4fc646aSLois Curfman McInnes $ of relative error in the function 271b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_umin <umin>, where umin is described above 272e0907662SLois Curfman McInnes $ -mat_fd_coloring_freq <freq> where <freq> is the frequency of 273e0907662SLois Curfman McInnes $ computing a new Jacobian 274ca71c51bSBarry Smith $ -mat_fd_coloring_view 275ca71c51bSBarry Smith $ -mat_fd_coloring_view_info 276ca71c51bSBarry Smith $ -mat_fd_coloring_view_draw 277639f9d9dSBarry Smith 278b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 279639f9d9dSBarry Smith @*/ 280639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 281639f9d9dSBarry Smith { 282005c665bSBarry Smith int ierr,flag,freq = 1; 283639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 284*3a40ed3dSBarry Smith 285*3a40ed3dSBarry Smith PetscFunctionBegin; 286639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 287639f9d9dSBarry Smith 288639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 289639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 290639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 291005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 292005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 293005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 294639f9d9dSBarry Smith if (flag) { 295639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 296639f9d9dSBarry Smith } 297*3a40ed3dSBarry Smith PetscFunctionReturn(0); 298639f9d9dSBarry Smith } 299639f9d9dSBarry Smith 3005615d1e5SSatish Balay #undef __FUNC__ 301d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 302639f9d9dSBarry Smith /*@ 303639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 304639f9d9dSBarry Smith using coloring. 305639f9d9dSBarry Smith 306639f9d9dSBarry Smith Input Parameter: 307639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 308639f9d9dSBarry Smith 309639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 310639f9d9dSBarry Smith @*/ 311639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 312639f9d9dSBarry Smith { 313*3a40ed3dSBarry Smith PetscFunctionBegin; 314639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 315639f9d9dSBarry Smith 316e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel); 317e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin); 318e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq); 319005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view\n"); 320005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n"); 321ae09f205SBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n"); 322*3a40ed3dSBarry Smith PetscFunctionReturn(0); 323005c665bSBarry Smith } 324005c665bSBarry Smith 325005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 326005c665bSBarry Smith { 327005c665bSBarry Smith int ierr,flg; 328005c665bSBarry Smith 329*3a40ed3dSBarry Smith PetscFunctionBegin; 330005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 331005c665bSBarry Smith if (flg) { 332f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 333005c665bSBarry Smith } 334ae09f205SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 335ae09f205SBarry Smith if (flg) { 336f8590f6eSBarry Smith ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 337f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 338f8590f6eSBarry Smith ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339ae09f205SBarry Smith } 340005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 341005c665bSBarry Smith if (flg) { 342005c665bSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 343005c665bSBarry Smith ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 344005c665bSBarry Smith } 345*3a40ed3dSBarry Smith PetscFunctionReturn(0); 346bbf0e169SBarry Smith } 347bbf0e169SBarry Smith 3485615d1e5SSatish Balay #undef __FUNC__ 3495615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 350bbf0e169SBarry Smith /*@C 351639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 352639f9d9dSBarry Smith computation of Jacobians. 353bbf0e169SBarry Smith 354639f9d9dSBarry Smith Input Parameters: 355639f9d9dSBarry Smith . mat - the matrix containing the nonzero structure of the Jacobian 356639f9d9dSBarry Smith . iscoloring - the coloring of the matrix 357bbf0e169SBarry Smith 358bbf0e169SBarry Smith Output Parameter: 359639f9d9dSBarry Smith . color - the new coloring context 360bbf0e169SBarry Smith 361b4fc646aSLois Curfman McInnes Options Database Keys: 362b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view 363b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_draw 364b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_info 365639f9d9dSBarry Smith 366639f9d9dSBarry Smith .seealso: MatFDColoringDestroy() 367bbf0e169SBarry Smith @*/ 368639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 369bbf0e169SBarry Smith { 370639f9d9dSBarry Smith MatFDColoring c; 371639f9d9dSBarry Smith MPI_Comm comm; 372639f9d9dSBarry Smith int ierr,M,N; 373639f9d9dSBarry Smith 374*3a40ed3dSBarry Smith PetscFunctionBegin; 375639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 376e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 377639f9d9dSBarry Smith 378639f9d9dSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 379d4bb536fSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 380639f9d9dSBarry Smith PLogObjectCreate(c); 381639f9d9dSBarry Smith 382639f9d9dSBarry Smith if (mat->ops.fdcoloringcreate) { 383639f9d9dSBarry Smith ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 384639f9d9dSBarry Smith } else { 385e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 386639f9d9dSBarry Smith } 387639f9d9dSBarry Smith 388639f9d9dSBarry Smith c->error_rel = 1.e-8; 389ae09f205SBarry Smith c->umin = 1.e-6; 390005c665bSBarry Smith c->freq = 1; 391005c665bSBarry Smith 392005c665bSBarry Smith ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 393639f9d9dSBarry Smith 394639f9d9dSBarry Smith *color = c; 395639f9d9dSBarry Smith 396*3a40ed3dSBarry Smith PetscFunctionReturn(0); 397bbf0e169SBarry Smith } 398bbf0e169SBarry Smith 3995615d1e5SSatish Balay #undef __FUNC__ 400d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 401bbf0e169SBarry Smith /*@C 402639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 403639f9d9dSBarry Smith via MatFDColoringCreate(). 404bbf0e169SBarry Smith 405b4fc646aSLois Curfman McInnes Input Parameter: 406639f9d9dSBarry Smith . c - coloring context 407bbf0e169SBarry Smith 408639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 409bbf0e169SBarry Smith @*/ 410639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 411bbf0e169SBarry Smith { 412263760aaSBarry Smith int i,ierr; 413bbf0e169SBarry Smith 414*3a40ed3dSBarry Smith PetscFunctionBegin; 415*3a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 416d4bb536fSBarry Smith 417639f9d9dSBarry Smith 418639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 419639f9d9dSBarry Smith if (c->columns[i]) PetscFree(c->columns[i]); 420639f9d9dSBarry Smith if (c->rows[i]) PetscFree(c->rows[i]); 421639f9d9dSBarry Smith if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 422bbf0e169SBarry Smith } 423639f9d9dSBarry Smith PetscFree(c->ncolumns); 424639f9d9dSBarry Smith PetscFree(c->columns); 425639f9d9dSBarry Smith PetscFree(c->nrows); 426639f9d9dSBarry Smith PetscFree(c->rows); 427639f9d9dSBarry Smith PetscFree(c->columnsforrow); 428639f9d9dSBarry Smith PetscFree(c->scale); 429005c665bSBarry Smith if (c->w1) { 430005c665bSBarry Smith ierr = VecDestroy(c->w1); CHKERRQ(ierr); 431005c665bSBarry Smith ierr = VecDestroy(c->w2); CHKERRQ(ierr); 432005c665bSBarry Smith ierr = VecDestroy(c->w3); CHKERRQ(ierr); 433005c665bSBarry Smith } 434639f9d9dSBarry Smith PLogObjectDestroy(c); 435639f9d9dSBarry Smith PetscHeaderDestroy(c); 436*3a40ed3dSBarry Smith PetscFunctionReturn(0); 437bbf0e169SBarry Smith } 43843a90d84SBarry Smith 439005c665bSBarry Smith #include "snes.h" 440005c665bSBarry Smith 4415615d1e5SSatish Balay #undef __FUNC__ 4425615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 44343a90d84SBarry Smith /*@ 444e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 445e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 44643a90d84SBarry Smith 44743a90d84SBarry Smith Input Parameters: 44843a90d84SBarry Smith . mat - location to store Jacobian 44943a90d84SBarry Smith . coloring - coloring context created with MatFDColoringCreate() 45043a90d84SBarry Smith . x1 - location at which Jacobian is to be computed 451005c665bSBarry Smith . sctx - optional context required by function (actually a SNES context) 45243a90d84SBarry Smith 4538bba8e72SBarry Smith Options Database Keys: 4548bba8e72SBarry Smith $ -mat_fd_coloring_freq <freq> 4558bba8e72SBarry Smith 45643a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 45743a90d84SBarry Smith 45843a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 45943a90d84SBarry Smith @*/ 460005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 46143a90d84SBarry Smith { 462e0907662SLois Curfman McInnes int k,fg,ierr,N,start,end,l,row,col,srow; 46343a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 46443a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 46543a90d84SBarry Smith MPI_Comm comm = coloring->comm; 466005c665bSBarry Smith Vec w1,w2,w3; 467840b8ebdSBarry Smith int (*f)(void *,Vec,Vec,void *) = ( int (*)(void *,Vec,Vec,void *))coloring->f; 468005c665bSBarry Smith void *fctx = coloring->fctx; 469005c665bSBarry Smith 470*3a40ed3dSBarry Smith PetscFunctionBegin; 471e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 472e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 473e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 474e0907662SLois Curfman McInnes 475005c665bSBarry Smith 476005c665bSBarry Smith if (!coloring->w1) { 477005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 478005c665bSBarry Smith PLogObjectParent(coloring,coloring->w1); 479005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 480005c665bSBarry Smith PLogObjectParent(coloring,coloring->w2); 481005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 482005c665bSBarry Smith PLogObjectParent(coloring,coloring->w3); 483005c665bSBarry Smith } 484005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 48543a90d84SBarry Smith 486e0907662SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 487e0907662SLois Curfman McInnes if (fg) { 488e0907662SLois Curfman McInnes PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 489e0907662SLois Curfman McInnes } else { 49043a90d84SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 491e0907662SLois Curfman McInnes } 49243a90d84SBarry Smith 49343a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 49443a90d84SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 49543a90d84SBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 49643a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 49743a90d84SBarry Smith 49843a90d84SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 49943a90d84SBarry Smith /* 50043a90d84SBarry Smith Loop over each color 50143a90d84SBarry Smith */ 50243a90d84SBarry Smith 50343a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 50443a90d84SBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 50543a90d84SBarry Smith /* 50643a90d84SBarry Smith Loop over each column associated with color adding the 50743a90d84SBarry Smith perturbation to the vector w3. 50843a90d84SBarry Smith */ 50943a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 51043a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 51143a90d84SBarry Smith dx = xx[col-start]; 512ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 513*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 514ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 515ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 51643a90d84SBarry Smith #else 517ae09f205SBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 518ae09f205SBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 51943a90d84SBarry Smith #endif 52043a90d84SBarry Smith dx *= epsilon; 52143a90d84SBarry Smith wscale[col] = 1.0/dx; 52243a90d84SBarry Smith VecSetValues(w3,1,&col,&dx,ADD_VALUES); 52343a90d84SBarry Smith } 52443a90d84SBarry Smith VecRestoreArray(x1,&xx); 52543a90d84SBarry Smith /* 526e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 52743a90d84SBarry Smith */ 52843a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 52943a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 53043a90d84SBarry Smith /* Communicate scale to all processors */ 531*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 53243a90d84SBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 53343a90d84SBarry Smith #else 53443a90d84SBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 53543a90d84SBarry Smith #endif 53643a90d84SBarry Smith /* 537e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 53843a90d84SBarry Smith */ 53943a90d84SBarry Smith VecGetArray(w2,&y); 54043a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 54143a90d84SBarry Smith row = coloring->rows[k][l]; 54243a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 54343a90d84SBarry Smith y[row] *= scale[col]; 54443a90d84SBarry Smith srow = row + start; 54543a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 54643a90d84SBarry Smith } 54743a90d84SBarry Smith VecRestoreArray(w2,&y); 54843a90d84SBarry Smith } 54943a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 55043a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 551*3a40ed3dSBarry Smith PetscFunctionReturn(0); 55243a90d84SBarry Smith } 553840b8ebdSBarry Smith 554840b8ebdSBarry Smith #include "ts.h" 555840b8ebdSBarry Smith 556840b8ebdSBarry Smith #undef __FUNC__ 557840b8ebdSBarry Smith #define __FUNC__ "MatFDColoringApplyTS" 558840b8ebdSBarry Smith /*@ 559840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 560840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 561840b8ebdSBarry Smith 562840b8ebdSBarry Smith Input Parameters: 563840b8ebdSBarry Smith . mat - location to store Jacobian 564840b8ebdSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 565840b8ebdSBarry Smith . x1 - location at which Jacobian is to be computed 566840b8ebdSBarry Smith . sctx - optional context required by function (actually a SNES context) 567840b8ebdSBarry Smith 568840b8ebdSBarry Smith Options Database Keys: 569840b8ebdSBarry Smith $ -mat_fd_coloring_freq <freq> 570840b8ebdSBarry Smith 571840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 572840b8ebdSBarry Smith 573840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 574840b8ebdSBarry Smith @*/ 575840b8ebdSBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,double t,Vec x1,MatStructure *flag,void *sctx) 576840b8ebdSBarry Smith { 577840b8ebdSBarry Smith int k,fg,ierr,N,start,end,l,row,col,srow; 578840b8ebdSBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 579840b8ebdSBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 580840b8ebdSBarry Smith MPI_Comm comm = coloring->comm; 581840b8ebdSBarry Smith Vec w1,w2,w3; 582840b8ebdSBarry Smith int (*f)(void *,double,Vec,Vec,void *) = ( int (*)(void *,double,Vec,Vec,void *))coloring->f; 583840b8ebdSBarry Smith void *fctx = coloring->fctx; 584840b8ebdSBarry Smith 585*3a40ed3dSBarry Smith PetscFunctionBegin; 586840b8ebdSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE); 587840b8ebdSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 588840b8ebdSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE); 589840b8ebdSBarry Smith 590840b8ebdSBarry Smith 591840b8ebdSBarry Smith if (!coloring->w1) { 592840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 593840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w1); 594840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 595840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w2); 596840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 597840b8ebdSBarry Smith PLogObjectParent(coloring,coloring->w3); 598840b8ebdSBarry Smith } 599840b8ebdSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 600840b8ebdSBarry Smith 601840b8ebdSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 602840b8ebdSBarry Smith if (fg) { 603840b8ebdSBarry Smith PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n"); 604840b8ebdSBarry Smith } else { 605840b8ebdSBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 606840b8ebdSBarry Smith } 607840b8ebdSBarry Smith 608840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 609840b8ebdSBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 610840b8ebdSBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 611840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx); CHKERRQ(ierr); 612840b8ebdSBarry Smith 613840b8ebdSBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 614840b8ebdSBarry Smith /* 615840b8ebdSBarry Smith Loop over each color 616840b8ebdSBarry Smith */ 617840b8ebdSBarry Smith 618840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 619840b8ebdSBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 620840b8ebdSBarry Smith /* 621840b8ebdSBarry Smith Loop over each column associated with color adding the 622840b8ebdSBarry Smith perturbation to the vector w3. 623840b8ebdSBarry Smith */ 624840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 625840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 626840b8ebdSBarry Smith dx = xx[col-start]; 627840b8ebdSBarry Smith if (dx == 0.0) dx = 1.0; 628*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 629840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 630840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 631840b8ebdSBarry Smith #else 632840b8ebdSBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 633840b8ebdSBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 634840b8ebdSBarry Smith #endif 635840b8ebdSBarry Smith dx *= epsilon; 636840b8ebdSBarry Smith wscale[col] = 1.0/dx; 637840b8ebdSBarry Smith VecSetValues(w3,1,&col,&dx,ADD_VALUES); 638840b8ebdSBarry Smith } 639840b8ebdSBarry Smith VecRestoreArray(x1,&xx); 640840b8ebdSBarry Smith /* 641840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 642840b8ebdSBarry Smith */ 643840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx); CHKERRQ(ierr); 644840b8ebdSBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 645840b8ebdSBarry Smith /* Communicate scale to all processors */ 646*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 647840b8ebdSBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 648840b8ebdSBarry Smith #else 649840b8ebdSBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 650840b8ebdSBarry Smith #endif 651840b8ebdSBarry Smith /* 652840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 653840b8ebdSBarry Smith */ 654840b8ebdSBarry Smith VecGetArray(w2,&y); 655840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 656840b8ebdSBarry Smith row = coloring->rows[k][l]; 657840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 658840b8ebdSBarry Smith y[row] *= scale[col]; 659840b8ebdSBarry Smith srow = row + start; 660840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 661840b8ebdSBarry Smith } 662840b8ebdSBarry Smith VecRestoreArray(w2,&y); 663840b8ebdSBarry Smith } 664840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 665840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 666*3a40ed3dSBarry Smith PetscFunctionReturn(0); 667840b8ebdSBarry Smith } 668