1be1d678aSKris Buschelman 2bbf0e169SBarry Smith /* 3639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 4639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 5bbf0e169SBarry Smith */ 6bbf0e169SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 8af0996ceSBarry Smith #include <petsc/private/isimpl.h> 9bbf0e169SBarry Smith 107087cfbeSBarry Smith PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F) 113a7fca6bSBarry Smith { 124e269d77SPeter Brune PetscErrorCode ierr; 136e111a19SKarl Rupp 143a7fca6bSBarry Smith PetscFunctionBegin; 154e269d77SPeter Brune if (F) { 164e269d77SPeter Brune ierr = VecCopy(F,fd->w1);CHKERRQ(ierr); 174e269d77SPeter Brune fd->fset = PETSC_TRUE; 184e269d77SPeter Brune } else { 194e269d77SPeter Brune fd->fset = PETSC_FALSE; 204e269d77SPeter Brune } 213a7fca6bSBarry Smith PetscFunctionReturn(0); 223a7fca6bSBarry Smith } 233a7fca6bSBarry Smith 249804daf3SBarry Smith #include <petscdraw.h> 256849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 269194eea9SBarry Smith { 279194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 28dfbe8321SBarry Smith PetscErrorCode ierr; 29b312708bSHong Zhang PetscInt i,j,nz,row; 309194eea9SBarry Smith PetscReal x,y; 31b312708bSHong Zhang MatEntry *Jentry=fd->matentry; 329194eea9SBarry Smith 339194eea9SBarry Smith PetscFunctionBegin; 349194eea9SBarry Smith /* loop over colors */ 35b312708bSHong Zhang nz = 0; 369194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 379194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 38b312708bSHong Zhang row = Jentry[nz].row; 39b312708bSHong Zhang y = fd->M - row - fd->rstart; 40b312708bSHong Zhang x = (PetscReal)Jentry[nz++].col; 41b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 429194eea9SBarry Smith } 439194eea9SBarry Smith } 449194eea9SBarry Smith PetscFunctionReturn(0); 459194eea9SBarry Smith } 469194eea9SBarry Smith 476849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 48005c665bSBarry Smith { 49dfbe8321SBarry Smith PetscErrorCode ierr; 50ace3abfcSBarry Smith PetscBool isnull; 51b0a32e0cSBarry Smith PetscDraw draw; 5254d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 53005c665bSBarry Smith 543a40ed3dSBarry Smith PetscFunctionBegin; 55b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 56832b7cebSLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 57832b7cebSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 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); 62832b7cebSLisandro Dalcin ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 63b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 64f77a5242SKarl Rupp ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr); 65832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 663a40ed3dSBarry Smith PetscFunctionReturn(0); 67005c665bSBarry Smith } 68005c665bSBarry Smith 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 @*/ 987087cfbeSBarry Smith PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 99bbf0e169SBarry Smith { 1006849ba73SBarry Smith PetscErrorCode ierr; 10138baddfdSBarry Smith PetscInt i,j; 102ace3abfcSBarry Smith PetscBool isdraw,iascii; 103f3ef73ceSBarry Smith PetscViewerFormat format; 104bbf0e169SBarry Smith 1053a40ed3dSBarry Smith PetscFunctionBegin; 1060700a824SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 1073050cee2SBarry Smith if (!viewer) { 108ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr); 1093050cee2SBarry Smith } 1100700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 111c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 112bbf0e169SBarry Smith 113251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 114251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1150f5bd95cSBarry Smith if (isdraw) { 116b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 11732077d6dSBarry Smith } else if (iascii) { 118dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr); 11957622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr); 12057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",(double)c->umin);CHKERRQ(ierr); 12177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 122ae09f205SBarry Smith 123b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 124fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 125b312708bSHong Zhang PetscInt row,col,nz; 126b312708bSHong Zhang nz = 0; 127b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 12877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 12977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 130b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 13177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 132639f9d9dSBarry Smith } 13377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 1345bdb020cSBarry Smith if (c->matentry) { 135b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 136b312708bSHong Zhang row = c->matentry[nz].row; 137b312708bSHong Zhang col = c->matentry[nz++].col; 138b312708bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",row,col);CHKERRQ(ierr); 139b4fc646aSLois Curfman McInnes } 140bbf0e169SBarry Smith } 141bbf0e169SBarry Smith } 1425bdb020cSBarry Smith } 143b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 144bbf0e169SBarry Smith } 1453a40ed3dSBarry Smith PetscFunctionReturn(0); 146639f9d9dSBarry Smith } 147639f9d9dSBarry Smith 148639f9d9dSBarry Smith /*@ 149b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 150b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 151639f9d9dSBarry Smith 1523f9fe445SBarry Smith Logically Collective on MatFDColoring 153ef5ee4d1SLois Curfman McInnes 154ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 155ef5ee4d1SLois Curfman McInnes .vb 15665f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 157*d461c19dSHong Zhang htype = 'ds': 158f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 159f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 160ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 161*d461c19dSHong Zhang 162*d461c19dSHong Zhang htype = 'wp': 163*d461c19dSHong Zhang h = error_rel * sqrt(1 + ||u||) 164ef5ee4d1SLois Curfman McInnes .ve 165639f9d9dSBarry Smith 166639f9d9dSBarry Smith Input Parameters: 167ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 168639f9d9dSBarry Smith . error_rel - relative error 169f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 170fee21e36SBarry Smith 17115091d37SBarry Smith Level: advanced 17215091d37SBarry Smith 173b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 174b4fc646aSLois Curfman McInnes 17595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 17695b89fc3SBarry Smith 177639f9d9dSBarry Smith @*/ 1787087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 179639f9d9dSBarry Smith { 1803a40ed3dSBarry Smith PetscFunctionBegin; 1810700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 182c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,error,2); 183c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,umin,3); 184639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 185639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1863a40ed3dSBarry Smith PetscFunctionReturn(0); 187639f9d9dSBarry Smith } 188639f9d9dSBarry Smith 189f86b9fbaSHong Zhang /*@ 190c8a9c622SHong Zhang MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix. 191005c665bSBarry Smith 192f86b9fbaSHong Zhang Logically Collective on MatFDColoring 193f86b9fbaSHong Zhang 194f86b9fbaSHong Zhang Input Parameters: 195f86b9fbaSHong Zhang + coloring - the coloring context 196f86b9fbaSHong Zhang . brows - number of rows in the block 197f86b9fbaSHong Zhang - bcols - number of columns in the block 198f86b9fbaSHong Zhang 199f86b9fbaSHong Zhang Level: intermediate 200f86b9fbaSHong Zhang 201f86b9fbaSHong Zhang .keywords: Mat, coloring 202f86b9fbaSHong Zhang 203f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 204f86b9fbaSHong Zhang 205f86b9fbaSHong Zhang @*/ 206f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols) 207f86b9fbaSHong Zhang { 208f86b9fbaSHong Zhang PetscFunctionBegin; 209f86b9fbaSHong Zhang PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 210f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,brows,2); 211f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,bcols,3); 212f86b9fbaSHong Zhang if (brows != PETSC_DEFAULT) matfd->brows = brows; 213f86b9fbaSHong Zhang if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 214f86b9fbaSHong Zhang PetscFunctionReturn(0); 215f86b9fbaSHong Zhang } 216f86b9fbaSHong Zhang 217f86b9fbaSHong Zhang /*@ 218f86b9fbaSHong Zhang MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 219f86b9fbaSHong Zhang 220f86b9fbaSHong Zhang Collective on Mat 221f86b9fbaSHong Zhang 222f86b9fbaSHong Zhang Input Parameters: 223f86b9fbaSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 224f86b9fbaSHong Zhang . iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 225f86b9fbaSHong Zhang - color - the matrix coloring context 226f86b9fbaSHong Zhang 227f86b9fbaSHong Zhang Level: beginner 228f86b9fbaSHong Zhang 229f86b9fbaSHong Zhang .keywords: MatFDColoring, setup 230f86b9fbaSHong Zhang 231f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy() 232f86b9fbaSHong Zhang @*/ 233f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color) 234f86b9fbaSHong Zhang { 235f86b9fbaSHong Zhang PetscErrorCode ierr; 236f86b9fbaSHong Zhang 237f86b9fbaSHong Zhang PetscFunctionBegin; 238c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 239c8a9c622SHong Zhang PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3); 240c8a9c622SHong Zhang if (color->setupcalled) PetscFunctionReturn(0); 241c8a9c622SHong Zhang 2420df34763SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 243f86b9fbaSHong Zhang if (mat->ops->fdcoloringsetup) { 244f86b9fbaSHong Zhang ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr); 245f86b9fbaSHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 246c8a9c622SHong Zhang 247c8a9c622SHong Zhang color->setupcalled = PETSC_TRUE; 2480df34763SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 249f86b9fbaSHong Zhang PetscFunctionReturn(0); 250f86b9fbaSHong Zhang } 251ff0cfa39SBarry Smith 2524a9d489dSBarry Smith /*@C 2534a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2544a9d489dSBarry Smith 2553f9fe445SBarry Smith Not Collective 2564a9d489dSBarry Smith 2574a9d489dSBarry Smith Input Parameters: 2584a9d489dSBarry Smith . coloring - the coloring context 2594a9d489dSBarry Smith 2604a9d489dSBarry Smith Output Parameters: 2614a9d489dSBarry Smith + f - the function 2624a9d489dSBarry Smith - fctx - the optional user-defined function context 2634a9d489dSBarry Smith 2644a9d489dSBarry Smith Level: intermediate 2654a9d489dSBarry Smith 2664a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 26795b89fc3SBarry Smith 26895b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 26995b89fc3SBarry Smith 2704a9d489dSBarry Smith @*/ 2717087cfbeSBarry Smith PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2724a9d489dSBarry Smith { 2734a9d489dSBarry Smith PetscFunctionBegin; 2740700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 2754a9d489dSBarry Smith if (f) *f = matfd->f; 2764a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2774a9d489dSBarry Smith PetscFunctionReturn(0); 2784a9d489dSBarry Smith } 2794a9d489dSBarry Smith 280d64ed03dSBarry Smith /*@C 281005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 282005c665bSBarry Smith 2833f9fe445SBarry Smith Logically Collective on MatFDColoring 284fee21e36SBarry Smith 285ef5ee4d1SLois Curfman McInnes Input Parameters: 286ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 287ef5ee4d1SLois Curfman McInnes . f - the function 288ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 289ef5ee4d1SLois Curfman McInnes 2907850c7c0SBarry Smith Calling sequence of (*f) function: 2917850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 292ab637aeaSJed Brown If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 29315091d37SBarry Smith 2947850c7c0SBarry Smith Level: advanced 2957850c7c0SBarry Smith 296ab637aeaSJed Brown Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 2978d359177SBarry Smith SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by 298ab637aeaSJed Brown calling MatFDColoringApply() 2997850c7c0SBarry Smith 3007850c7c0SBarry Smith Fortran Notes: 3017850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 302ab637aeaSJed Brown be used without SNES or within the SNES solvers. 303f881d145SBarry Smith 304b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 30595b89fc3SBarry Smith 30695b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 30795b89fc3SBarry Smith 308005c665bSBarry Smith @*/ 3097087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 310005c665bSBarry Smith { 3113a40ed3dSBarry Smith PetscFunctionBegin; 3120700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 313005c665bSBarry Smith matfd->f = f; 314005c665bSBarry Smith matfd->fctx = fctx; 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 316005c665bSBarry Smith } 317005c665bSBarry Smith 318639f9d9dSBarry Smith /*@ 319b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 320639f9d9dSBarry Smith the options database. 321639f9d9dSBarry Smith 322fee21e36SBarry Smith Collective on MatFDColoring 323fee21e36SBarry Smith 32465f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 325ef5ee4d1SLois Curfman McInnes .vb 32665f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 327f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 328f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 329ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 330ef5ee4d1SLois Curfman McInnes .ve 331ef5ee4d1SLois Curfman McInnes 332ef5ee4d1SLois Curfman McInnes Input Parameter: 333ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 334ef5ee4d1SLois Curfman McInnes 335b4fc646aSLois Curfman McInnes Options Database Keys: 3365620d6dcSBarry Smith + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 337f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 3383ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 339ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 340e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info 341e350db43SBarry Smith - -mat_fd_coloring_view draw - Activates drawing 342639f9d9dSBarry Smith 34315091d37SBarry Smith Level: intermediate 34415091d37SBarry Smith 345b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 346d1c39f55SBarry Smith 347d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 348d1c39f55SBarry Smith 349639f9d9dSBarry Smith @*/ 3507087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 351639f9d9dSBarry Smith { 352dfbe8321SBarry Smith PetscErrorCode ierr; 353ace3abfcSBarry Smith PetscBool flg; 354efb30889SBarry Smith char value[3]; 3553a40ed3dSBarry Smith 3563a40ed3dSBarry Smith PetscFunctionBegin; 3570700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 358639f9d9dSBarry Smith 3593194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr); 36087828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 36187828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 3621bc75a8dSBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr); 363efb30889SBarry Smith if (flg) { 364efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 365efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 366e32f2f54SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 367efb30889SBarry Smith } 368f86b9fbaSHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr); 36993dfae19SHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr); 37093dfae19SHong Zhang if (flg && matfd->bcols > matfd->ncolors) { 37193dfae19SHong Zhang /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */ 37293dfae19SHong Zhang matfd->bcols = matfd->ncolors; 37393dfae19SHong Zhang } 374f86b9fbaSHong Zhang 3755d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 3760633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr); 377b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 379005c665bSBarry Smith } 380005c665bSBarry Smith 38161ab5df0SBarry Smith /*@C 38261ab5df0SBarry Smith MatFDColoringSetType - Sets the approach for computing the finite difference parameter 38361ab5df0SBarry Smith 38461ab5df0SBarry Smith Collective on MatFDColoring 38561ab5df0SBarry Smith 38661ab5df0SBarry Smith Input Parameters: 38761ab5df0SBarry Smith + coloring - the coloring context 38861ab5df0SBarry Smith - type - either MATMFFD_WP or MATMFFD_DS 38961ab5df0SBarry Smith 39061ab5df0SBarry Smith Options Database Keys: 39161ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds" 39261ab5df0SBarry Smith 39361ab5df0SBarry Smith Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries 39461ab5df0SBarry Smith but the process of computing the entries is the same as as with the MatMFFD operation so we should reuse the names instead of 39561ab5df0SBarry Smith introducing another one. 39661ab5df0SBarry Smith 39761ab5df0SBarry Smith Level: intermediate 39861ab5df0SBarry Smith 39961ab5df0SBarry Smith .keywords: Mat, finite differences, parameters 40061ab5df0SBarry Smith 40161ab5df0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 40261ab5df0SBarry Smith 40361ab5df0SBarry Smith @*/ 40461ab5df0SBarry Smith PetscErrorCode MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type) 40561ab5df0SBarry Smith { 40661ab5df0SBarry Smith PetscFunctionBegin; 40761ab5df0SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 40861ab5df0SBarry Smith /* 40961ab5df0SBarry Smith It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype 41061ab5df0SBarry Smith and this function is being provided as patch for a release so it shouldn't change the implementaton 41161ab5df0SBarry Smith */ 41261ab5df0SBarry Smith if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp"; 41361ab5df0SBarry Smith else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds"; 41461ab5df0SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type); 41561ab5df0SBarry Smith PetscFunctionReturn(0); 41661ab5df0SBarry Smith } 41761ab5df0SBarry Smith 418146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[]) 419005c665bSBarry Smith { 420dfbe8321SBarry Smith PetscErrorCode ierr; 421e350db43SBarry Smith PetscBool flg; 4223050cee2SBarry Smith PetscViewer viewer; 423cffb1e40SBarry Smith PetscViewerFormat format; 424005c665bSBarry Smith 4253a40ed3dSBarry Smith PetscFunctionBegin; 426146574abSBarry Smith if (prefix) { 427146574abSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 428146574abSBarry Smith } else { 429ce94432eSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 430146574abSBarry Smith } 431005c665bSBarry Smith if (flg) { 432cffb1e40SBarry Smith ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 4333050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 434cffb1e40SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 435cffb1e40SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 436005c665bSBarry Smith } 4373a40ed3dSBarry Smith PetscFunctionReturn(0); 438bbf0e169SBarry Smith } 439bbf0e169SBarry Smith 44005869f15SSatish Balay /*@ 441639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 442639f9d9dSBarry Smith computation of Jacobians. 443bbf0e169SBarry Smith 444ef5ee4d1SLois Curfman McInnes Collective on Mat 445ef5ee4d1SLois Curfman McInnes 446639f9d9dSBarry Smith Input Parameters: 447ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 4487086a01eSPeter Brune - iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring() 449bbf0e169SBarry Smith 450bbf0e169SBarry Smith Output Parameter: 451639f9d9dSBarry Smith . color - the new coloring context 452bbf0e169SBarry Smith 45315091d37SBarry Smith Level: intermediate 45415091d37SBarry Smith 4558d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(), 456d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 4577086a01eSPeter Brune MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring() 458bbf0e169SBarry Smith @*/ 4597087cfbeSBarry Smith PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 460bbf0e169SBarry Smith { 461639f9d9dSBarry Smith MatFDColoring c; 462639f9d9dSBarry Smith MPI_Comm comm; 463dfbe8321SBarry Smith PetscErrorCode ierr; 46438baddfdSBarry Smith PetscInt M,N; 465639f9d9dSBarry Smith 4663a40ed3dSBarry Smith PetscFunctionBegin; 467c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 468f141eedfSHong Zhang if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();"); 469d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 470639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 471ce94432eSBarry Smith if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices"); 472f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 47373107ff1SLisandro Dalcin ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 474639f9d9dSBarry Smith 475b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 476b8f8c88eSHong Zhang 47793dfae19SHong Zhang if (mat->ops->fdcoloringcreate) { 47893dfae19SHong Zhang ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 47993dfae19SHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 48093dfae19SHong Zhang 4812a7a6963SBarry Smith ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr); 4823bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr); 483b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 4843bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr); 485b8f8c88eSHong Zhang 48677d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 48777d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 48849b058dcSBarry Smith c->currentcolor = -1; 489efb30889SBarry Smith c->htype = "wp"; 4904e269d77SPeter Brune c->fset = PETSC_FALSE; 491c8a9c622SHong Zhang c->setupcalled = PETSC_FALSE; 492639f9d9dSBarry Smith 493639f9d9dSBarry Smith *color = c; 4944e269d77SPeter Brune ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr); 495d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4963a40ed3dSBarry Smith PetscFunctionReturn(0); 497bbf0e169SBarry Smith } 498bbf0e169SBarry Smith 49905869f15SSatish Balay /*@ 500639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 501639f9d9dSBarry Smith via MatFDColoringCreate(). 502bbf0e169SBarry Smith 503ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 504ef5ee4d1SLois Curfman McInnes 505b4fc646aSLois Curfman McInnes Input Parameter: 506639f9d9dSBarry Smith . c - coloring context 507bbf0e169SBarry Smith 50815091d37SBarry Smith Level: intermediate 50915091d37SBarry Smith 510639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 511bbf0e169SBarry Smith @*/ 5126bf464f9SBarry Smith PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 513bbf0e169SBarry Smith { 5146849ba73SBarry Smith PetscErrorCode ierr; 51538baddfdSBarry Smith PetscInt i; 5160df34763SHong Zhang MatFDColoring color = *c; 517bbf0e169SBarry Smith 5183a40ed3dSBarry Smith PetscFunctionBegin; 5196bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 5200df34763SHong Zhang if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);} 521d4bb536fSBarry Smith 5220df34763SHong Zhang for (i=0; i<color->ncolors; i++) { 5230df34763SHong Zhang ierr = PetscFree(color->columns[i]);CHKERRQ(ierr); 524bbf0e169SBarry Smith } 5250df34763SHong Zhang ierr = PetscFree(color->ncolumns);CHKERRQ(ierr); 5260df34763SHong Zhang ierr = PetscFree(color->columns);CHKERRQ(ierr); 5270df34763SHong Zhang ierr = PetscFree(color->nrows);CHKERRQ(ierr); 5280df34763SHong Zhang if (color->htype[0] == 'w') { 5290df34763SHong Zhang ierr = PetscFree(color->matentry2);CHKERRQ(ierr); 5300df34763SHong Zhang } else { 5310df34763SHong Zhang ierr = PetscFree(color->matentry);CHKERRQ(ierr); 5320df34763SHong Zhang } 5330df34763SHong Zhang ierr = PetscFree(color->dy);CHKERRQ(ierr); 5340df34763SHong Zhang if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);} 5350df34763SHong Zhang ierr = VecDestroy(&color->w1);CHKERRQ(ierr); 5360df34763SHong Zhang ierr = VecDestroy(&color->w2);CHKERRQ(ierr); 5370df34763SHong Zhang ierr = VecDestroy(&color->w3);CHKERRQ(ierr); 538d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 5393a40ed3dSBarry Smith PetscFunctionReturn(0); 540bbf0e169SBarry Smith } 54143a90d84SBarry Smith 54249b058dcSBarry Smith /*@C 54349b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 54449b058dcSBarry Smith that are currently being perturbed. 54549b058dcSBarry Smith 54649b058dcSBarry Smith Not Collective 54749b058dcSBarry Smith 54849b058dcSBarry Smith Input Parameters: 54949b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 55049b058dcSBarry Smith 55149b058dcSBarry Smith Output Parameters: 55249b058dcSBarry Smith + n - the number of local columns being perturbed 55349b058dcSBarry Smith - cols - the column indices, in global numbering 55449b058dcSBarry Smith 55521fcc2ddSBarry Smith Level: advanced 55649b058dcSBarry Smith 557edaa7c33SBarry Smith Fortran Note: 558edaa7c33SBarry Smith This routine has a different interface for Fortran 55921fcc2ddSBarry Smith $ #include <petsc/finclude/petscmat.h> 56021fcc2ddSBarry Smith $ use petscmat 561edaa7c33SBarry Smith $ PetscInt, pointer :: array(:) 562edaa7c33SBarry Smith $ PetscErrorCode ierr 563edaa7c33SBarry Smith $ MatFDColoring i 564edaa7c33SBarry Smith $ call MatFDColoringGetPerturbedColumnsF90(i,array,ierr) 565edaa7c33SBarry Smith $ use the entries of array ... 566edaa7c33SBarry Smith $ call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr) 567edaa7c33SBarry Smith 56849b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 56949b058dcSBarry Smith 57049b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 57149b058dcSBarry Smith @*/ 572edaa7c33SBarry Smith PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[]) 57349b058dcSBarry Smith { 57449b058dcSBarry Smith PetscFunctionBegin; 57549b058dcSBarry Smith if (coloring->currentcolor >= 0) { 57649b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 57749b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 57849b058dcSBarry Smith } else { 57949b058dcSBarry Smith *n = 0; 58049b058dcSBarry Smith } 58149b058dcSBarry Smith PetscFunctionReturn(0); 58249b058dcSBarry Smith } 58349b058dcSBarry Smith 58443a90d84SBarry Smith /*@ 585e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 586e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 58743a90d84SBarry Smith 588fee21e36SBarry Smith Collective on MatFDColoring 589fee21e36SBarry Smith 590ef5ee4d1SLois Curfman McInnes Input Parameters: 591ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 592ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 593ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 5947850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 595ef5ee4d1SLois Curfman McInnes 5968bba8e72SBarry Smith Options Database Keys: 597ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 598b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 599e350db43SBarry Smith . -mat_fd_coloring_view draw - Activates drawing of coloring 600e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 6018bba8e72SBarry Smith 60215091d37SBarry Smith Level: intermediate 60398d79826SSatish Balay 6047850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 60543a90d84SBarry Smith 60643a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 60743a90d84SBarry Smith @*/ 608d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 60943a90d84SBarry Smith { 6103acb8795SBarry Smith PetscErrorCode ierr; 611684f2004SHong Zhang PetscBool flg = PETSC_FALSE; 6123acb8795SBarry Smith 6133acb8795SBarry Smith PetscFunctionBegin; 6140700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 6150700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 6160700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 617ce94432eSBarry Smith if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 618e32f2f54SBarry Smith if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 619c8a9c622SHong Zhang if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()"); 620684f2004SHong Zhang 621684f2004SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 622c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 623684f2004SHong Zhang if (flg) { 624684f2004SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 625684f2004SHong Zhang } else { 626684f2004SHong Zhang PetscBool assembled; 627684f2004SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 628684f2004SHong Zhang if (assembled) { 629684f2004SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 630684f2004SHong Zhang } 631684f2004SHong Zhang } 632684f2004SHong Zhang 6335922145eSHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 634d1e9a80fSBarry Smith ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr); 6355922145eSHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 6365bdb020cSBarry Smith if (!coloring->viewed) { 6375bdb020cSBarry Smith ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr); 6385bdb020cSBarry Smith coloring->viewed = PETSC_TRUE; 6395bdb020cSBarry Smith } 6403acb8795SBarry Smith PetscFunctionReturn(0); 6413acb8795SBarry Smith } 642