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 104a2ae208SSatish Balay #undef __FUNCT__ 113a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF" 127087cfbeSBarry Smith PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F) 133a7fca6bSBarry Smith { 144e269d77SPeter Brune PetscErrorCode ierr; 156e111a19SKarl Rupp 163a7fca6bSBarry Smith PetscFunctionBegin; 174e269d77SPeter Brune if (F) { 184e269d77SPeter Brune ierr = VecCopy(F,fd->w1);CHKERRQ(ierr); 194e269d77SPeter Brune fd->fset = PETSC_TRUE; 204e269d77SPeter Brune } else { 214e269d77SPeter Brune fd->fset = PETSC_FALSE; 224e269d77SPeter Brune } 233a7fca6bSBarry Smith PetscFunctionReturn(0); 243a7fca6bSBarry Smith } 253a7fca6bSBarry Smith 269804daf3SBarry Smith #include <petscdraw.h> 273a7fca6bSBarry Smith #undef __FUNCT__ 284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 296849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 309194eea9SBarry Smith { 319194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 32dfbe8321SBarry Smith PetscErrorCode ierr; 33b312708bSHong Zhang PetscInt i,j,nz,row; 349194eea9SBarry Smith PetscReal x,y; 35b312708bSHong Zhang MatEntry *Jentry=fd->matentry; 369194eea9SBarry Smith 379194eea9SBarry Smith PetscFunctionBegin; 389194eea9SBarry Smith /* loop over colors */ 39b312708bSHong Zhang nz = 0; 409194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 419194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 42b312708bSHong Zhang row = Jentry[nz].row; 43b312708bSHong Zhang y = fd->M - row - fd->rstart; 44b312708bSHong Zhang x = (PetscReal)Jentry[nz++].col; 45b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 469194eea9SBarry Smith } 479194eea9SBarry Smith } 489194eea9SBarry Smith PetscFunctionReturn(0); 499194eea9SBarry Smith } 509194eea9SBarry Smith 514a2ae208SSatish Balay #undef __FUNCT__ 524a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 536849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 54005c665bSBarry Smith { 55dfbe8321SBarry Smith PetscErrorCode ierr; 56ace3abfcSBarry Smith PetscBool isnull; 57b0a32e0cSBarry Smith PetscDraw draw; 5854d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 59005c665bSBarry Smith 603a40ed3dSBarry Smith PetscFunctionBegin; 61b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 62832b7cebSLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 63832b7cebSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 64005c665bSBarry Smith 65005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 66005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 67b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 68832b7cebSLisandro Dalcin ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 69b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 70f77a5242SKarl Rupp ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr); 71832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 723a40ed3dSBarry Smith PetscFunctionReturn(0); 73005c665bSBarry Smith } 74005c665bSBarry Smith 754a2ae208SSatish Balay #undef __FUNCT__ 764a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 77bbf0e169SBarry Smith /*@C 78639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 79bbf0e169SBarry Smith 804c49b128SBarry Smith Collective on MatFDColoring 81fee21e36SBarry Smith 82ef5ee4d1SLois Curfman McInnes Input Parameters: 83ef5ee4d1SLois Curfman McInnes + c - the coloring context 84ef5ee4d1SLois Curfman McInnes - viewer - visualization context 85ef5ee4d1SLois Curfman McInnes 8615091d37SBarry Smith Level: intermediate 8715091d37SBarry Smith 88b4fc646aSLois Curfman McInnes Notes: 89b4fc646aSLois Curfman McInnes The available visualization contexts include 90b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 91b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 92ef5ee4d1SLois Curfman McInnes output where only the first processor opens 93ef5ee4d1SLois Curfman McInnes the file. All other processors send their 94ef5ee4d1SLois Curfman McInnes data to the first processor to print. 95b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 96bbf0e169SBarry Smith 979194eea9SBarry Smith Notes: 989194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 99b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 1009194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 1019194eea9SBarry Smith 102639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 103005c665bSBarry Smith 104b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 105bbf0e169SBarry Smith @*/ 1067087cfbeSBarry Smith PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 107bbf0e169SBarry Smith { 1086849ba73SBarry Smith PetscErrorCode ierr; 10938baddfdSBarry Smith PetscInt i,j; 110ace3abfcSBarry Smith PetscBool isdraw,iascii; 111f3ef73ceSBarry Smith PetscViewerFormat format; 112bbf0e169SBarry Smith 1133a40ed3dSBarry Smith PetscFunctionBegin; 1140700a824SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 1153050cee2SBarry Smith if (!viewer) { 116ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr); 1173050cee2SBarry Smith } 1180700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 119c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 120bbf0e169SBarry Smith 121251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 122251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1230f5bd95cSBarry Smith if (isdraw) { 124b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 12532077d6dSBarry Smith } else if (iascii) { 126dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr); 12757622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",(double)c->error_rel);CHKERRQ(ierr); 12857622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",(double)c->umin);CHKERRQ(ierr); 12977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 130ae09f205SBarry Smith 131b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 132fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 133b312708bSHong Zhang PetscInt row,col,nz; 134b312708bSHong Zhang nz = 0; 135b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 13677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 13777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 138b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 13977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 140639f9d9dSBarry Smith } 14177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 142*5bdb020cSBarry Smith if (c->matentry) { 143b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 144b312708bSHong Zhang row = c->matentry[nz].row; 145b312708bSHong Zhang col = c->matentry[nz++].col; 146b312708bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",row,col);CHKERRQ(ierr); 147b4fc646aSLois Curfman McInnes } 148bbf0e169SBarry Smith } 149bbf0e169SBarry Smith } 150*5bdb020cSBarry Smith } 151b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 152bbf0e169SBarry Smith } 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 154639f9d9dSBarry Smith } 155639f9d9dSBarry Smith 1564a2ae208SSatish Balay #undef __FUNCT__ 1574a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 158639f9d9dSBarry Smith /*@ 159b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 160b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 161639f9d9dSBarry Smith 1623f9fe445SBarry Smith Logically Collective on MatFDColoring 163ef5ee4d1SLois Curfman McInnes 164ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 165ef5ee4d1SLois Curfman McInnes .vb 16665f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 167f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 168f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 169ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 170ef5ee4d1SLois Curfman McInnes .ve 171639f9d9dSBarry Smith 172639f9d9dSBarry Smith Input Parameters: 173ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 174639f9d9dSBarry Smith . error_rel - relative error 175f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 176fee21e36SBarry Smith 17715091d37SBarry Smith Level: advanced 17815091d37SBarry Smith 179b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 180b4fc646aSLois Curfman McInnes 18195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 18295b89fc3SBarry Smith 183639f9d9dSBarry Smith @*/ 1847087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 185639f9d9dSBarry Smith { 1863a40ed3dSBarry Smith PetscFunctionBegin; 1870700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 188c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,error,2); 189c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,umin,3); 190639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 191639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 193639f9d9dSBarry Smith } 194639f9d9dSBarry Smith 195f86b9fbaSHong Zhang #undef __FUNCT__ 196f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetBlockSize" 197f86b9fbaSHong Zhang /*@ 198c8a9c622SHong Zhang MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix. 199005c665bSBarry Smith 200f86b9fbaSHong Zhang Logically Collective on MatFDColoring 201f86b9fbaSHong Zhang 202f86b9fbaSHong Zhang Input Parameters: 203f86b9fbaSHong Zhang + coloring - the coloring context 204f86b9fbaSHong Zhang . brows - number of rows in the block 205f86b9fbaSHong Zhang - bcols - number of columns in the block 206f86b9fbaSHong Zhang 207f86b9fbaSHong Zhang Level: intermediate 208f86b9fbaSHong Zhang 209f86b9fbaSHong Zhang .keywords: Mat, coloring 210f86b9fbaSHong Zhang 211f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 212f86b9fbaSHong Zhang 213f86b9fbaSHong Zhang @*/ 214f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols) 215f86b9fbaSHong Zhang { 216f86b9fbaSHong Zhang PetscFunctionBegin; 217f86b9fbaSHong Zhang PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 218f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,brows,2); 219f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,bcols,3); 220f86b9fbaSHong Zhang if (brows != PETSC_DEFAULT) matfd->brows = brows; 221f86b9fbaSHong Zhang if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 222f86b9fbaSHong Zhang PetscFunctionReturn(0); 223f86b9fbaSHong Zhang } 224f86b9fbaSHong Zhang 225f86b9fbaSHong Zhang #undef __FUNCT__ 226f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp" 227f86b9fbaSHong Zhang /*@ 228f86b9fbaSHong Zhang MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 229f86b9fbaSHong Zhang 230f86b9fbaSHong Zhang Collective on Mat 231f86b9fbaSHong Zhang 232f86b9fbaSHong Zhang Input Parameters: 233f86b9fbaSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 234f86b9fbaSHong Zhang . iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 235f86b9fbaSHong Zhang - color - the matrix coloring context 236f86b9fbaSHong Zhang 237f86b9fbaSHong Zhang Level: beginner 238f86b9fbaSHong Zhang 239f86b9fbaSHong Zhang .keywords: MatFDColoring, setup 240f86b9fbaSHong Zhang 241f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy() 242f86b9fbaSHong Zhang @*/ 243f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color) 244f86b9fbaSHong Zhang { 245f86b9fbaSHong Zhang PetscErrorCode ierr; 246f86b9fbaSHong Zhang 247f86b9fbaSHong Zhang PetscFunctionBegin; 248c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 249c8a9c622SHong Zhang PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3); 250c8a9c622SHong Zhang if (color->setupcalled) PetscFunctionReturn(0); 251c8a9c622SHong Zhang 2520df34763SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 253f86b9fbaSHong Zhang if (mat->ops->fdcoloringsetup) { 254f86b9fbaSHong Zhang ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr); 255f86b9fbaSHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 256c8a9c622SHong Zhang 257c8a9c622SHong Zhang color->setupcalled = PETSC_TRUE; 2580df34763SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0);CHKERRQ(ierr); 259f86b9fbaSHong Zhang PetscFunctionReturn(0); 260f86b9fbaSHong Zhang } 261ff0cfa39SBarry Smith 2624a2ae208SSatish Balay #undef __FUNCT__ 2634a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction" 2644a9d489dSBarry Smith /*@C 2654a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2664a9d489dSBarry Smith 2673f9fe445SBarry Smith Not Collective 2684a9d489dSBarry Smith 2694a9d489dSBarry Smith Input Parameters: 2704a9d489dSBarry Smith . coloring - the coloring context 2714a9d489dSBarry Smith 2724a9d489dSBarry Smith Output Parameters: 2734a9d489dSBarry Smith + f - the function 2744a9d489dSBarry Smith - fctx - the optional user-defined function context 2754a9d489dSBarry Smith 2764a9d489dSBarry Smith Level: intermediate 2774a9d489dSBarry Smith 2784a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 27995b89fc3SBarry Smith 28095b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 28195b89fc3SBarry Smith 2824a9d489dSBarry Smith @*/ 2837087cfbeSBarry Smith PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2844a9d489dSBarry Smith { 2854a9d489dSBarry Smith PetscFunctionBegin; 2860700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 2874a9d489dSBarry Smith if (f) *f = matfd->f; 2884a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2894a9d489dSBarry Smith PetscFunctionReturn(0); 2904a9d489dSBarry Smith } 2914a9d489dSBarry Smith 2924a9d489dSBarry Smith #undef __FUNCT__ 2934a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 294d64ed03dSBarry Smith /*@C 295005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 296005c665bSBarry Smith 2973f9fe445SBarry Smith Logically Collective on MatFDColoring 298fee21e36SBarry Smith 299ef5ee4d1SLois Curfman McInnes Input Parameters: 300ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 301ef5ee4d1SLois Curfman McInnes . f - the function 302ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 303ef5ee4d1SLois Curfman McInnes 3047850c7c0SBarry Smith Calling sequence of (*f) function: 3057850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 306ab637aeaSJed Brown If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 30715091d37SBarry Smith 3087850c7c0SBarry Smith Level: advanced 3097850c7c0SBarry Smith 310ab637aeaSJed Brown Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 3118d359177SBarry Smith SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by 312ab637aeaSJed Brown calling MatFDColoringApply() 3137850c7c0SBarry Smith 3147850c7c0SBarry Smith Fortran Notes: 3157850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 316ab637aeaSJed Brown be used without SNES or within the SNES solvers. 317f881d145SBarry Smith 318b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 31995b89fc3SBarry Smith 32095b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 32195b89fc3SBarry Smith 322005c665bSBarry Smith @*/ 3237087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 324005c665bSBarry Smith { 3253a40ed3dSBarry Smith PetscFunctionBegin; 3260700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 327005c665bSBarry Smith matfd->f = f; 328005c665bSBarry Smith matfd->fctx = fctx; 3293a40ed3dSBarry Smith PetscFunctionReturn(0); 330005c665bSBarry Smith } 331005c665bSBarry Smith 3324a2ae208SSatish Balay #undef __FUNCT__ 3334a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 334639f9d9dSBarry Smith /*@ 335b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 336639f9d9dSBarry Smith the options database. 337639f9d9dSBarry Smith 338fee21e36SBarry Smith Collective on MatFDColoring 339fee21e36SBarry Smith 34065f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 341ef5ee4d1SLois Curfman McInnes .vb 34265f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 343f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 344f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 345ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 346ef5ee4d1SLois Curfman McInnes .ve 347ef5ee4d1SLois Curfman McInnes 348ef5ee4d1SLois Curfman McInnes Input Parameter: 349ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 350ef5ee4d1SLois Curfman McInnes 351b4fc646aSLois Curfman McInnes Options Database Keys: 3525620d6dcSBarry Smith + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 353f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 3543ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 355ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 356e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info 357e350db43SBarry Smith - -mat_fd_coloring_view draw - Activates drawing 358639f9d9dSBarry Smith 35915091d37SBarry Smith Level: intermediate 36015091d37SBarry Smith 361b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 362d1c39f55SBarry Smith 363d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 364d1c39f55SBarry Smith 365639f9d9dSBarry Smith @*/ 3667087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 367639f9d9dSBarry Smith { 368dfbe8321SBarry Smith PetscErrorCode ierr; 369ace3abfcSBarry Smith PetscBool flg; 370efb30889SBarry Smith char value[3]; 3713a40ed3dSBarry Smith 3723a40ed3dSBarry Smith PetscFunctionBegin; 3730700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 374639f9d9dSBarry Smith 3753194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr); 37687828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 37787828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 3781bc75a8dSBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr); 379efb30889SBarry Smith if (flg) { 380efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 381efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 382e32f2f54SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 383efb30889SBarry Smith } 384f86b9fbaSHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr); 38593dfae19SHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg);CHKERRQ(ierr); 38693dfae19SHong Zhang if (flg && matfd->bcols > matfd->ncolors) { 38793dfae19SHong Zhang /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */ 38893dfae19SHong Zhang matfd->bcols = matfd->ncolors; 38993dfae19SHong Zhang } 390f86b9fbaSHong Zhang 3915d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 3920633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)matfd);CHKERRQ(ierr); 393b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3943a40ed3dSBarry Smith PetscFunctionReturn(0); 395005c665bSBarry Smith } 396005c665bSBarry Smith 3974a2ae208SSatish Balay #undef __FUNCT__ 39861ab5df0SBarry Smith #define __FUNCT__ "MatFDColoringSetType" 39961ab5df0SBarry Smith /*@C 40061ab5df0SBarry Smith MatFDColoringSetType - Sets the approach for computing the finite difference parameter 40161ab5df0SBarry Smith 40261ab5df0SBarry Smith Collective on MatFDColoring 40361ab5df0SBarry Smith 40461ab5df0SBarry Smith Input Parameters: 40561ab5df0SBarry Smith + coloring - the coloring context 40661ab5df0SBarry Smith - type - either MATMFFD_WP or MATMFFD_DS 40761ab5df0SBarry Smith 40861ab5df0SBarry Smith Options Database Keys: 40961ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds" 41061ab5df0SBarry Smith 41161ab5df0SBarry Smith Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries 41261ab5df0SBarry 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 41361ab5df0SBarry Smith introducing another one. 41461ab5df0SBarry Smith 41561ab5df0SBarry Smith Level: intermediate 41661ab5df0SBarry Smith 41761ab5df0SBarry Smith .keywords: Mat, finite differences, parameters 41861ab5df0SBarry Smith 41961ab5df0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 42061ab5df0SBarry Smith 42161ab5df0SBarry Smith @*/ 42261ab5df0SBarry Smith PetscErrorCode MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type) 42361ab5df0SBarry Smith { 42461ab5df0SBarry Smith PetscFunctionBegin; 42561ab5df0SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 42661ab5df0SBarry Smith /* 42761ab5df0SBarry Smith It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype 42861ab5df0SBarry Smith and this function is being provided as patch for a release so it shouldn't change the implementaton 42961ab5df0SBarry Smith */ 43061ab5df0SBarry Smith if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp"; 43161ab5df0SBarry Smith else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds"; 43261ab5df0SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type); 43361ab5df0SBarry Smith PetscFunctionReturn(0); 43461ab5df0SBarry Smith } 43561ab5df0SBarry Smith 43661ab5df0SBarry Smith #undef __FUNCT__ 437e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions" 438146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[]) 439005c665bSBarry Smith { 440dfbe8321SBarry Smith PetscErrorCode ierr; 441e350db43SBarry Smith PetscBool flg; 4423050cee2SBarry Smith PetscViewer viewer; 443cffb1e40SBarry Smith PetscViewerFormat format; 444005c665bSBarry Smith 4453a40ed3dSBarry Smith PetscFunctionBegin; 446146574abSBarry Smith if (prefix) { 447146574abSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 448146574abSBarry Smith } else { 449ce94432eSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 450146574abSBarry Smith } 451005c665bSBarry Smith if (flg) { 452cffb1e40SBarry Smith ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 4533050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 454cffb1e40SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 455cffb1e40SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 456005c665bSBarry Smith } 4573a40ed3dSBarry Smith PetscFunctionReturn(0); 458bbf0e169SBarry Smith } 459bbf0e169SBarry Smith 4604a2ae208SSatish Balay #undef __FUNCT__ 4614a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 46205869f15SSatish Balay /*@ 463639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 464639f9d9dSBarry Smith computation of Jacobians. 465bbf0e169SBarry Smith 466ef5ee4d1SLois Curfman McInnes Collective on Mat 467ef5ee4d1SLois Curfman McInnes 468639f9d9dSBarry Smith Input Parameters: 469ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 4707086a01eSPeter Brune - iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring() 471bbf0e169SBarry Smith 472bbf0e169SBarry Smith Output Parameter: 473639f9d9dSBarry Smith . color - the new coloring context 474bbf0e169SBarry Smith 47515091d37SBarry Smith Level: intermediate 47615091d37SBarry Smith 4778d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(), 478d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 4797086a01eSPeter Brune MatFDColoringView(), MatFDColoringSetParameters(), MatColoringCreate(), DMCreateColoring() 480bbf0e169SBarry Smith @*/ 4817087cfbeSBarry Smith PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 482bbf0e169SBarry Smith { 483639f9d9dSBarry Smith MatFDColoring c; 484639f9d9dSBarry Smith MPI_Comm comm; 485dfbe8321SBarry Smith PetscErrorCode ierr; 48638baddfdSBarry Smith PetscInt M,N; 487639f9d9dSBarry Smith 4883a40ed3dSBarry Smith PetscFunctionBegin; 489c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 490f141eedfSHong Zhang if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();"); 491d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 492639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 493ce94432eSBarry Smith if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices"); 494f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 49573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 496639f9d9dSBarry Smith 497b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 498b8f8c88eSHong Zhang 49993dfae19SHong Zhang if (mat->ops->fdcoloringcreate) { 50093dfae19SHong Zhang ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 50193dfae19SHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 50293dfae19SHong Zhang 5032a7a6963SBarry Smith ierr = MatCreateVecs(mat,NULL,&c->w1);CHKERRQ(ierr); 5043bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr); 505b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 5063bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr); 507b8f8c88eSHong Zhang 50877d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 50977d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 51049b058dcSBarry Smith c->currentcolor = -1; 511efb30889SBarry Smith c->htype = "wp"; 5124e269d77SPeter Brune c->fset = PETSC_FALSE; 513c8a9c622SHong Zhang c->setupcalled = PETSC_FALSE; 514639f9d9dSBarry Smith 515639f9d9dSBarry Smith *color = c; 5164e269d77SPeter Brune ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr); 517d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 5183a40ed3dSBarry Smith PetscFunctionReturn(0); 519bbf0e169SBarry Smith } 520bbf0e169SBarry Smith 5214a2ae208SSatish Balay #undef __FUNCT__ 5224a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 52305869f15SSatish Balay /*@ 524639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 525639f9d9dSBarry Smith via MatFDColoringCreate(). 526bbf0e169SBarry Smith 527ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 528ef5ee4d1SLois Curfman McInnes 529b4fc646aSLois Curfman McInnes Input Parameter: 530639f9d9dSBarry Smith . c - coloring context 531bbf0e169SBarry Smith 53215091d37SBarry Smith Level: intermediate 53315091d37SBarry Smith 534639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 535bbf0e169SBarry Smith @*/ 5366bf464f9SBarry Smith PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 537bbf0e169SBarry Smith { 5386849ba73SBarry Smith PetscErrorCode ierr; 53938baddfdSBarry Smith PetscInt i; 5400df34763SHong Zhang MatFDColoring color = *c; 541bbf0e169SBarry Smith 5423a40ed3dSBarry Smith PetscFunctionBegin; 5436bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 5440df34763SHong Zhang if (--((PetscObject)color)->refct > 0) {*c = 0; PetscFunctionReturn(0);} 545d4bb536fSBarry Smith 5460df34763SHong Zhang for (i=0; i<color->ncolors; i++) { 5470df34763SHong Zhang ierr = PetscFree(color->columns[i]);CHKERRQ(ierr); 548bbf0e169SBarry Smith } 5490df34763SHong Zhang ierr = PetscFree(color->ncolumns);CHKERRQ(ierr); 5500df34763SHong Zhang ierr = PetscFree(color->columns);CHKERRQ(ierr); 5510df34763SHong Zhang ierr = PetscFree(color->nrows);CHKERRQ(ierr); 5520df34763SHong Zhang if (color->htype[0] == 'w') { 5530df34763SHong Zhang ierr = PetscFree(color->matentry2);CHKERRQ(ierr); 5540df34763SHong Zhang } else { 5550df34763SHong Zhang ierr = PetscFree(color->matentry);CHKERRQ(ierr); 5560df34763SHong Zhang } 5570df34763SHong Zhang ierr = PetscFree(color->dy);CHKERRQ(ierr); 5580df34763SHong Zhang if (color->vscale) {ierr = VecDestroy(&color->vscale);CHKERRQ(ierr);} 5590df34763SHong Zhang ierr = VecDestroy(&color->w1);CHKERRQ(ierr); 5600df34763SHong Zhang ierr = VecDestroy(&color->w2);CHKERRQ(ierr); 5610df34763SHong Zhang ierr = VecDestroy(&color->w3);CHKERRQ(ierr); 562d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 5633a40ed3dSBarry Smith PetscFunctionReturn(0); 564bbf0e169SBarry Smith } 56543a90d84SBarry Smith 5664a2ae208SSatish Balay #undef __FUNCT__ 56749b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 56849b058dcSBarry Smith /*@C 56949b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 57049b058dcSBarry Smith that are currently being perturbed. 57149b058dcSBarry Smith 57249b058dcSBarry Smith Not Collective 57349b058dcSBarry Smith 57449b058dcSBarry Smith Input Parameters: 57549b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 57649b058dcSBarry Smith 57749b058dcSBarry Smith Output Parameters: 57849b058dcSBarry Smith + n - the number of local columns being perturbed 57949b058dcSBarry Smith - cols - the column indices, in global numbering 58049b058dcSBarry Smith 58149b058dcSBarry Smith Level: intermediate 58249b058dcSBarry Smith 58349b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 58449b058dcSBarry Smith 58549b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 58649b058dcSBarry Smith @*/ 5877087cfbeSBarry Smith PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 58849b058dcSBarry Smith { 58949b058dcSBarry Smith PetscFunctionBegin; 59049b058dcSBarry Smith if (coloring->currentcolor >= 0) { 59149b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 59249b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 59349b058dcSBarry Smith } else { 59449b058dcSBarry Smith *n = 0; 59549b058dcSBarry Smith } 59649b058dcSBarry Smith PetscFunctionReturn(0); 59749b058dcSBarry Smith } 59849b058dcSBarry Smith 59949b058dcSBarry Smith #undef __FUNCT__ 6004a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 60143a90d84SBarry Smith /*@ 602e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 603e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 60443a90d84SBarry Smith 605fee21e36SBarry Smith Collective on MatFDColoring 606fee21e36SBarry Smith 607ef5ee4d1SLois Curfman McInnes Input Parameters: 608ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 609ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 610ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 6117850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 612ef5ee4d1SLois Curfman McInnes 6138bba8e72SBarry Smith Options Database Keys: 614ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 615b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 616e350db43SBarry Smith . -mat_fd_coloring_view draw - Activates drawing of coloring 617e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 6188bba8e72SBarry Smith 61915091d37SBarry Smith Level: intermediate 62098d79826SSatish Balay 6217850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 62243a90d84SBarry Smith 62343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 62443a90d84SBarry Smith @*/ 625d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 62643a90d84SBarry Smith { 6273acb8795SBarry Smith PetscErrorCode ierr; 628684f2004SHong Zhang PetscBool flg = PETSC_FALSE; 6293acb8795SBarry Smith 6303acb8795SBarry Smith PetscFunctionBegin; 6310700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 6320700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 6330700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 634ce94432eSBarry Smith if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 635e32f2f54SBarry Smith if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 636c8a9c622SHong Zhang if (!coloring->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()"); 637684f2004SHong Zhang 638684f2004SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 639c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)coloring)->options,NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 640684f2004SHong Zhang if (flg) { 641684f2004SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 642684f2004SHong Zhang } else { 643684f2004SHong Zhang PetscBool assembled; 644684f2004SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 645684f2004SHong Zhang if (assembled) { 646684f2004SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 647684f2004SHong Zhang } 648684f2004SHong Zhang } 649684f2004SHong Zhang 6505922145eSHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 651d1e9a80fSBarry Smith ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,sctx);CHKERRQ(ierr); 6525922145eSHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 653*5bdb020cSBarry Smith if (!coloring->viewed) { 654*5bdb020cSBarry Smith ierr = MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view");CHKERRQ(ierr); 655*5bdb020cSBarry Smith coloring->viewed = PETSC_TRUE; 656*5bdb020cSBarry Smith } 6573acb8795SBarry Smith PetscFunctionReturn(0); 6583acb8795SBarry Smith } 659