1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3bbf0e169SBarry Smith /* 4639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 5639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 6bbf0e169SBarry Smith */ 7bbf0e169SBarry Smith 8e090d566SSatish Balay #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 9bbf0e169SBarry Smith 108ba1e511SMatthew Knepley /* Logging support */ 11be1d678aSKris Buschelman PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE = 0; 128ba1e511SMatthew Knepley 134a2ae208SSatish Balay #undef __FUNCT__ 143a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF" 15be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F) 163a7fca6bSBarry Smith { 173a7fca6bSBarry Smith PetscFunctionBegin; 183a7fca6bSBarry Smith fd->F = F; 193a7fca6bSBarry Smith PetscFunctionReturn(0); 203a7fca6bSBarry Smith } 213a7fca6bSBarry Smith 223a7fca6bSBarry Smith #undef __FUNCT__ 234a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 246849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 259194eea9SBarry Smith { 269194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 27dfbe8321SBarry Smith PetscErrorCode ierr; 2838baddfdSBarry Smith PetscInt i,j; 299194eea9SBarry Smith PetscReal x,y; 309194eea9SBarry Smith 319194eea9SBarry Smith PetscFunctionBegin; 329194eea9SBarry Smith 339194eea9SBarry Smith /* loop over colors */ 349194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 359194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 369194eea9SBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 379194eea9SBarry Smith x = fd->columnsforrow[i][j]; 38b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 399194eea9SBarry Smith } 409194eea9SBarry Smith } 419194eea9SBarry Smith PetscFunctionReturn(0); 429194eea9SBarry Smith } 439194eea9SBarry Smith 444a2ae208SSatish Balay #undef __FUNCT__ 454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 466849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 47005c665bSBarry Smith { 48dfbe8321SBarry Smith PetscErrorCode ierr; 49005c665bSBarry Smith PetscTruth isnull; 50b0a32e0cSBarry Smith PetscDraw draw; 5154d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 52005c665bSBarry Smith 533a40ed3dSBarry Smith PetscFunctionBegin; 54b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 55b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 569194eea9SBarry Smith 579194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 58005c665bSBarry Smith 59005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 60005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 61b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 62b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 639194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 65005c665bSBarry Smith } 66005c665bSBarry Smith 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 69bbf0e169SBarry Smith /*@C 70639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 71bbf0e169SBarry Smith 724c49b128SBarry Smith Collective on MatFDColoring 73fee21e36SBarry Smith 74ef5ee4d1SLois Curfman McInnes Input Parameters: 75ef5ee4d1SLois Curfman McInnes + c - the coloring context 76ef5ee4d1SLois Curfman McInnes - viewer - visualization context 77ef5ee4d1SLois Curfman McInnes 7815091d37SBarry Smith Level: intermediate 7915091d37SBarry Smith 80b4fc646aSLois Curfman McInnes Notes: 81b4fc646aSLois Curfman McInnes The available visualization contexts include 82b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 83b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 84ef5ee4d1SLois Curfman McInnes output where only the first processor opens 85ef5ee4d1SLois Curfman McInnes the file. All other processors send their 86ef5ee4d1SLois Curfman McInnes data to the first processor to print. 87b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 88bbf0e169SBarry Smith 899194eea9SBarry Smith Notes: 909194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 91b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 929194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 939194eea9SBarry Smith 94639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 95005c665bSBarry Smith 96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 97bbf0e169SBarry Smith @*/ 98be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer) 99bbf0e169SBarry Smith { 1006849ba73SBarry Smith PetscErrorCode ierr; 10138baddfdSBarry Smith PetscInt i,j; 10232077d6dSBarry Smith PetscTruth isdraw,iascii; 103f3ef73ceSBarry Smith PetscViewerFormat format; 104bbf0e169SBarry Smith 1053a40ed3dSBarry Smith PetscFunctionBegin; 1064482741eSBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1); 107b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 1084482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 109c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 110bbf0e169SBarry Smith 111fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 11232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1130f5bd95cSBarry Smith if (isdraw) { 114b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 11532077d6dSBarry Smith } else if (iascii) { 116b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 117a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 118a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 11977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 120ae09f205SBarry Smith 121b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 122fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 123b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 12477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 12577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 126b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 12777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 128639f9d9dSBarry Smith } 12977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 130b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 13177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 132b4fc646aSLois Curfman McInnes } 133bbf0e169SBarry Smith } 134bbf0e169SBarry Smith } 135b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1365cd90555SBarry Smith } else { 13779a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 138bbf0e169SBarry Smith } 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 140639f9d9dSBarry Smith } 141639f9d9dSBarry Smith 1424a2ae208SSatish Balay #undef __FUNCT__ 1434a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 144639f9d9dSBarry Smith /*@ 145b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 146b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 147639f9d9dSBarry Smith 148ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 149ef5ee4d1SLois Curfman McInnes 150ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 151ef5ee4d1SLois Curfman McInnes .vb 15265f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 153f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 154f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 155ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 156ef5ee4d1SLois Curfman McInnes .ve 157639f9d9dSBarry Smith 158639f9d9dSBarry Smith Input Parameters: 159ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 160639f9d9dSBarry Smith . error_rel - relative error 161f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 162fee21e36SBarry Smith 16315091d37SBarry Smith Level: advanced 16415091d37SBarry Smith 165b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 166b4fc646aSLois Curfman McInnes 167b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 168639f9d9dSBarry Smith @*/ 169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 170639f9d9dSBarry Smith { 1713a40ed3dSBarry Smith PetscFunctionBegin; 1724482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 173639f9d9dSBarry Smith 174639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 175639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1763a40ed3dSBarry Smith PetscFunctionReturn(0); 177639f9d9dSBarry Smith } 178639f9d9dSBarry Smith 1794a2ae208SSatish Balay #undef __FUNCT__ 1804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency" 181005c665bSBarry Smith /*@ 182e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 183e0907662SLois Curfman McInnes matrices. 184005c665bSBarry Smith 185fee21e36SBarry Smith Collective on MatFDColoring 186fee21e36SBarry Smith 187ef5ee4d1SLois Curfman McInnes Input Parameters: 188ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 189ef5ee4d1SLois Curfman McInnes - freq - frequency (default is 1) 190ef5ee4d1SLois Curfman McInnes 19115091d37SBarry Smith Options Database Keys: 19215091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 19315091d37SBarry Smith 19415091d37SBarry Smith Level: advanced 19515091d37SBarry Smith 196e0907662SLois Curfman McInnes Notes: 197e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 198e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 199e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 200e0907662SLois Curfman McInnes <freq> nonlinear iterations. 201e0907662SLois Curfman McInnes 202b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 203ef5ee4d1SLois Curfman McInnes 20440964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute() 205005c665bSBarry Smith @*/ 206be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring matfd,PetscInt freq) 207005c665bSBarry Smith { 2083a40ed3dSBarry Smith PetscFunctionBegin; 2094482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 210005c665bSBarry Smith 211005c665bSBarry Smith matfd->freq = freq; 2123a40ed3dSBarry Smith PetscFunctionReturn(0); 213005c665bSBarry Smith } 214005c665bSBarry Smith 2154a2ae208SSatish Balay #undef __FUNCT__ 2164a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency" 217ff0cfa39SBarry Smith /*@ 218ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 219ff0cfa39SBarry Smith matrices. 220ff0cfa39SBarry Smith 221ef5ee4d1SLois Curfman McInnes Not Collective 222ef5ee4d1SLois Curfman McInnes 223ff0cfa39SBarry Smith Input Parameters: 224ff0cfa39SBarry Smith . coloring - the coloring context 225ff0cfa39SBarry Smith 226ff0cfa39SBarry Smith Output Parameters: 227ff0cfa39SBarry Smith . freq - frequency (default is 1) 228ff0cfa39SBarry Smith 22915091d37SBarry Smith Options Database Keys: 23015091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 23115091d37SBarry Smith 23215091d37SBarry Smith Level: advanced 23315091d37SBarry Smith 234ff0cfa39SBarry Smith Notes: 235ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 236ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 237ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 238ff0cfa39SBarry Smith <freq> nonlinear iterations. 239ff0cfa39SBarry Smith 240ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 241ef5ee4d1SLois Curfman McInnes 242ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency() 243ff0cfa39SBarry Smith @*/ 244be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring matfd,PetscInt *freq) 245ff0cfa39SBarry Smith { 2463a40ed3dSBarry Smith PetscFunctionBegin; 2474482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 248ff0cfa39SBarry Smith *freq = matfd->freq; 2493a40ed3dSBarry Smith PetscFunctionReturn(0); 250ff0cfa39SBarry Smith } 251ff0cfa39SBarry Smith 2524a2ae208SSatish Balay #undef __FUNCT__ 2534a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction" 2544a9d489dSBarry Smith /*@C 2554a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2564a9d489dSBarry Smith 2574a9d489dSBarry Smith Collective on MatFDColoring 2584a9d489dSBarry Smith 2594a9d489dSBarry Smith Input Parameters: 2604a9d489dSBarry Smith . coloring - the coloring context 2614a9d489dSBarry Smith 2624a9d489dSBarry Smith Output Parameters: 2634a9d489dSBarry Smith + f - the function 2644a9d489dSBarry Smith - fctx - the optional user-defined function context 2654a9d489dSBarry Smith 2664a9d489dSBarry Smith Level: intermediate 2674a9d489dSBarry Smith 2684a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 2694a9d489dSBarry Smith @*/ 2704a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2714a9d489dSBarry Smith { 2724a9d489dSBarry Smith PetscFunctionBegin; 2734a9d489dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 2744a9d489dSBarry Smith if (f) *f = matfd->f; 2754a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2764a9d489dSBarry Smith PetscFunctionReturn(0); 2774a9d489dSBarry Smith } 2784a9d489dSBarry Smith 2794a9d489dSBarry Smith #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 281d64ed03dSBarry Smith /*@C 282005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 283005c665bSBarry Smith 284fee21e36SBarry Smith Collective on MatFDColoring 285fee21e36SBarry Smith 286ef5ee4d1SLois Curfman McInnes Input Parameters: 287ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 288ef5ee4d1SLois Curfman McInnes . f - the function 289ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 290ef5ee4d1SLois Curfman McInnes 29115091d37SBarry Smith Level: intermediate 29215091d37SBarry Smith 293f881d145SBarry Smith Notes: 294f881d145SBarry Smith In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 295f881d145SBarry Smith be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 296f881d145SBarry Smith with the TS solvers. 297f881d145SBarry Smith 298b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 299005c665bSBarry Smith @*/ 300be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 301005c665bSBarry Smith { 3023a40ed3dSBarry Smith PetscFunctionBegin; 3034482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 304005c665bSBarry Smith matfd->f = f; 305005c665bSBarry Smith matfd->fctx = fctx; 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 307005c665bSBarry Smith } 308005c665bSBarry Smith 3094a2ae208SSatish Balay #undef __FUNCT__ 3104a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 311639f9d9dSBarry Smith /*@ 312b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 313639f9d9dSBarry Smith the options database. 314639f9d9dSBarry Smith 315fee21e36SBarry Smith Collective on MatFDColoring 316fee21e36SBarry Smith 31765f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 318ef5ee4d1SLois Curfman McInnes .vb 31965f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 320f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 321f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 322ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 323ef5ee4d1SLois Curfman McInnes .ve 324ef5ee4d1SLois Curfman McInnes 325ef5ee4d1SLois Curfman McInnes Input Parameter: 326ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 327ef5ee4d1SLois Curfman McInnes 328b4fc646aSLois Curfman McInnes Options Database Keys: 329d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 330ef5ee4d1SLois Curfman McInnes of relative error in the function) 331f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 332ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 333efb30889SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS) 334ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 335ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 336ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 337639f9d9dSBarry Smith 33815091d37SBarry Smith Level: intermediate 33915091d37SBarry Smith 340b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 341d1c39f55SBarry Smith 342d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 343d1c39f55SBarry Smith 344639f9d9dSBarry Smith @*/ 345be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd) 346639f9d9dSBarry Smith { 347dfbe8321SBarry Smith PetscErrorCode ierr; 348efb30889SBarry Smith PetscTruth flg; 349efb30889SBarry Smith char value[3]; 3503a40ed3dSBarry Smith 3513a40ed3dSBarry Smith PetscFunctionBegin; 3524482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 353639f9d9dSBarry Smith 354b0a32e0cSBarry Smith ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 35587828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 35687828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 357b0a32e0cSBarry Smith ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 358efb30889SBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr); 359efb30889SBarry Smith if (flg) { 360efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 361efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 362efb30889SBarry Smith else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 363efb30889SBarry Smith } 364186905e3SBarry Smith /* not used here; but so they are presented in the GUI */ 365b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 366b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 367b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 368b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3693a40ed3dSBarry Smith PetscFunctionReturn(0); 370005c665bSBarry Smith } 371005c665bSBarry Smith 3724a2ae208SSatish Balay #undef __FUNCT__ 3734a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private" 374dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd) 375005c665bSBarry Smith { 376dfbe8321SBarry Smith PetscErrorCode ierr; 377f1af5d2fSBarry Smith PetscTruth flg; 378005c665bSBarry Smith 3793a40ed3dSBarry Smith PetscFunctionBegin; 380b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 381005c665bSBarry Smith if (flg) { 382b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 383005c665bSBarry Smith } 384b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 385ae09f205SBarry Smith if (flg) { 386fb9695e5SSatish Balay ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 387b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 388b0a32e0cSBarry Smith ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 389ae09f205SBarry Smith } 390b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 391005c665bSBarry Smith if (flg) { 392b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 393b0a32e0cSBarry Smith ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 394005c665bSBarry Smith } 3953a40ed3dSBarry Smith PetscFunctionReturn(0); 396bbf0e169SBarry Smith } 397bbf0e169SBarry Smith 3984a2ae208SSatish Balay #undef __FUNCT__ 3994a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 40005869f15SSatish Balay /*@ 401639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 402639f9d9dSBarry Smith computation of Jacobians. 403bbf0e169SBarry Smith 404ef5ee4d1SLois Curfman McInnes Collective on Mat 405ef5ee4d1SLois Curfman McInnes 406639f9d9dSBarry Smith Input Parameters: 407ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 408ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 409bbf0e169SBarry Smith 410bbf0e169SBarry Smith Output Parameter: 411639f9d9dSBarry Smith . color - the new coloring context 412bbf0e169SBarry Smith 41315091d37SBarry Smith Level: intermediate 41415091d37SBarry Smith 4156831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 416d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 417d1c39f55SBarry Smith MatFDColoringSetFrequency(), MatFDColoringSetRecompute(), MatFDColoringView(), 418d1c39f55SBarry Smith MatFDColoringSetParameters() 419bbf0e169SBarry Smith @*/ 420be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 421bbf0e169SBarry Smith { 422639f9d9dSBarry Smith MatFDColoring c; 423639f9d9dSBarry Smith MPI_Comm comm; 424dfbe8321SBarry Smith PetscErrorCode ierr; 42538baddfdSBarry Smith PetscInt M,N; 426b8f8c88eSHong Zhang PetscMPIInt size; 427639f9d9dSBarry Smith 4283a40ed3dSBarry Smith PetscFunctionBegin; 429d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 430639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 43129bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 432639f9d9dSBarry Smith 433f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 43452e6d16bSBarry Smith ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 435639f9d9dSBarry Smith 436b8f8c88eSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 437b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 438b8f8c88eSHong Zhang if (size == 1) c->ctype = iscoloring->ctype = IS_COLORING_LOCAL; 439b8f8c88eSHong Zhang 440f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 441f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 442639f9d9dSBarry Smith } else { 44329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 444639f9d9dSBarry Smith } 445639f9d9dSBarry Smith 446b8f8c88eSHong Zhang ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr); 447b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr); 448b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 449b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr); 450b8f8c88eSHong Zhang 45177d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 45277d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 453005c665bSBarry Smith c->freq = 1; 45440964ad5SBarry Smith c->usersetsrecompute = PETSC_FALSE; 45540964ad5SBarry Smith c->recompute = PETSC_FALSE; 45649b058dcSBarry Smith c->currentcolor = -1; 457efb30889SBarry Smith c->htype = "wp"; 458639f9d9dSBarry Smith 459639f9d9dSBarry Smith *color = c; 460d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4613a40ed3dSBarry Smith PetscFunctionReturn(0); 462bbf0e169SBarry Smith } 463bbf0e169SBarry Smith 4644a2ae208SSatish Balay #undef __FUNCT__ 4654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 46605869f15SSatish Balay /*@ 467639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 468639f9d9dSBarry Smith via MatFDColoringCreate(). 469bbf0e169SBarry Smith 470ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 471ef5ee4d1SLois Curfman McInnes 472b4fc646aSLois Curfman McInnes Input Parameter: 473639f9d9dSBarry Smith . c - coloring context 474bbf0e169SBarry Smith 47515091d37SBarry Smith Level: intermediate 47615091d37SBarry Smith 477639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 478bbf0e169SBarry Smith @*/ 479be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c) 480bbf0e169SBarry Smith { 4816849ba73SBarry Smith PetscErrorCode ierr; 48238baddfdSBarry Smith PetscInt i; 483bbf0e169SBarry Smith 4843a40ed3dSBarry Smith PetscFunctionBegin; 4853a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 486d4bb536fSBarry Smith 487639f9d9dSBarry Smith for (i=0; i<c->ncolors; i++) { 48805b42c5fSBarry Smith ierr = PetscFree(c->columns[i]);CHKERRQ(ierr); 48905b42c5fSBarry Smith ierr = PetscFree(c->rows[i]);CHKERRQ(ierr); 49005b42c5fSBarry Smith ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr); 49105b42c5fSBarry Smith if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 492bbf0e169SBarry Smith } 493606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 494606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 495606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 496606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 497606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 49805b42c5fSBarry Smith ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr); 4994f113abfSBarry Smith if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 500005c665bSBarry Smith if (c->w1) { 501005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 502005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 503b8f8c88eSHong Zhang } 504b8f8c88eSHong Zhang if (c->w3){ 505005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 506005c665bSBarry Smith } 507d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 5083a40ed3dSBarry Smith PetscFunctionReturn(0); 509bbf0e169SBarry Smith } 51043a90d84SBarry Smith 5114a2ae208SSatish Balay #undef __FUNCT__ 51249b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 51349b058dcSBarry Smith /*@C 51449b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 51549b058dcSBarry Smith that are currently being perturbed. 51649b058dcSBarry Smith 51749b058dcSBarry Smith Not Collective 51849b058dcSBarry Smith 51949b058dcSBarry Smith Input Parameters: 52049b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 52149b058dcSBarry Smith 52249b058dcSBarry Smith Output Parameters: 52349b058dcSBarry Smith + n - the number of local columns being perturbed 52449b058dcSBarry Smith - cols - the column indices, in global numbering 52549b058dcSBarry Smith 52649b058dcSBarry Smith Level: intermediate 52749b058dcSBarry Smith 52849b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 52949b058dcSBarry Smith 53049b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 53149b058dcSBarry Smith @*/ 532be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 53349b058dcSBarry Smith { 53449b058dcSBarry Smith PetscFunctionBegin; 53549b058dcSBarry Smith if (coloring->currentcolor >= 0) { 53649b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 53749b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 53849b058dcSBarry Smith } else { 53949b058dcSBarry Smith *n = 0; 54049b058dcSBarry Smith } 54149b058dcSBarry Smith PetscFunctionReturn(0); 54249b058dcSBarry Smith } 54349b058dcSBarry Smith 544b8f8c88eSHong Zhang #include "petscda.h" /*I "petscda.h" I*/ 545b8f8c88eSHong Zhang 54649b058dcSBarry Smith #undef __FUNCT__ 5474a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 54843a90d84SBarry Smith /*@ 549e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 550e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 55143a90d84SBarry Smith 552fee21e36SBarry Smith Collective on MatFDColoring 553fee21e36SBarry Smith 554ef5ee4d1SLois Curfman McInnes Input Parameters: 555ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 556ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 557ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 558ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 559ef5ee4d1SLois Curfman McInnes 5608bba8e72SBarry Smith Options Database Keys: 561b9722508SLois Curfman McInnes + -mat_fd_coloring_freq <freq> - Sets coloring frequency 562efb30889SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATSNESMF_WP or MATSNESMF_DS) 563b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 564b9722508SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 565b9722508SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 5668bba8e72SBarry Smith 56715091d37SBarry Smith Level: intermediate 56815091d37SBarry Smith 56943a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 57043a90d84SBarry Smith 57143a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 57243a90d84SBarry Smith @*/ 573be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 57443a90d84SBarry Smith { 5756849ba73SBarry Smith PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; 5766849ba73SBarry Smith PetscErrorCode ierr; 577b8f8c88eSHong Zhang PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 578efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 579ea709b57SSatish Balay PetscScalar *vscale_array; 580efb30889SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 581b8f8c88eSHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 582005c665bSBarry Smith void *fctx = coloring->fctx; 583f1af5d2fSBarry Smith PetscTruth flg; 584*705d48f7SSatish Balay PetscInt ctype=coloring->ctype,N,col_start=0,col_end=0; 585b8f8c88eSHong Zhang Vec x1_tmp; 586a2e34c3dSBarry Smith 5873a40ed3dSBarry Smith PetscFunctionBegin; 5884482741eSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE,1); 5894482741eSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 5904482741eSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE,3); 591feba9168SBarry Smith if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 592e0907662SLois Curfman McInnes 59340964ad5SBarry Smith if (coloring->usersetsrecompute) { 59440964ad5SBarry Smith if (!coloring->recompute) { 59540964ad5SBarry Smith *flag = SAME_PRECONDITIONER; 596ae15b995SBarry Smith ierr = PetscInfo(J,"Skipping Jacobian, since user called MatFDColorSetRecompute()\n");CHKERRQ(ierr); 59740964ad5SBarry Smith PetscFunctionReturn(0); 59840964ad5SBarry Smith } else { 59940964ad5SBarry Smith coloring->recompute = PETSC_FALSE; 60040964ad5SBarry Smith } 60140964ad5SBarry Smith } 60240964ad5SBarry Smith 603d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 6044c49b128SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 605b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 606f1af5d2fSBarry Smith if (flg) { 607ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 608e0907662SLois Curfman McInnes } else { 6090b9b6f31SBarry Smith PetscTruth assembled; 6100b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 6110b9b6f31SBarry Smith if (assembled) { 61243a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 613e0907662SLois Curfman McInnes } 6140b9b6f31SBarry Smith } 61543a90d84SBarry Smith 616b8f8c88eSHong Zhang x1_tmp = x1; 617dac7f5fdSBarry Smith if (!coloring->vscale){ 618b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 619b8f8c88eSHong Zhang } 620be2fbe1fSBarry Smith 6213a7fca6bSBarry Smith /* 6223a7fca6bSBarry Smith This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets 6233a7fca6bSBarry Smith coloring->F for the coarser grids from the finest 6243a7fca6bSBarry Smith */ 6253a7fca6bSBarry Smith if (coloring->F) { 6263a7fca6bSBarry Smith ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 6273a7fca6bSBarry Smith ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 6283a7fca6bSBarry Smith if (m1 != m2) { 6293a7fca6bSBarry Smith coloring->F = 0; 6303a7fca6bSBarry Smith } 6313a7fca6bSBarry Smith } 6323a7fca6bSBarry Smith 633b8f8c88eSHong Zhang if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 634b8f8c88eSHong Zhang ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 635b8f8c88eSHong Zhang } 636b8f8c88eSHong Zhang ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 637b8f8c88eSHong Zhang 638b8f8c88eSHong Zhang /* Set w1 = F(x1) */ 639be2fbe1fSBarry Smith if (coloring->F) { 640be2fbe1fSBarry Smith w1 = coloring->F; /* use already computed value of function */ 641be2fbe1fSBarry Smith coloring->F = 0; 642be2fbe1fSBarry Smith } else { 64366f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 644b8f8c88eSHong Zhang ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 64566f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 646be2fbe1fSBarry Smith } 64743a90d84SBarry Smith 648b8f8c88eSHong Zhang if (!coloring->w3) { 649b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 650b8f8c88eSHong Zhang ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 651efb30889SBarry Smith } 652b8f8c88eSHong Zhang w3 = coloring->w3; 653efb30889SBarry Smith 654b8f8c88eSHong Zhang /* Compute all the local scale factors, including ghost points */ 655b8f8c88eSHong Zhang ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 656b8f8c88eSHong Zhang ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 657b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 658b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED){ 659b8f8c88eSHong Zhang col_start = 0; col_end = N; 660b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_LOCAL){ 661b8f8c88eSHong Zhang xx = xx - start; 662b8f8c88eSHong Zhang vscale_array = vscale_array - start; 663b8f8c88eSHong Zhang col_start = start; col_end = N + start; 664b8f8c88eSHong Zhang } 665b8f8c88eSHong Zhang for (col=col_start; col<col_end; col++){ 666b8f8c88eSHong Zhang /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 667efb30889SBarry Smith if (coloring->htype[0] == 'w') { 668efb30889SBarry Smith dx = 1.0 + unorm; 669efb30889SBarry Smith } else { 6709782cf98SBarry Smith dx = xx[col]; 671efb30889SBarry Smith } 6725b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 6739782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX) 6749782cf98SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 6759782cf98SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 6769782cf98SBarry Smith #else 6779782cf98SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 6789782cf98SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 6799782cf98SBarry Smith #endif 6809782cf98SBarry Smith dx *= epsilon; 68130b34957SBarry Smith vscale_array[col] = 1.0/dx; 6829782cf98SBarry Smith } 683b8f8c88eSHong Zhang if (ctype == IS_COLORING_LOCAL) vscale_array = vscale_array + start; 684efb30889SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 685b8f8c88eSHong Zhang if (ctype == IS_COLORING_LOCAL){ 68630b34957SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 68730b34957SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 688b8f8c88eSHong Zhang } 6899782cf98SBarry Smith 690b8f8c88eSHong Zhang if (coloring->vscaleforrow) { 691b8f8c88eSHong Zhang vscaleforrow = coloring->vscaleforrow; 692b8f8c88eSHong Zhang } else { 693b8f8c88eSHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 694b8f8c88eSHong Zhang } 695b0a32e0cSBarry Smith 6969782cf98SBarry Smith /* 69743a90d84SBarry Smith Loop over each color 69843a90d84SBarry Smith */ 699b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 70043a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 70149b058dcSBarry Smith coloring->currentcolor = k; 702b8f8c88eSHong Zhang ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 703b8f8c88eSHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 704b8f8c88eSHong Zhang if (ctype == IS_COLORING_LOCAL) w3_array = w3_array - start; 70543a90d84SBarry Smith /* 706b8f8c88eSHong Zhang Loop over each column associated with color 707b8f8c88eSHong Zhang adding the perturbation to the vector w3. 70843a90d84SBarry Smith */ 70943a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 710b8f8c88eSHong Zhang col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 711efb30889SBarry Smith if (coloring->htype[0] == 'w') { 712efb30889SBarry Smith dx = 1.0 + unorm; 713efb30889SBarry Smith } else { 71442460c72SBarry Smith dx = xx[col]; 715efb30889SBarry Smith } 7165b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 717aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 718ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 719ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 72043a90d84SBarry Smith #else 721329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 722329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 72343a90d84SBarry Smith #endif 72443a90d84SBarry Smith dx *= epsilon; 7251302d50aSBarry Smith if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 72642460c72SBarry Smith w3_array[col] += dx; 72743a90d84SBarry Smith } 728b8f8c88eSHong Zhang if (ctype == IS_COLORING_LOCAL) w3_array = w3_array + start; 729b8f8c88eSHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 7303b28642cSBarry Smith 73143a90d84SBarry Smith /* 732b8f8c88eSHong Zhang Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 733b8f8c88eSHong Zhang w2 = F(x1 + dx) - F(x1) 73443a90d84SBarry Smith */ 73566f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 73643a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 73766f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 738efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 7399782cf98SBarry Smith 74043a90d84SBarry Smith /* 741e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 74243a90d84SBarry Smith */ 7433b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 74443a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 745b8f8c88eSHong Zhang row = coloring->rows[k][l]; /* local row index */ 746b8f8c88eSHong Zhang col = coloring->columnsforrow[k][l]; /* global column index */ 7475904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 74843a90d84SBarry Smith srow = row + start; 74943a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 75043a90d84SBarry Smith } 7513b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 752aeb7e98aSMatthew Knepley } /* endof for each color */ 753b8f8c88eSHong Zhang if (ctype == IS_COLORING_LOCAL) xx = xx + start; 75430b34957SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 755b8f8c88eSHong Zhang ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 756b8f8c88eSHong Zhang 757b8f8c88eSHong Zhang coloring->currentcolor = -1; 75843a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75943a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 760d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 761a2e34c3dSBarry Smith 762b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 763a2e34c3dSBarry Smith if (flg) { 764a2e34c3dSBarry Smith ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 765a2e34c3dSBarry Smith } 766b9722508SLois Curfman McInnes ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 7673a40ed3dSBarry Smith PetscFunctionReturn(0); 76843a90d84SBarry Smith } 769840b8ebdSBarry Smith 7704a2ae208SSatish Balay #undef __FUNCT__ 7714a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS" 772840b8ebdSBarry Smith /*@ 773840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 774840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 775840b8ebdSBarry Smith 776fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 777fee21e36SBarry Smith 778ef5ee4d1SLois Curfman McInnes Input Parameters: 7793b28642cSBarry Smith + mat - location to store Jacobian 780ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 781ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 782ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 783ef5ee4d1SLois Curfman McInnes 784840b8ebdSBarry Smith Options Database Keys: 785ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 786840b8ebdSBarry Smith 78715091d37SBarry Smith Level: intermediate 78815091d37SBarry Smith 789840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 790840b8ebdSBarry Smith 791840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 792840b8ebdSBarry Smith @*/ 793be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 794840b8ebdSBarry Smith { 7956849ba73SBarry Smith PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f; 7966849ba73SBarry Smith PetscErrorCode ierr; 79738baddfdSBarry Smith PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow; 798efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 799ea709b57SSatish Balay PetscScalar *vscale_array; 800329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 801ced766e8SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 802840b8ebdSBarry Smith void *fctx = coloring->fctx; 803f1af5d2fSBarry Smith PetscTruth flg; 804840b8ebdSBarry Smith 8053a40ed3dSBarry Smith PetscFunctionBegin; 8064482741eSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE,1); 8074482741eSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 8084482741eSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE,4); 809840b8ebdSBarry Smith 810d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 811ced766e8SHong Zhang if (!coloring->w3) { 812840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 81352e6d16bSBarry Smith ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 814840b8ebdSBarry Smith } 815ced766e8SHong Zhang w3 = coloring->w3; 816840b8ebdSBarry Smith 8175904e1b1SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 818b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 819f1af5d2fSBarry Smith if (flg) { 820ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 821840b8ebdSBarry Smith } else { 8220b9b6f31SBarry Smith PetscTruth assembled; 8230b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 8240b9b6f31SBarry Smith if (assembled) { 825840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 826840b8ebdSBarry Smith } 8270b9b6f31SBarry Smith } 828840b8ebdSBarry Smith 829840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 830840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 83166f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 832840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 83366f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 834840b8ebdSBarry Smith 835840b8ebdSBarry Smith /* 8365904e1b1SBarry Smith Compute all the scale factors and share with other processors 837840b8ebdSBarry Smith */ 8385904e1b1SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 8395904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 840840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 841840b8ebdSBarry Smith /* 842840b8ebdSBarry Smith Loop over each column associated with color adding the 843840b8ebdSBarry Smith perturbation to the vector w3. 844840b8ebdSBarry Smith */ 845840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 846840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 8475904e1b1SBarry Smith dx = xx[col]; 8485b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 849aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 850840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 851840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 852840b8ebdSBarry Smith #else 853329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 854329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 855840b8ebdSBarry Smith #endif 856840b8ebdSBarry Smith dx *= epsilon; 8575904e1b1SBarry Smith vscale_array[col] = 1.0/dx; 858840b8ebdSBarry Smith } 8595904e1b1SBarry Smith } 8605904e1b1SBarry Smith vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8615904e1b1SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8625904e1b1SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8635904e1b1SBarry Smith 8645904e1b1SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 8655904e1b1SBarry Smith else vscaleforrow = coloring->columnsforrow; 8665904e1b1SBarry Smith 8675904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8685904e1b1SBarry Smith /* 8695904e1b1SBarry Smith Loop over each color 8705904e1b1SBarry Smith */ 8715904e1b1SBarry Smith for (k=0; k<coloring->ncolors; k++) { 8725904e1b1SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 8735904e1b1SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 8745904e1b1SBarry Smith /* 8755904e1b1SBarry Smith Loop over each column associated with color adding the 8765904e1b1SBarry Smith perturbation to the vector w3. 8775904e1b1SBarry Smith */ 8785904e1b1SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 8795904e1b1SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 8805904e1b1SBarry Smith dx = xx[col]; 8815b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 8825904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 8835904e1b1SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 8845904e1b1SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 8855904e1b1SBarry Smith #else 8865904e1b1SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 8875904e1b1SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 8885904e1b1SBarry Smith #endif 8895904e1b1SBarry Smith dx *= epsilon; 8905904e1b1SBarry Smith w3_array[col] += dx; 8915904e1b1SBarry Smith } 8925904e1b1SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 8935904e1b1SBarry Smith 894840b8ebdSBarry Smith /* 895840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 896840b8ebdSBarry Smith */ 89766f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 898840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 89966f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 900efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 9015904e1b1SBarry Smith 902840b8ebdSBarry Smith /* 903840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 904840b8ebdSBarry Smith */ 9053b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 906840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 907840b8ebdSBarry Smith row = coloring->rows[k][l]; 908840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 9095904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 910840b8ebdSBarry Smith srow = row + start; 911840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 912840b8ebdSBarry Smith } 9133b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 914840b8ebdSBarry Smith } 9155904e1b1SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 9165904e1b1SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 917840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 918840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 919d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 9203a40ed3dSBarry Smith PetscFunctionReturn(0); 921840b8ebdSBarry Smith } 9223b28642cSBarry Smith 9233b28642cSBarry Smith 9244a2ae208SSatish Balay #undef __FUNCT__ 9254a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()" 92640964ad5SBarry Smith /*@C 92740964ad5SBarry Smith MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 92840964ad5SBarry Smith is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 92940964ad5SBarry Smith no additional Jacobian's will be computed (the same one will be used) until this is 93040964ad5SBarry Smith called again. 93140964ad5SBarry Smith 93240964ad5SBarry Smith Collective on MatFDColoring 93340964ad5SBarry Smith 93440964ad5SBarry Smith Input Parameters: 93540964ad5SBarry Smith . c - the coloring context 93640964ad5SBarry Smith 93740964ad5SBarry Smith Level: intermediate 93840964ad5SBarry Smith 93940964ad5SBarry Smith Notes: The MatFDColoringSetFrequency() is ignored once this is called 94040964ad5SBarry Smith 94140964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 94240964ad5SBarry Smith 94340964ad5SBarry Smith .keywords: Mat, finite differences, coloring 94440964ad5SBarry Smith @*/ 945be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring c) 94640964ad5SBarry Smith { 94740964ad5SBarry Smith PetscFunctionBegin; 9484482741eSBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1); 94940964ad5SBarry Smith c->usersetsrecompute = PETSC_TRUE; 95040964ad5SBarry Smith c->recompute = PETSC_TRUE; 95140964ad5SBarry Smith PetscFunctionReturn(0); 95240964ad5SBarry Smith } 95340964ad5SBarry Smith 9543b28642cSBarry Smith 955