1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*f8590f6eSBarry Smith static char vcid[] = "$Id: fdmatrix.c,v 1.17 1997/09/27 23:29:32 curfman Exp bsmith $"; 3bbf0e169SBarry Smith #endif 4bbf0e169SBarry Smith 5bbf0e169SBarry Smith /* 6639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 7639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 8bbf0e169SBarry Smith */ 9bbf0e169SBarry Smith 10bbf0e169SBarry Smith #include "petsc.h" 11bbf0e169SBarry Smith #include "src/mat/matimpl.h" /*I "mat.h" I*/ 12bbf0e169SBarry Smith #include "src/vec/vecimpl.h" 13bbf0e169SBarry Smith #include "pinclude/pviewer.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 25005c665bSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 26005c665bSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 27005c665bSBarry Smith ierr = DrawSyncClear(draw); CHKERRQ(ierr); 28005c665bSBarry Smith 29005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 30005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 31005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 32005c665bSBarry Smith 33005c665bSBarry Smith /* loop over colors */ 34005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 35005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 36005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 37005c665bSBarry Smith x = fd->columnsforrow[i][j]; 38005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 39005c665bSBarry Smith } 40005c665bSBarry Smith } 41005c665bSBarry Smith ierr = DrawSyncFlush(draw); CHKERRQ(ierr); 42005c665bSBarry Smith ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr); 43005c665bSBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 44005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 45005c665bSBarry Smith ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0); 46005c665bSBarry Smith while (button != BUTTON_RIGHT) { 47005c665bSBarry Smith ierr = DrawSyncClear(draw); CHKERRQ(ierr); 48005c665bSBarry Smith if (button == BUTTON_LEFT) scale = .5; 49005c665bSBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 50005c665bSBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 51005c665bSBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 52005c665bSBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 53005c665bSBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 54005c665bSBarry Smith w *= scale; h *= scale; 55005c665bSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 56005c665bSBarry Smith /* loop over colors */ 57005c665bSBarry Smith for (i=0; i<fd->ncolors; i++ ) { 58005c665bSBarry Smith for ( j=0; j<fd->nrows[i]; j++ ) { 59005c665bSBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 60005c665bSBarry Smith x = fd->columnsforrow[i][j]; 61005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 62005c665bSBarry Smith } 63005c665bSBarry Smith } 64005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 65005c665bSBarry Smith ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0); 66005c665bSBarry Smith } 67005c665bSBarry Smith 68005c665bSBarry Smith return 0; 69005c665bSBarry Smith } 70005c665bSBarry Smith 71005c665bSBarry Smith #undef __FUNC__ 72d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView" 73bbf0e169SBarry Smith /*@C 74639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 75bbf0e169SBarry Smith 76b4fc646aSLois Curfman McInnes Input Parameters: 77b4fc646aSLois Curfman McInnes . c - the coloring context 78b4fc646aSLois Curfman McInnes . viewer - visualization context 79b4fc646aSLois Curfman McInnes 80b4fc646aSLois Curfman McInnes Notes: 81b4fc646aSLois Curfman McInnes The available visualization contexts include 82b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_SELF - standard output (default) 83b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_WORLD - synchronized standard 84b4fc646aSLois Curfman McInnes $ output where only the first processor opens 85b4fc646aSLois Curfman McInnes $ the file. All other processors send their 86b4fc646aSLois Curfman McInnes $ data to the first processor to print. 87b4fc646aSLois Curfman McInnes $ VIEWER_DRAWX_WORLD - graphical display of nonzero structure 88bbf0e169SBarry Smith 89639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 90005c665bSBarry Smith 91b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 92bbf0e169SBarry Smith @*/ 93b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer) 94bbf0e169SBarry Smith { 95005c665bSBarry Smith ViewerType vtype; 96639f9d9dSBarry Smith int i,j,format,ierr; 97b4fc646aSLois Curfman McInnes FILE *fd; 98bbf0e169SBarry Smith 99b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 100b4fc646aSLois Curfman McInnes if (viewer) {PetscValidHeader(viewer);} 101b4fc646aSLois Curfman McInnes else {viewer = VIEWER_STDOUT_SELF;} 102bbf0e169SBarry Smith 103005c665bSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 104005c665bSBarry Smith if (vtype == DRAW_VIEWER) { 105b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr); 106005c665bSBarry Smith return 0; 107005c665bSBarry Smith } 108b4fc646aSLois Curfman McInnes else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 109b4fc646aSLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 110b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n"); 111b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Error tolerance=%g\n",c->error_rel); 112b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Umin=%g\n",c->umin); 113b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of colors=%d\n",c->ncolors); 114ae09f205SBarry Smith 115ae09f205SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 116ae09f205SBarry Smith if (format != VIEWER_FORMAT_ASCII_INFO) { 117b4fc646aSLois Curfman McInnes for ( i=0; i<c->ncolors; i++ ) { 118b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Information for color %d\n",i); 119b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of columns %d\n",c->ncolumns[i]); 120b4fc646aSLois Curfman McInnes for ( j=0; j<c->ncolumns[i]; j++ ) { 121b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d\n",c->columns[i][j]); 122639f9d9dSBarry Smith } 123b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of rows %d\n",c->nrows[i]); 124b4fc646aSLois Curfman McInnes for ( j=0; j<c->nrows[i]; j++ ) { 125b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]); 126b4fc646aSLois Curfman McInnes } 127bbf0e169SBarry Smith } 128bbf0e169SBarry Smith } 129bbf0e169SBarry Smith } 130639f9d9dSBarry Smith return 0; 131639f9d9dSBarry Smith } 132639f9d9dSBarry Smith 1335615d1e5SSatish Balay #undef __FUNC__ 1345615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters" 135639f9d9dSBarry Smith /*@ 136b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 137b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 138639f9d9dSBarry Smith 139639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 140639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 141639f9d9dSBarry Smith $ = error_rel*umin else 142639f9d9dSBarry Smith $ 143639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 144639f9d9dSBarry Smith 145639f9d9dSBarry Smith Input Parameters: 146b4fc646aSLois Curfman McInnes . coloring - the coloring context 147639f9d9dSBarry Smith . error_rel - relative error 148639f9d9dSBarry Smith . umin - minimum allowable u-value 149639f9d9dSBarry Smith 150b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 151b4fc646aSLois Curfman McInnes 152b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 153639f9d9dSBarry Smith @*/ 154639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 155639f9d9dSBarry Smith { 156639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 157639f9d9dSBarry Smith 158639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 159639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 160639f9d9dSBarry Smith return 0; 161639f9d9dSBarry Smith } 162639f9d9dSBarry Smith 1635615d1e5SSatish Balay #undef __FUNC__ 164005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 165005c665bSBarry Smith /*@ 166e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 167e0907662SLois Curfman McInnes matrices. 168005c665bSBarry Smith 169005c665bSBarry Smith Input Parameters: 170b4fc646aSLois Curfman McInnes . coloring - the coloring context 171e0907662SLois Curfman McInnes . freq - frequency (default is 1) 172e0907662SLois Curfman McInnes 173e0907662SLois Curfman McInnes Notes: 174e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 175e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 176e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 177e0907662SLois Curfman McInnes <freq> nonlinear iterations. 178e0907662SLois Curfman McInnes 179e0907662SLois Curfman McInnes Options Database Keys: 180e0907662SLois Curfman McInnes $ -mat_fd_coloring_freq <freq> 181005c665bSBarry Smith 182b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 183005c665bSBarry Smith @*/ 184005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 185005c665bSBarry Smith { 186005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 187005c665bSBarry Smith 188005c665bSBarry Smith matfd->freq = freq; 189005c665bSBarry Smith return 0; 190005c665bSBarry Smith } 191005c665bSBarry Smith 192005c665bSBarry Smith #undef __FUNC__ 193005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 194005c665bSBarry Smith /*@ 195005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 196005c665bSBarry Smith 197005c665bSBarry Smith Input Parameters: 198b4fc646aSLois Curfman McInnes . coloring - the coloring context 199005c665bSBarry Smith . f - the function 200005c665bSBarry Smith . fctx - the function context 201005c665bSBarry Smith 202b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 203005c665bSBarry Smith @*/ 204005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx) 205005c665bSBarry Smith { 206005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 207005c665bSBarry Smith 208005c665bSBarry Smith matfd->f = f; 209005c665bSBarry Smith matfd->fctx = fctx; 210005c665bSBarry Smith 211005c665bSBarry Smith return 0; 212005c665bSBarry Smith } 213005c665bSBarry Smith 214005c665bSBarry Smith #undef __FUNC__ 215d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 216639f9d9dSBarry Smith /*@ 217b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 218639f9d9dSBarry Smith the options database. 219639f9d9dSBarry Smith 220b4fc646aSLois Curfman McInnes The Jacobian is estimated with the differencing approximation 221639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 222639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 223ae09f205SBarry Smith $ = error_rel*umin else 224639f9d9dSBarry Smith $ 225639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 226639f9d9dSBarry Smith 227639f9d9dSBarry Smith Input Parameters: 228b4fc646aSLois Curfman McInnes . coloring - the coloring context 229639f9d9dSBarry Smith 230b4fc646aSLois Curfman McInnes Options Database Keys: 231b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_error <err>, where <err> is the square root 232b4fc646aSLois Curfman McInnes $ of relative error in the function 233b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_umin <umin>, where umin is described above 234e0907662SLois Curfman McInnes $ -mat_fd_coloring_freq <freq> where <freq> is the frequency of 235e0907662SLois Curfman McInnes $ computing a new Jacobian 236639f9d9dSBarry Smith 237b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 238639f9d9dSBarry Smith @*/ 239639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 240639f9d9dSBarry Smith { 241005c665bSBarry Smith int ierr,flag,freq = 1; 242639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 243639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 244639f9d9dSBarry Smith 245639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 246639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 247639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 248005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 249005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 250005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 251639f9d9dSBarry Smith if (flag) { 252639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 253639f9d9dSBarry Smith } 254639f9d9dSBarry Smith return 0; 255639f9d9dSBarry Smith } 256639f9d9dSBarry Smith 2575615d1e5SSatish Balay #undef __FUNC__ 258d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 259639f9d9dSBarry Smith /*@ 260639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 261639f9d9dSBarry Smith using coloring. 262639f9d9dSBarry Smith 263639f9d9dSBarry Smith Input Parameter: 264639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 265639f9d9dSBarry Smith 266639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 267639f9d9dSBarry Smith @*/ 268639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 269639f9d9dSBarry Smith { 270639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 271639f9d9dSBarry Smith 272e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel); 273e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin); 274e0907662SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq); 275005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view\n"); 276005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n"); 277ae09f205SBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n"); 278005c665bSBarry Smith return 0; 279005c665bSBarry Smith } 280005c665bSBarry Smith 281005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 282005c665bSBarry Smith { 283005c665bSBarry Smith int ierr,flg; 284005c665bSBarry Smith 285005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 286005c665bSBarry Smith if (flg) { 287*f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 288005c665bSBarry Smith } 289ae09f205SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 290ae09f205SBarry Smith if (flg) { 291*f8590f6eSBarry Smith ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 292*f8590f6eSBarry Smith ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm)); CHKERRQ(ierr); 293*f8590f6eSBarry Smith ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 294ae09f205SBarry Smith } 295005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 296005c665bSBarry Smith if (flg) { 297005c665bSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 298005c665bSBarry Smith ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 299005c665bSBarry Smith } 300bbf0e169SBarry Smith return 0; 301bbf0e169SBarry Smith } 302bbf0e169SBarry Smith 3035615d1e5SSatish Balay #undef __FUNC__ 3045615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 305bbf0e169SBarry Smith /*@C 306639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 307639f9d9dSBarry Smith computation of Jacobians. 308bbf0e169SBarry Smith 309639f9d9dSBarry Smith Input Parameters: 310639f9d9dSBarry Smith . mat - the matrix containing the nonzero structure of the Jacobian 311639f9d9dSBarry Smith . iscoloring - the coloring of the matrix 312bbf0e169SBarry Smith 313bbf0e169SBarry Smith Output Parameter: 314639f9d9dSBarry Smith . color - the new coloring context 315bbf0e169SBarry Smith 316b4fc646aSLois Curfman McInnes Options Database Keys: 317b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view 318b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_draw 319b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_info 320639f9d9dSBarry Smith 321639f9d9dSBarry Smith .seealso: MatFDColoringDestroy() 322bbf0e169SBarry Smith @*/ 323639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 324bbf0e169SBarry Smith { 325639f9d9dSBarry Smith MatFDColoring c; 326639f9d9dSBarry Smith MPI_Comm comm; 327639f9d9dSBarry Smith int ierr,M,N; 328639f9d9dSBarry Smith 329639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 330e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 331639f9d9dSBarry Smith 332639f9d9dSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 333d4bb536fSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 334639f9d9dSBarry Smith PLogObjectCreate(c); 335639f9d9dSBarry Smith 336639f9d9dSBarry Smith if (mat->ops.fdcoloringcreate) { 337639f9d9dSBarry Smith ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 338639f9d9dSBarry Smith } else { 339e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 340639f9d9dSBarry Smith } 341639f9d9dSBarry Smith 342639f9d9dSBarry Smith c->error_rel = 1.e-8; 343ae09f205SBarry Smith c->umin = 1.e-6; 344005c665bSBarry Smith c->freq = 1; 345005c665bSBarry Smith 346005c665bSBarry Smith ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 347639f9d9dSBarry Smith 348639f9d9dSBarry Smith *color = c; 349639f9d9dSBarry Smith 350bbf0e169SBarry Smith return 0; 351bbf0e169SBarry Smith } 352bbf0e169SBarry Smith 3535615d1e5SSatish Balay #undef __FUNC__ 354d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 355bbf0e169SBarry Smith /*@C 356639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 357639f9d9dSBarry Smith via MatFDColoringCreate(). 358bbf0e169SBarry Smith 359b4fc646aSLois Curfman McInnes Input Parameter: 360639f9d9dSBarry Smith . c - coloring context 361bbf0e169SBarry Smith 362639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 363bbf0e169SBarry Smith @*/ 364639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 365bbf0e169SBarry Smith { 366639f9d9dSBarry Smith int i,ierr,flag; 367bbf0e169SBarry Smith 368d4bb536fSBarry Smith if (--c->refct > 0) return 0; 369d4bb536fSBarry Smith 370639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag); 371639f9d9dSBarry Smith if (flag) { 372*f8590f6eSBarry Smith ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr); 373bbf0e169SBarry Smith } 374639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag); 375639f9d9dSBarry Smith if (flag) { 376*f8590f6eSBarry Smith ierr = ViewerPushFormat(VIEWER_STDOUT_(c->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 377*f8590f6eSBarry Smith ierr = MatFDColoringView(c,VIEWER_STDOUT_(c->comm));CHKERRQ(ierr); 378*f8590f6eSBarry Smith ierr = ViewerPopFormat(VIEWER_STDOUT_(c->comm)); 379bbf0e169SBarry Smith } 380639f9d9dSBarry Smith 381639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 382639f9d9dSBarry Smith if (c->columns[i]) PetscFree(c->columns[i]); 383639f9d9dSBarry Smith if (c->rows[i]) PetscFree(c->rows[i]); 384639f9d9dSBarry Smith if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 385bbf0e169SBarry Smith } 386639f9d9dSBarry Smith PetscFree(c->ncolumns); 387639f9d9dSBarry Smith PetscFree(c->columns); 388639f9d9dSBarry Smith PetscFree(c->nrows); 389639f9d9dSBarry Smith PetscFree(c->rows); 390639f9d9dSBarry Smith PetscFree(c->columnsforrow); 391639f9d9dSBarry Smith PetscFree(c->scale); 392005c665bSBarry Smith if (c->w1) { 393005c665bSBarry Smith ierr = VecDestroy(c->w1); CHKERRQ(ierr); 394005c665bSBarry Smith ierr = VecDestroy(c->w2); CHKERRQ(ierr); 395005c665bSBarry Smith ierr = VecDestroy(c->w3); CHKERRQ(ierr); 396005c665bSBarry Smith } 397639f9d9dSBarry Smith PLogObjectDestroy(c); 398639f9d9dSBarry Smith PetscHeaderDestroy(c); 399bbf0e169SBarry Smith return 0; 400bbf0e169SBarry Smith } 40143a90d84SBarry Smith 402005c665bSBarry Smith #include "snes.h" 403005c665bSBarry Smith 4045615d1e5SSatish Balay #undef __FUNC__ 4055615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 40643a90d84SBarry Smith /*@ 407e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 408e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 40943a90d84SBarry Smith 41043a90d84SBarry Smith Input Parameters: 41143a90d84SBarry Smith . mat - location to store Jacobian 41243a90d84SBarry Smith . coloring - coloring context created with MatFDColoringCreate() 41343a90d84SBarry Smith . x1 - location at which Jacobian is to be computed 414005c665bSBarry Smith . sctx - optional context required by function (actually a SNES context) 41543a90d84SBarry Smith 41643a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 41743a90d84SBarry Smith 41843a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 41943a90d84SBarry Smith @*/ 420005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 42143a90d84SBarry Smith { 422e0907662SLois Curfman McInnes int k,fg,ierr,N,start,end,l,row,col,srow; 42343a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 42443a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 42543a90d84SBarry Smith MPI_Comm comm = coloring->comm; 426005c665bSBarry Smith Vec w1,w2,w3; 427dd9fa9a5SLois Curfman McInnes int (*f)(void *,Vec,Vec,void *) = coloring->f; 428005c665bSBarry Smith void *fctx = coloring->fctx; 429005c665bSBarry Smith 430e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 431e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 432e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 433e0907662SLois Curfman McInnes 434d4bb536fSBarry Smith /* 435005c665bSBarry Smith ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr); 436005c665bSBarry Smith if ((freq > 1) && ((it % freq) != 1)) { 437005c665bSBarry Smith PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq); 438005c665bSBarry Smith *flag = SAME_PRECONDITIONER; 439005c665bSBarry Smith return 0; 440005c665bSBarry Smith } else { 441005c665bSBarry Smith PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq); 442005c665bSBarry Smith *flag = SAME_NONZERO_PATTERN; 443d4bb536fSBarry Smith }*/ 444005c665bSBarry Smith 445005c665bSBarry Smith if (!coloring->w1) { 446005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 447005c665bSBarry Smith PLogObjectParent(coloring,coloring->w1); 448005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 449005c665bSBarry Smith PLogObjectParent(coloring,coloring->w2); 450005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 451005c665bSBarry Smith PLogObjectParent(coloring,coloring->w3); 452005c665bSBarry Smith } 453005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 45443a90d84SBarry Smith 455e0907662SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&fg); CHKERRQ(ierr); 456e0907662SLois Curfman McInnes if (fg) { 457e0907662SLois Curfman McInnes PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 458e0907662SLois Curfman McInnes } else { 45943a90d84SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 460e0907662SLois Curfman McInnes } 46143a90d84SBarry Smith 46243a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 46343a90d84SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 46443a90d84SBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 46543a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 46643a90d84SBarry Smith 46743a90d84SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 46843a90d84SBarry Smith /* 46943a90d84SBarry Smith Loop over each color 47043a90d84SBarry Smith */ 47143a90d84SBarry Smith 47243a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 47343a90d84SBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 47443a90d84SBarry Smith /* 47543a90d84SBarry Smith Loop over each column associated with color adding the 47643a90d84SBarry Smith perturbation to the vector w3. 47743a90d84SBarry Smith */ 47843a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 47943a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 48043a90d84SBarry Smith dx = xx[col-start]; 481ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 48243a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 483ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 484ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 48543a90d84SBarry Smith #else 486ae09f205SBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 487ae09f205SBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 48843a90d84SBarry Smith #endif 48943a90d84SBarry Smith dx *= epsilon; 49043a90d84SBarry Smith wscale[col] = 1.0/dx; 49143a90d84SBarry Smith VecSetValues(w3,1,&col,&dx,ADD_VALUES); 49243a90d84SBarry Smith } 49343a90d84SBarry Smith VecRestoreArray(x1,&xx); 49443a90d84SBarry Smith /* 495e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 49643a90d84SBarry Smith */ 49743a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 49843a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 49943a90d84SBarry Smith /* Communicate scale to all processors */ 50043a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 50143a90d84SBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 50243a90d84SBarry Smith #else 50343a90d84SBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 50443a90d84SBarry Smith #endif 50543a90d84SBarry Smith /* 506e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 50743a90d84SBarry Smith */ 50843a90d84SBarry Smith VecGetArray(w2,&y); 50943a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 51043a90d84SBarry Smith row = coloring->rows[k][l]; 51143a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 51243a90d84SBarry Smith y[row] *= scale[col]; 51343a90d84SBarry Smith srow = row + start; 51443a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 51543a90d84SBarry Smith } 51643a90d84SBarry Smith VecRestoreArray(w2,&y); 51743a90d84SBarry Smith } 51843a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 51943a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 52043a90d84SBarry Smith return 0; 52143a90d84SBarry Smith } 522