1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*dd9fa9a5SLois Curfman McInnes static char vcid[] = "$Id: fdmatrix.c,v 1.13 1997/08/22 15:13:06 bsmith Exp curfman $"; 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 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]; 39005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 40005c665bSBarry Smith } 41005c665bSBarry Smith } 42005c665bSBarry Smith ierr = DrawSyncFlush(draw); CHKERRQ(ierr); 43005c665bSBarry Smith ierr = DrawGetPause(draw,&pause); CHKERRQ(ierr); 44005c665bSBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 45005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 46005c665bSBarry Smith ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0); 47005c665bSBarry Smith while (button != BUTTON_RIGHT) { 48005c665bSBarry Smith ierr = DrawSyncClear(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]; 62005c665bSBarry Smith DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1); 63005c665bSBarry Smith } 64005c665bSBarry Smith } 65005c665bSBarry Smith ierr = DrawCheckResizedWindow(draw); 66005c665bSBarry Smith ierr = DrawSyncGetMouseButton(draw,&button,&xc,&yc,0,0); 67005c665bSBarry Smith } 68005c665bSBarry Smith 69005c665bSBarry Smith return 0; 70005c665bSBarry Smith } 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 78639f9d9dSBarry Smith Input Parameter: 79639f9d9dSBarry Smith . color - the coloring context 80bbf0e169SBarry Smith 81639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 82005c665bSBarry Smith 83bbf0e169SBarry Smith @*/ 84639f9d9dSBarry Smith int MatFDColoringView(MatFDColoring color,Viewer viewer) 85bbf0e169SBarry Smith { 86005c665bSBarry Smith ViewerType vtype; 87639f9d9dSBarry Smith int i,j,format,ierr; 88bbf0e169SBarry Smith 89639f9d9dSBarry Smith if (!viewer) viewer = VIEWER_STDOUT_WORLD; 90bbf0e169SBarry Smith 91005c665bSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 92005c665bSBarry Smith if (vtype == DRAW_VIEWER) { 93005c665bSBarry Smith ierr = MatFDColoringView_Draw(color,viewer); CHKERRQ(ierr); 94005c665bSBarry Smith return 0; 95005c665bSBarry Smith } 96005c665bSBarry Smith 97bbf0e169SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 98639f9d9dSBarry Smith 99639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 100639f9d9dSBarry Smith printf("MatFDColoring Object:\n"); 101639f9d9dSBarry Smith printf(" Error tolerance %g\n",color->error_rel); 102639f9d9dSBarry Smith printf(" umin %g\n",color->umin); 103639f9d9dSBarry Smith } else { 104639f9d9dSBarry Smith printf("MatFDColoring Object:\n"); 105639f9d9dSBarry Smith printf(" Error tolerance %g\n",color->error_rel); 106639f9d9dSBarry Smith printf(" umin %g\n",color->umin); 107639f9d9dSBarry Smith printf("Number of colors %d\n",color->ncolors); 108639f9d9dSBarry Smith for ( i=0; i<color->ncolors; i++ ) { 109639f9d9dSBarry Smith printf("Information for color %d\n",i); 110639f9d9dSBarry Smith printf("Number of columns %d\n",color->ncolumns[i]); 111639f9d9dSBarry Smith for ( j=0; j<color->ncolumns[i]; j++ ) { 112639f9d9dSBarry Smith printf(" %d\n",color->columns[i][j]); 113639f9d9dSBarry Smith } 114639f9d9dSBarry Smith printf("Number of rows %d\n",color->nrows[i]); 115639f9d9dSBarry Smith for ( j=0; j<color->nrows[i]; j++ ) { 116639f9d9dSBarry Smith printf(" %d %d \n",color->rows[i][j],color->columnsforrow[i][j]); 117bbf0e169SBarry Smith } 118bbf0e169SBarry Smith } 119bbf0e169SBarry Smith } 120639f9d9dSBarry Smith return 0; 121639f9d9dSBarry Smith } 122639f9d9dSBarry Smith 1235615d1e5SSatish Balay #undef __FUNC__ 1245615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters" 125639f9d9dSBarry Smith /*@ 126639f9d9dSBarry Smith MatFDColoringSetParameters - Sets the parameters for the approximation of 127639f9d9dSBarry Smith Jacobian using finite differences. 128639f9d9dSBarry Smith 129639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 130639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 131639f9d9dSBarry Smith $ = error_rel*umin else 132639f9d9dSBarry Smith $ 133639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 134639f9d9dSBarry Smith 135639f9d9dSBarry Smith Input Parameters: 136639f9d9dSBarry Smith . coloring - the jacobian coloring context 137639f9d9dSBarry Smith . error_rel - relative error 138639f9d9dSBarry Smith . umin - minimum allowable u-value 139639f9d9dSBarry Smith 140639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters 141639f9d9dSBarry Smith @*/ 142639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 143639f9d9dSBarry Smith { 144639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 145639f9d9dSBarry Smith 146639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 147639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 148639f9d9dSBarry Smith return 0; 149639f9d9dSBarry Smith } 150639f9d9dSBarry Smith 1515615d1e5SSatish Balay #undef __FUNC__ 152005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 153005c665bSBarry Smith /*@ 154005c665bSBarry Smith MatFDColoringSetFrequency - Sets the frequency at which new Jacobians 155005c665bSBarry Smith are computed. 156005c665bSBarry Smith 157005c665bSBarry Smith Input Parameters: 158005c665bSBarry Smith . coloring - the jacobian coloring context 159005c665bSBarry Smith . freq - frequency 160005c665bSBarry Smith 161005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters 162005c665bSBarry Smith @*/ 163005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 164005c665bSBarry Smith { 165005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 166005c665bSBarry Smith 167005c665bSBarry Smith matfd->freq = freq; 168005c665bSBarry Smith return 0; 169005c665bSBarry Smith } 170005c665bSBarry Smith 171005c665bSBarry Smith #undef __FUNC__ 172005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 173005c665bSBarry Smith /*@ 174005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 175005c665bSBarry Smith 176005c665bSBarry Smith Input Parameters: 177005c665bSBarry Smith . coloring - the jacobian coloring context 178005c665bSBarry Smith . f - the function 179005c665bSBarry Smith . fctx - the function context 180005c665bSBarry Smith 181005c665bSBarry Smith 182005c665bSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters, function 183005c665bSBarry Smith @*/ 184005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx) 185005c665bSBarry Smith { 186005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 187005c665bSBarry Smith 188005c665bSBarry Smith matfd->f = f; 189005c665bSBarry Smith matfd->fctx = fctx; 190005c665bSBarry Smith 191005c665bSBarry Smith return 0; 192005c665bSBarry Smith } 193005c665bSBarry Smith 194005c665bSBarry Smith #undef __FUNC__ 195d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 196639f9d9dSBarry Smith /*@ 197639f9d9dSBarry Smith MatFDColoringSetFromOptions - Set coloring finite difference parameters from 198639f9d9dSBarry Smith the options database. 199639f9d9dSBarry Smith 200639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 201639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 202639f9d9dSBarry Smith $ = error_rel*.1 else 203639f9d9dSBarry Smith $ 204639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 205639f9d9dSBarry Smith 206639f9d9dSBarry Smith Input Parameters: 207639f9d9dSBarry Smith . coloring - the jacobian coloring context 208639f9d9dSBarry Smith 209639f9d9dSBarry Smith Options Database: 210639f9d9dSBarry Smith . -mat_fd_coloring_error square root of relative error in function 211639f9d9dSBarry Smith . -mat_fd_coloring_umin see above 212005c665bSBarry Smith . -mat_fd_coloring_freq frequency at which Jacobian is computed 213639f9d9dSBarry Smith 214639f9d9dSBarry Smith .keywords: SNES, Jacobian, finite differences, parameters 215639f9d9dSBarry Smith @*/ 216639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 217639f9d9dSBarry Smith { 218005c665bSBarry Smith int ierr,flag,freq = 1; 219639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 220639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 221639f9d9dSBarry Smith 222639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 223639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 224639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 225005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 226005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 227005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 228639f9d9dSBarry Smith if (flag) { 229639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 230639f9d9dSBarry Smith } 231639f9d9dSBarry Smith return 0; 232639f9d9dSBarry Smith } 233639f9d9dSBarry Smith 2345615d1e5SSatish Balay #undef __FUNC__ 235d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 236639f9d9dSBarry Smith /*@ 237639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 238639f9d9dSBarry Smith using coloring. 239639f9d9dSBarry Smith 240639f9d9dSBarry Smith Input Parameter: 241639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 242639f9d9dSBarry Smith 243639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 244639f9d9dSBarry Smith @*/ 245639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 246639f9d9dSBarry Smith { 247639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 248639f9d9dSBarry Smith 2490f665d81SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to 1.e-8\n"); 2500f665d81SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to 1.e-8\n"); 251005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults 1\n"); 252005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view\n"); 253005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n"); 254005c665bSBarry Smith return 0; 255005c665bSBarry Smith } 256005c665bSBarry Smith 257005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 258005c665bSBarry Smith { 259005c665bSBarry Smith int ierr,flg; 260005c665bSBarry Smith 261005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 262005c665bSBarry Smith if (flg) { 263005c665bSBarry Smith Viewer viewer; 264005c665bSBarry Smith ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr); 265005c665bSBarry Smith ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr); 266005c665bSBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 267005c665bSBarry Smith } 268005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 269005c665bSBarry Smith if (flg) { 270005c665bSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 271005c665bSBarry Smith ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 272005c665bSBarry Smith } 273bbf0e169SBarry Smith return 0; 274bbf0e169SBarry Smith } 275bbf0e169SBarry Smith 2765615d1e5SSatish Balay #undef __FUNC__ 2775615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 278bbf0e169SBarry Smith /*@C 279639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 280639f9d9dSBarry Smith computation of Jacobians. 281bbf0e169SBarry Smith 282639f9d9dSBarry Smith Input Parameters: 283639f9d9dSBarry Smith . mat - the matrix containing the nonzero structure of the Jacobian 284639f9d9dSBarry Smith . iscoloring - the coloring of the matrix 285bbf0e169SBarry Smith 286bbf0e169SBarry Smith Output Parameter: 287639f9d9dSBarry Smith . color - the new coloring context 288bbf0e169SBarry Smith 289005c665bSBarry Smith Options Database: 290005c665bSBarry Smith . -mat_fd_coloring_error square root of relative error in function 291005c665bSBarry Smith . -mat_fd_coloring_umin see above 292005c665bSBarry Smith . -mat_fd_coloring_view 293005c665bSBarry Smith . -mat_fd_coloring_view_draw 294639f9d9dSBarry Smith 295639f9d9dSBarry Smith .seealso: MatFDColoringDestroy() 296bbf0e169SBarry Smith @*/ 297639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 298bbf0e169SBarry Smith { 299639f9d9dSBarry Smith MatFDColoring c; 300639f9d9dSBarry Smith MPI_Comm comm; 301639f9d9dSBarry Smith int ierr,M,N; 302639f9d9dSBarry Smith 303639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 304e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 305639f9d9dSBarry Smith 306639f9d9dSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 307d4bb536fSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 308639f9d9dSBarry Smith PLogObjectCreate(c); 309639f9d9dSBarry Smith 310639f9d9dSBarry Smith if (mat->ops.fdcoloringcreate) { 311639f9d9dSBarry Smith ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 312639f9d9dSBarry Smith } else { 313e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 314639f9d9dSBarry Smith } 315639f9d9dSBarry Smith 316639f9d9dSBarry Smith c->error_rel = 1.e-8; 317639f9d9dSBarry Smith c->umin = 1.e-5; 318005c665bSBarry Smith c->freq = 1; 319005c665bSBarry Smith 320005c665bSBarry Smith ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 321639f9d9dSBarry Smith 322639f9d9dSBarry Smith *color = c; 323639f9d9dSBarry Smith 324bbf0e169SBarry Smith return 0; 325bbf0e169SBarry Smith } 326bbf0e169SBarry Smith 3275615d1e5SSatish Balay #undef __FUNC__ 328d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 329bbf0e169SBarry Smith /*@C 330639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 331639f9d9dSBarry Smith via MatFDColoringCreate(). 332bbf0e169SBarry Smith 333639f9d9dSBarry Smith Input Paramter: 334639f9d9dSBarry Smith . c - coloring context 335bbf0e169SBarry Smith 336639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 337bbf0e169SBarry Smith @*/ 338639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 339bbf0e169SBarry Smith { 340639f9d9dSBarry Smith int i,ierr,flag; 341bbf0e169SBarry Smith 342d4bb536fSBarry Smith if (--c->refct > 0) return 0; 343d4bb536fSBarry Smith 344639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag); 345639f9d9dSBarry Smith if (flag) { 346bbf0e169SBarry Smith Viewer viewer; 347639f9d9dSBarry Smith ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr); 348639f9d9dSBarry Smith ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr); 349bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 350bbf0e169SBarry Smith } 351639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag); 352639f9d9dSBarry Smith if (flag) { 353bbf0e169SBarry Smith Viewer viewer; 354639f9d9dSBarry Smith ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr); 355639f9d9dSBarry Smith ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr); 356639f9d9dSBarry Smith ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr); 357639f9d9dSBarry Smith ierr = ViewerPopFormat(viewer); 358bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 359bbf0e169SBarry Smith } 360639f9d9dSBarry Smith 361639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 362639f9d9dSBarry Smith if (c->columns[i]) PetscFree(c->columns[i]); 363639f9d9dSBarry Smith if (c->rows[i]) PetscFree(c->rows[i]); 364639f9d9dSBarry Smith if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 365bbf0e169SBarry Smith } 366639f9d9dSBarry Smith PetscFree(c->ncolumns); 367639f9d9dSBarry Smith PetscFree(c->columns); 368639f9d9dSBarry Smith PetscFree(c->nrows); 369639f9d9dSBarry Smith PetscFree(c->rows); 370639f9d9dSBarry Smith PetscFree(c->columnsforrow); 371639f9d9dSBarry Smith PetscFree(c->scale); 372005c665bSBarry Smith if (c->w1) { 373005c665bSBarry Smith ierr = VecDestroy(c->w1); CHKERRQ(ierr); 374005c665bSBarry Smith ierr = VecDestroy(c->w2); CHKERRQ(ierr); 375005c665bSBarry Smith ierr = VecDestroy(c->w3); CHKERRQ(ierr); 376005c665bSBarry Smith } 377639f9d9dSBarry Smith PLogObjectDestroy(c); 378639f9d9dSBarry Smith PetscHeaderDestroy(c); 379bbf0e169SBarry Smith return 0; 380bbf0e169SBarry Smith } 38143a90d84SBarry Smith 382005c665bSBarry Smith #include "snes.h" 383005c665bSBarry Smith 3845615d1e5SSatish Balay #undef __FUNC__ 3855615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 38643a90d84SBarry Smith /*@ 38743a90d84SBarry Smith MatFDColoringApply - Given a matrix for which a MatFDColoring has been created, 38843a90d84SBarry Smith computes the Jacobian for a function via finite differences. 38943a90d84SBarry Smith 39043a90d84SBarry Smith Input Parameters: 39143a90d84SBarry Smith . mat - location to store Jacobian 39243a90d84SBarry Smith . coloring - coloring context created with MatFDColoringCreate() 39343a90d84SBarry Smith . x1 - location at which Jacobian is to be computed 394005c665bSBarry Smith . sctx - optional context required by function (actually a SNES context) 39543a90d84SBarry Smith 39643a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 39743a90d84SBarry Smith 39843a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 39943a90d84SBarry Smith @*/ 400005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 40143a90d84SBarry Smith { 402*dd9fa9a5SLois Curfman McInnes int k, ierr,N,start,end,l,row,col,srow; 40343a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 40443a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 40543a90d84SBarry Smith MPI_Comm comm = coloring->comm; 406005c665bSBarry Smith Vec w1,w2,w3; 407*dd9fa9a5SLois Curfman McInnes int (*f)(void *,Vec,Vec,void *) = coloring->f; 408005c665bSBarry Smith void *fctx = coloring->fctx; 409005c665bSBarry Smith 410d4bb536fSBarry Smith /* 411005c665bSBarry Smith ierr = SNESGetIterationNumber((SNES)sctx,&it); CHKERRQ(ierr); 412005c665bSBarry Smith if ((freq > 1) && ((it % freq) != 1)) { 413005c665bSBarry Smith PLogInfo(coloring,"MatFDColoringApply:Skipping Jacobian, iteration %d, freq %d\n",it,freq); 414005c665bSBarry Smith *flag = SAME_PRECONDITIONER; 415005c665bSBarry Smith return 0; 416005c665bSBarry Smith } else { 417005c665bSBarry Smith PLogInfo(coloring,"MatFDColoringApply:Computing Jacobian, iteration %d, freq %d\n",it,freq); 418005c665bSBarry Smith *flag = SAME_NONZERO_PATTERN; 419d4bb536fSBarry Smith }*/ 420005c665bSBarry Smith 421005c665bSBarry Smith if (!coloring->w1) { 422005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1); CHKERRQ(ierr); 423005c665bSBarry Smith PLogObjectParent(coloring,coloring->w1); 424005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2); CHKERRQ(ierr); 425005c665bSBarry Smith PLogObjectParent(coloring,coloring->w2); 426005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3); CHKERRQ(ierr); 427005c665bSBarry Smith PLogObjectParent(coloring,coloring->w3); 428005c665bSBarry Smith } 429005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 43043a90d84SBarry Smith 43143a90d84SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 43243a90d84SBarry Smith 43343a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 43443a90d84SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 43543a90d84SBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 43643a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 43743a90d84SBarry Smith 43843a90d84SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 43943a90d84SBarry Smith /* 44043a90d84SBarry Smith Loop over each color 44143a90d84SBarry Smith */ 44243a90d84SBarry Smith 44343a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 44443a90d84SBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 44543a90d84SBarry Smith /* 44643a90d84SBarry Smith Loop over each column associated with color adding the 44743a90d84SBarry Smith perturbation to the vector w3. 44843a90d84SBarry Smith */ 44943a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 45043a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 45143a90d84SBarry Smith dx = xx[col-start]; 45243a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 45343a90d84SBarry Smith if (dx < umin && dx >= 0.0) dx = .1; 45443a90d84SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -.1; 45543a90d84SBarry Smith #else 45643a90d84SBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = .1; 45743a90d84SBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1; 45843a90d84SBarry Smith #endif 45943a90d84SBarry Smith dx *= epsilon; 46043a90d84SBarry Smith wscale[col] = 1.0/dx; 46143a90d84SBarry Smith VecSetValues(w3,1,&col,&dx,ADD_VALUES); 46243a90d84SBarry Smith } 46343a90d84SBarry Smith VecRestoreArray(x1,&xx); 46443a90d84SBarry Smith /* 46543a90d84SBarry Smith Evaluate function at x1 + dx (here dx is a vector, of perturbations) 46643a90d84SBarry Smith */ 46743a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 46843a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 46943a90d84SBarry Smith /* Communicate scale to all processors */ 47043a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 47143a90d84SBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 47243a90d84SBarry Smith #else 47343a90d84SBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 47443a90d84SBarry Smith #endif 47543a90d84SBarry Smith /* 47643a90d84SBarry Smith Loop over rows of vector putting results into Jacobian matrix 47743a90d84SBarry Smith */ 47843a90d84SBarry Smith VecGetArray(w2,&y); 47943a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 48043a90d84SBarry Smith row = coloring->rows[k][l]; 48143a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 48243a90d84SBarry Smith y[row] *= scale[col]; 48343a90d84SBarry Smith srow = row + start; 48443a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 48543a90d84SBarry Smith } 48643a90d84SBarry Smith VecRestoreArray(w2,&y); 48743a90d84SBarry Smith } 48843a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 48943a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 49043a90d84SBarry Smith return 0; 49143a90d84SBarry Smith } 492