1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*b4fc646aSLois Curfman McInnes static char vcid[] = "$Id: fdmatrix.c,v 1.15 1997/09/19 19:26:53 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 #undef __FUNC__ 73d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringView" 74bbf0e169SBarry Smith /*@C 75639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 76bbf0e169SBarry Smith 77*b4fc646aSLois Curfman McInnes Input Parameters: 78*b4fc646aSLois Curfman McInnes . c - the coloring context 79*b4fc646aSLois Curfman McInnes . viewer - visualization context 80*b4fc646aSLois Curfman McInnes 81*b4fc646aSLois Curfman McInnes Notes: 82*b4fc646aSLois Curfman McInnes The available visualization contexts include 83*b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_SELF - standard output (default) 84*b4fc646aSLois Curfman McInnes $ VIEWER_STDOUT_WORLD - synchronized standard 85*b4fc646aSLois Curfman McInnes $ output where only the first processor opens 86*b4fc646aSLois Curfman McInnes $ the file. All other processors send their 87*b4fc646aSLois Curfman McInnes $ data to the first processor to print. 88*b4fc646aSLois Curfman McInnes $ VIEWER_DRAWX_WORLD - graphical display of nonzero structure 89bbf0e169SBarry Smith 90639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 91005c665bSBarry Smith 92*b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 93bbf0e169SBarry Smith @*/ 94*b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer) 95bbf0e169SBarry Smith { 96005c665bSBarry Smith ViewerType vtype; 97639f9d9dSBarry Smith int i,j,format,ierr; 98*b4fc646aSLois Curfman McInnes FILE *fd; 99bbf0e169SBarry Smith 100*b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 101*b4fc646aSLois Curfman McInnes if (viewer) {PetscValidHeader(viewer);} 102*b4fc646aSLois Curfman McInnes else {viewer = VIEWER_STDOUT_SELF;} 103bbf0e169SBarry Smith 104005c665bSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 105005c665bSBarry Smith if (vtype == DRAW_VIEWER) { 106*b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer); CHKERRQ(ierr); 107005c665bSBarry Smith return 0; 108005c665bSBarry Smith } 109*b4fc646aSLois Curfman McInnes else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 110*b4fc646aSLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 111*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd,"MatFDColoring Object:\n"); 112*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Error tolerance=%g\n",c->error_rel); 113*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Umin=%g\n",c->umin); 114*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," 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) { 118*b4fc646aSLois Curfman McInnes for ( i=0; i<c->ncolors; i++ ) { 119*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Information for color %d\n",i); 120*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of columns %d\n",c->ncolumns[i]); 121*b4fc646aSLois Curfman McInnes for ( j=0; j<c->ncolumns[i]; j++ ) { 122*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d\n",c->columns[i][j]); 123639f9d9dSBarry Smith } 124*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," Number of rows %d\n",c->nrows[i]); 125*b4fc646aSLois Curfman McInnes for ( j=0; j<c->nrows[i]; j++ ) { 126*b4fc646aSLois Curfman McInnes PetscFPrintf(c->comm,fd," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]); 127*b4fc646aSLois Curfman McInnes } 128bbf0e169SBarry Smith } 129bbf0e169SBarry Smith } 130bbf0e169SBarry Smith } 131639f9d9dSBarry Smith return 0; 132639f9d9dSBarry Smith } 133639f9d9dSBarry Smith 1345615d1e5SSatish Balay #undef __FUNC__ 1355615d1e5SSatish Balay #define __FUNC__ "MatFDColoringSetParameters" 136639f9d9dSBarry Smith /*@ 137*b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 138*b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 139639f9d9dSBarry Smith 140639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 141639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 142639f9d9dSBarry Smith $ = error_rel*umin else 143639f9d9dSBarry Smith $ 144639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 145639f9d9dSBarry Smith 146639f9d9dSBarry Smith Input Parameters: 147*b4fc646aSLois Curfman McInnes . coloring - the coloring context 148639f9d9dSBarry Smith . error_rel - relative error 149639f9d9dSBarry Smith . umin - minimum allowable u-value 150639f9d9dSBarry Smith 151*b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 152*b4fc646aSLois Curfman McInnes 153*b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 154639f9d9dSBarry Smith @*/ 155639f9d9dSBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,double error,double umin) 156639f9d9dSBarry Smith { 157639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 158639f9d9dSBarry Smith 159639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 160639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 161639f9d9dSBarry Smith return 0; 162639f9d9dSBarry Smith } 163639f9d9dSBarry Smith 1645615d1e5SSatish Balay #undef __FUNC__ 165005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 166005c665bSBarry Smith /*@ 167*b4fc646aSLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency at which new Jacobian 168*b4fc646aSLois Curfman McInnes matrices are computed. 169005c665bSBarry Smith 170005c665bSBarry Smith Input Parameters: 171*b4fc646aSLois Curfman McInnes . coloring - the coloring context 172005c665bSBarry Smith . freq - frequency 173005c665bSBarry Smith 174*b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 175005c665bSBarry Smith @*/ 176005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 177005c665bSBarry Smith { 178005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 179005c665bSBarry Smith 180005c665bSBarry Smith matfd->freq = freq; 181005c665bSBarry Smith return 0; 182005c665bSBarry Smith } 183005c665bSBarry Smith 184005c665bSBarry Smith #undef __FUNC__ 185005c665bSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 186005c665bSBarry Smith /*@ 187005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 188005c665bSBarry Smith 189005c665bSBarry Smith Input Parameters: 190*b4fc646aSLois Curfman McInnes . coloring - the coloring context 191005c665bSBarry Smith . f - the function 192005c665bSBarry Smith . fctx - the function context 193005c665bSBarry Smith 194*b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 195005c665bSBarry Smith @*/ 196005c665bSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void *,Vec,Vec,void *),void *fctx) 197005c665bSBarry Smith { 198005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 199005c665bSBarry Smith 200005c665bSBarry Smith matfd->f = f; 201005c665bSBarry Smith matfd->fctx = fctx; 202005c665bSBarry Smith 203005c665bSBarry Smith return 0; 204005c665bSBarry Smith } 205005c665bSBarry Smith 206005c665bSBarry Smith #undef __FUNC__ 207d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 208639f9d9dSBarry Smith /*@ 209*b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 210639f9d9dSBarry Smith the options database. 211639f9d9dSBarry Smith 212*b4fc646aSLois Curfman McInnes The Jacobian is estimated with the differencing approximation 213639f9d9dSBarry Smith $ J(u)_{:,i} = [J(u+h*dx_{i}) - J(u)]/h where 214639f9d9dSBarry Smith $ h = error_rel*u[i] if u[i] > umin 215ae09f205SBarry Smith $ = error_rel*umin else 216639f9d9dSBarry Smith $ 217639f9d9dSBarry Smith $ dx_{i} = (0, ... 1, .... 0) 218639f9d9dSBarry Smith 219639f9d9dSBarry Smith Input Parameters: 220*b4fc646aSLois Curfman McInnes . coloring - the coloring context 221639f9d9dSBarry Smith 222*b4fc646aSLois Curfman McInnes Options Database Keys: 223*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_error <err>, where <err> is the square root 224*b4fc646aSLois Curfman McInnes $ of relative error in the function 225*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_freq <freq> where <freq> is the frequency at 226*b4fc646aSLois Curfman McInnes $ which Jacobian is computed 227*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_umin <umin>, where umin is described above 228639f9d9dSBarry Smith 229*b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 230639f9d9dSBarry Smith @*/ 231639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 232639f9d9dSBarry Smith { 233005c665bSBarry Smith int ierr,flag,freq = 1; 234639f9d9dSBarry Smith double error = PETSC_DEFAULT,umin = PETSC_DEFAULT; 235639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 236639f9d9dSBarry Smith 237639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,&flag);CHKERRQ(ierr); 238639f9d9dSBarry Smith ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,&flag);CHKERRQ(ierr); 239639f9d9dSBarry Smith ierr = MatFDColoringSetParameters(matfd,error,umin); CHKERRQ(ierr); 240005c665bSBarry Smith ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,&flag);CHKERRQ(ierr); 241005c665bSBarry Smith ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr); 242005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flag); CHKERRQ(ierr); 243639f9d9dSBarry Smith if (flag) { 244639f9d9dSBarry Smith ierr = MatFDColoringPrintHelp(matfd); CHKERRQ(ierr); 245639f9d9dSBarry Smith } 246639f9d9dSBarry Smith return 0; 247639f9d9dSBarry Smith } 248639f9d9dSBarry Smith 2495615d1e5SSatish Balay #undef __FUNC__ 250d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringPrintHelp" 251639f9d9dSBarry Smith /*@ 252639f9d9dSBarry Smith MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations 253639f9d9dSBarry Smith using coloring. 254639f9d9dSBarry Smith 255639f9d9dSBarry Smith Input Parameter: 256639f9d9dSBarry Smith . fdcoloring - the MatFDColoring context 257639f9d9dSBarry Smith 258639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions() 259639f9d9dSBarry Smith @*/ 260639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd) 261639f9d9dSBarry Smith { 262639f9d9dSBarry Smith PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE); 263639f9d9dSBarry Smith 2640f665d81SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to 1.e-8\n"); 2650f665d81SLois Curfman McInnes PetscPrintf(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to 1.e-8\n"); 266005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults 1\n"); 267005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view\n"); 268005c665bSBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_draw\n"); 269ae09f205SBarry Smith PetscPrintf(fd->comm,"-mat_fd_coloring_view_info\n"); 270005c665bSBarry Smith return 0; 271005c665bSBarry Smith } 272005c665bSBarry Smith 273005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 274005c665bSBarry Smith { 275005c665bSBarry Smith int ierr,flg; 276005c665bSBarry Smith 277005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg); CHKERRQ(ierr); 278005c665bSBarry Smith if (flg) { 279005c665bSBarry Smith Viewer viewer; 280005c665bSBarry Smith ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr); 281005c665bSBarry Smith ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr); 282005c665bSBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 283005c665bSBarry Smith } 284ae09f205SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg); CHKERRQ(ierr); 285ae09f205SBarry Smith if (flg) { 286ae09f205SBarry Smith Viewer viewer; 287ae09f205SBarry Smith ierr = ViewerFileOpenASCII(fd->comm,"stdout",&viewer);CHKERRQ(ierr); 288ae09f205SBarry Smith ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 289ae09f205SBarry Smith ierr = MatFDColoringView(fd,viewer); CHKERRQ(ierr); 290ae09f205SBarry Smith ierr = ViewerPopFormat(viewer);CHKERRQ(ierr); 291ae09f205SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 292ae09f205SBarry Smith } 293005c665bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg); CHKERRQ(ierr); 294005c665bSBarry Smith if (flg) { 295005c665bSBarry Smith ierr = MatFDColoringView(fd,VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 296005c665bSBarry Smith ierr = ViewerFlush(VIEWER_DRAWX_(fd->comm)); CHKERRQ(ierr); 297005c665bSBarry Smith } 298bbf0e169SBarry Smith return 0; 299bbf0e169SBarry Smith } 300bbf0e169SBarry Smith 3015615d1e5SSatish Balay #undef __FUNC__ 3025615d1e5SSatish Balay #define __FUNC__ "MatFDColoringCreate" 303bbf0e169SBarry Smith /*@C 304639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 305639f9d9dSBarry Smith computation of Jacobians. 306bbf0e169SBarry Smith 307639f9d9dSBarry Smith Input Parameters: 308639f9d9dSBarry Smith . mat - the matrix containing the nonzero structure of the Jacobian 309639f9d9dSBarry Smith . iscoloring - the coloring of the matrix 310bbf0e169SBarry Smith 311bbf0e169SBarry Smith Output Parameter: 312639f9d9dSBarry Smith . color - the new coloring context 313bbf0e169SBarry Smith 314*b4fc646aSLois Curfman McInnes Options Database Keys: 315*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view 316*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_draw 317*b4fc646aSLois Curfman McInnes $ -mat_fd_coloring_view_info 318639f9d9dSBarry Smith 319639f9d9dSBarry Smith .seealso: MatFDColoringDestroy() 320bbf0e169SBarry Smith @*/ 321639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 322bbf0e169SBarry Smith { 323639f9d9dSBarry Smith MatFDColoring c; 324639f9d9dSBarry Smith MPI_Comm comm; 325639f9d9dSBarry Smith int ierr,M,N; 326639f9d9dSBarry Smith 327639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N); CHKERRQ(ierr); 328e3372554SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices"); 329639f9d9dSBarry Smith 330639f9d9dSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 331d4bb536fSBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,MAT_FDCOLORING_COOKIE,0,comm,MatFDColoringDestroy,MatFDColoringView); 332639f9d9dSBarry Smith PLogObjectCreate(c); 333639f9d9dSBarry Smith 334639f9d9dSBarry Smith if (mat->ops.fdcoloringcreate) { 335639f9d9dSBarry Smith ierr = (*mat->ops.fdcoloringcreate)(mat,iscoloring,c); CHKERRQ(ierr); 336639f9d9dSBarry Smith } else { 337e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type"); 338639f9d9dSBarry Smith } 339639f9d9dSBarry Smith 340639f9d9dSBarry Smith c->error_rel = 1.e-8; 341ae09f205SBarry Smith c->umin = 1.e-6; 342005c665bSBarry Smith c->freq = 1; 343005c665bSBarry Smith 344005c665bSBarry Smith ierr = MatFDColoringView_Private(c); CHKERRQ(ierr); 345639f9d9dSBarry Smith 346639f9d9dSBarry Smith *color = c; 347639f9d9dSBarry Smith 348bbf0e169SBarry Smith return 0; 349bbf0e169SBarry Smith } 350bbf0e169SBarry Smith 3515615d1e5SSatish Balay #undef __FUNC__ 352d4bb536fSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 353bbf0e169SBarry Smith /*@C 354639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 355639f9d9dSBarry Smith via MatFDColoringCreate(). 356bbf0e169SBarry Smith 357*b4fc646aSLois Curfman McInnes Input Parameter: 358639f9d9dSBarry Smith . c - coloring context 359bbf0e169SBarry Smith 360639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 361bbf0e169SBarry Smith @*/ 362639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 363bbf0e169SBarry Smith { 364639f9d9dSBarry Smith int i,ierr,flag; 365bbf0e169SBarry Smith 366d4bb536fSBarry Smith if (--c->refct > 0) return 0; 367d4bb536fSBarry Smith 368639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view",&flag); 369639f9d9dSBarry Smith if (flag) { 370bbf0e169SBarry Smith Viewer viewer; 371639f9d9dSBarry Smith ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr); 372639f9d9dSBarry Smith ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr); 373bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 374bbf0e169SBarry Smith } 375639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matfdcoloring_view_info",&flag); 376639f9d9dSBarry Smith if (flag) { 377bbf0e169SBarry Smith Viewer viewer; 378639f9d9dSBarry Smith ierr = ViewerFileOpenASCII(c->comm,"stdout",&viewer);CHKERRQ(ierr); 379639f9d9dSBarry Smith ierr = ViewerPushFormat(viewer,VIEWER_FORMAT_ASCII_INFO, PETSC_NULL); CHKERRQ(ierr); 380639f9d9dSBarry Smith ierr = MatFDColoringView(c,viewer);CHKERRQ(ierr); 381639f9d9dSBarry Smith ierr = ViewerPopFormat(viewer); 382bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 383bbf0e169SBarry Smith } 384639f9d9dSBarry Smith 385639f9d9dSBarry Smith for ( i=0; i<c->ncolors; i++ ) { 386639f9d9dSBarry Smith if (c->columns[i]) PetscFree(c->columns[i]); 387639f9d9dSBarry Smith if (c->rows[i]) PetscFree(c->rows[i]); 388639f9d9dSBarry Smith if (c->columnsforrow[i]) PetscFree(c->columnsforrow[i]); 389bbf0e169SBarry Smith } 390639f9d9dSBarry Smith PetscFree(c->ncolumns); 391639f9d9dSBarry Smith PetscFree(c->columns); 392639f9d9dSBarry Smith PetscFree(c->nrows); 393639f9d9dSBarry Smith PetscFree(c->rows); 394639f9d9dSBarry Smith PetscFree(c->columnsforrow); 395639f9d9dSBarry Smith PetscFree(c->scale); 396005c665bSBarry Smith if (c->w1) { 397005c665bSBarry Smith ierr = VecDestroy(c->w1); CHKERRQ(ierr); 398005c665bSBarry Smith ierr = VecDestroy(c->w2); CHKERRQ(ierr); 399005c665bSBarry Smith ierr = VecDestroy(c->w3); CHKERRQ(ierr); 400005c665bSBarry Smith } 401639f9d9dSBarry Smith PLogObjectDestroy(c); 402639f9d9dSBarry Smith PetscHeaderDestroy(c); 403bbf0e169SBarry Smith return 0; 404bbf0e169SBarry Smith } 40543a90d84SBarry Smith 406005c665bSBarry Smith #include "snes.h" 407005c665bSBarry Smith 4085615d1e5SSatish Balay #undef __FUNC__ 4095615d1e5SSatish Balay #define __FUNC__ "MatFDColoringApply" 41043a90d84SBarry Smith /*@ 41143a90d84SBarry Smith MatFDColoringApply - Given a matrix for which a MatFDColoring has been created, 41243a90d84SBarry Smith computes the Jacobian for a function via finite differences. 41343a90d84SBarry Smith 41443a90d84SBarry Smith Input Parameters: 41543a90d84SBarry Smith . mat - location to store Jacobian 41643a90d84SBarry Smith . coloring - coloring context created with MatFDColoringCreate() 41743a90d84SBarry Smith . x1 - location at which Jacobian is to be computed 418005c665bSBarry Smith . sctx - optional context required by function (actually a SNES context) 41943a90d84SBarry Smith 42043a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 42143a90d84SBarry Smith 42243a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 42343a90d84SBarry Smith @*/ 424005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 42543a90d84SBarry Smith { 426dd9fa9a5SLois Curfman McInnes int k, ierr,N,start,end,l,row,col,srow; 42743a90d84SBarry Smith Scalar dx, mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale; 42843a90d84SBarry Smith double epsilon = coloring->error_rel, umin = coloring->umin; 42943a90d84SBarry Smith MPI_Comm comm = coloring->comm; 430005c665bSBarry Smith Vec w1,w2,w3; 431dd9fa9a5SLois Curfman McInnes int (*f)(void *,Vec,Vec,void *) = coloring->f; 432005c665bSBarry Smith void *fctx = coloring->fctx; 433005c665bSBarry Smith 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 45543a90d84SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 45643a90d84SBarry Smith 45743a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 45843a90d84SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 45943a90d84SBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 46043a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx); CHKERRQ(ierr); 46143a90d84SBarry Smith 46243a90d84SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 46343a90d84SBarry Smith /* 46443a90d84SBarry Smith Loop over each color 46543a90d84SBarry Smith */ 46643a90d84SBarry Smith 46743a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 46843a90d84SBarry Smith ierr = VecCopy(x1,w3); CHKERRQ(ierr); 46943a90d84SBarry Smith /* 47043a90d84SBarry Smith Loop over each column associated with color adding the 47143a90d84SBarry Smith perturbation to the vector w3. 47243a90d84SBarry Smith */ 47343a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 47443a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 47543a90d84SBarry Smith dx = xx[col-start]; 476ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 47743a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 478ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 479ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 48043a90d84SBarry Smith #else 481ae09f205SBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = umin; 482ae09f205SBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -umin; 48343a90d84SBarry Smith #endif 48443a90d84SBarry Smith dx *= epsilon; 48543a90d84SBarry Smith wscale[col] = 1.0/dx; 48643a90d84SBarry Smith VecSetValues(w3,1,&col,&dx,ADD_VALUES); 48743a90d84SBarry Smith } 48843a90d84SBarry Smith VecRestoreArray(x1,&xx); 48943a90d84SBarry Smith /* 49043a90d84SBarry Smith Evaluate function at x1 + dx (here dx is a vector, of perturbations) 49143a90d84SBarry Smith */ 49243a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx); CHKERRQ(ierr); 49343a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2); CHKERRQ(ierr); 49443a90d84SBarry Smith /* Communicate scale to all processors */ 49543a90d84SBarry Smith #if !defined(PETSC_COMPLEX) 49643a90d84SBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 49743a90d84SBarry Smith #else 49843a90d84SBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 49943a90d84SBarry Smith #endif 50043a90d84SBarry Smith /* 50143a90d84SBarry Smith Loop over rows of vector putting results into Jacobian matrix 50243a90d84SBarry Smith */ 50343a90d84SBarry Smith VecGetArray(w2,&y); 50443a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 50543a90d84SBarry Smith row = coloring->rows[k][l]; 50643a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 50743a90d84SBarry Smith y[row] *= scale[col]; 50843a90d84SBarry Smith srow = row + start; 50943a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 51043a90d84SBarry Smith } 51143a90d84SBarry Smith VecRestoreArray(w2,&y); 51243a90d84SBarry Smith } 51343a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 51443a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 51543a90d84SBarry Smith return 0; 51643a90d84SBarry Smith } 517