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 7b45d2f2cSJed Brown #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 8bbf0e169SBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 103a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF" 117087cfbeSBarry Smith PetscErrorCode MatFDColoringSetF(MatFDColoring fd,Vec F) 123a7fca6bSBarry Smith { 134e269d77SPeter Brune PetscErrorCode ierr; 146e111a19SKarl Rupp 153a7fca6bSBarry Smith PetscFunctionBegin; 164e269d77SPeter Brune if (F) { 174e269d77SPeter Brune ierr = VecCopy(F,fd->w1);CHKERRQ(ierr); 184e269d77SPeter Brune fd->fset = PETSC_TRUE; 194e269d77SPeter Brune } else { 204e269d77SPeter Brune fd->fset = PETSC_FALSE; 214e269d77SPeter Brune } 223a7fca6bSBarry Smith PetscFunctionReturn(0); 233a7fca6bSBarry Smith } 243a7fca6bSBarry Smith 259804daf3SBarry Smith #include <petscdraw.h> 263a7fca6bSBarry Smith #undef __FUNCT__ 274a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 286849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 299194eea9SBarry Smith { 309194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 31dfbe8321SBarry Smith PetscErrorCode ierr; 32b312708bSHong Zhang PetscInt i,j,nz,row; 339194eea9SBarry Smith PetscReal x,y; 34b312708bSHong Zhang MatEntry *Jentry=fd->matentry; 359194eea9SBarry Smith 369194eea9SBarry Smith PetscFunctionBegin; 379194eea9SBarry Smith /* loop over colors */ 38b312708bSHong Zhang nz = 0; 399194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 409194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 41b312708bSHong Zhang row = Jentry[nz].row; 42b312708bSHong Zhang y = fd->M - row - fd->rstart; 43b312708bSHong Zhang x = (PetscReal)Jentry[nz++].col; 44b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 459194eea9SBarry Smith } 469194eea9SBarry Smith } 479194eea9SBarry Smith PetscFunctionReturn(0); 489194eea9SBarry Smith } 499194eea9SBarry Smith 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 526849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 53005c665bSBarry Smith { 54dfbe8321SBarry Smith PetscErrorCode ierr; 55ace3abfcSBarry Smith PetscBool isnull; 56b0a32e0cSBarry Smith PetscDraw draw; 5754d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 58005c665bSBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 60b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 61b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 629194eea9SBarry Smith 639194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 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); 68b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 69f77a5242SKarl Rupp ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL);CHKERRQ(ierr); 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71005c665bSBarry Smith } 72005c665bSBarry Smith 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 75bbf0e169SBarry Smith /*@C 76639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 77bbf0e169SBarry Smith 784c49b128SBarry Smith Collective on MatFDColoring 79fee21e36SBarry Smith 80ef5ee4d1SLois Curfman McInnes Input Parameters: 81ef5ee4d1SLois Curfman McInnes + c - the coloring context 82ef5ee4d1SLois Curfman McInnes - viewer - visualization context 83ef5ee4d1SLois Curfman McInnes 8415091d37SBarry Smith Level: intermediate 8515091d37SBarry Smith 86b4fc646aSLois Curfman McInnes Notes: 87b4fc646aSLois Curfman McInnes The available visualization contexts include 88b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 89b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 90ef5ee4d1SLois Curfman McInnes output where only the first processor opens 91ef5ee4d1SLois Curfman McInnes the file. All other processors send their 92ef5ee4d1SLois Curfman McInnes data to the first processor to print. 93b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 94bbf0e169SBarry Smith 959194eea9SBarry Smith Notes: 969194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 97b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 989194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 999194eea9SBarry Smith 100639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 101005c665bSBarry Smith 102b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 103bbf0e169SBarry Smith @*/ 1047087cfbeSBarry Smith PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 105bbf0e169SBarry Smith { 1066849ba73SBarry Smith PetscErrorCode ierr; 10738baddfdSBarry Smith PetscInt i,j; 108ace3abfcSBarry Smith PetscBool isdraw,iascii; 109f3ef73ceSBarry Smith PetscViewerFormat format; 110bbf0e169SBarry Smith 1113a40ed3dSBarry Smith PetscFunctionBegin; 1120700a824SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 1133050cee2SBarry Smith if (!viewer) { 114ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr); 1153050cee2SBarry Smith } 1160700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 117c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 118bbf0e169SBarry Smith 119251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 120251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1210f5bd95cSBarry Smith if (isdraw) { 122b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 12332077d6dSBarry Smith } else if (iascii) { 124dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer);CHKERRQ(ierr); 125a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 126a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 12777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 128ae09f205SBarry Smith 129b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 130fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 131b312708bSHong Zhang PetscInt row,col,nz; 132b312708bSHong Zhang nz = 0; 133b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 13477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 13577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 136b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 13777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 138639f9d9dSBarry Smith } 13977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 140b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 141b312708bSHong Zhang row = c->matentry[nz].row; 142b312708bSHong Zhang col = c->matentry[nz++].col; 143b312708bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",row,col);CHKERRQ(ierr); 144b4fc646aSLois Curfman McInnes } 145bbf0e169SBarry Smith } 146bbf0e169SBarry Smith } 147b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 148bbf0e169SBarry Smith } 1493a40ed3dSBarry Smith PetscFunctionReturn(0); 150639f9d9dSBarry Smith } 151639f9d9dSBarry Smith 1524a2ae208SSatish Balay #undef __FUNCT__ 1534a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 154639f9d9dSBarry Smith /*@ 155b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 156b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 157639f9d9dSBarry Smith 1583f9fe445SBarry Smith Logically Collective on MatFDColoring 159ef5ee4d1SLois Curfman McInnes 160ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 161ef5ee4d1SLois Curfman McInnes .vb 16265f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 163f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 164f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 165ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 166ef5ee4d1SLois Curfman McInnes .ve 167639f9d9dSBarry Smith 168639f9d9dSBarry Smith Input Parameters: 169ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 170639f9d9dSBarry Smith . error_rel - relative error 171f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 172fee21e36SBarry Smith 17315091d37SBarry Smith Level: advanced 17415091d37SBarry Smith 175b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 176b4fc646aSLois Curfman McInnes 17795b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 17895b89fc3SBarry Smith 179639f9d9dSBarry Smith @*/ 1807087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 181639f9d9dSBarry Smith { 1823a40ed3dSBarry Smith PetscFunctionBegin; 1830700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 184c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,error,2); 185c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,umin,3); 186639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 187639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1883a40ed3dSBarry Smith PetscFunctionReturn(0); 189639f9d9dSBarry Smith } 190639f9d9dSBarry Smith 191*f86b9fbaSHong Zhang #undef __FUNCT__ 192*f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetBlockSize" 193*f86b9fbaSHong Zhang /*@ 194*f86b9fbaSHong Zhang MatFDColoringSetBlockSize - Sets block size for efficient 195*f86b9fbaSHong Zhang inserting entries of Jacobian matrix. 196005c665bSBarry Smith 197*f86b9fbaSHong Zhang Logically Collective on MatFDColoring 198*f86b9fbaSHong Zhang 199*f86b9fbaSHong Zhang Input Parameters: 200*f86b9fbaSHong Zhang + coloring - the coloring context 201*f86b9fbaSHong Zhang . brows - number of rows in the block 202*f86b9fbaSHong Zhang - bcols - number of columns in the block 203*f86b9fbaSHong Zhang 204*f86b9fbaSHong Zhang Level: intermediate 205*f86b9fbaSHong Zhang 206*f86b9fbaSHong Zhang .keywords: Mat, coloring 207*f86b9fbaSHong Zhang 208*f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 209*f86b9fbaSHong Zhang 210*f86b9fbaSHong Zhang @*/ 211*f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols) 212*f86b9fbaSHong Zhang { 213*f86b9fbaSHong Zhang PetscFunctionBegin; 214*f86b9fbaSHong Zhang PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 215*f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,brows,2); 216*f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,bcols,3); 217*f86b9fbaSHong Zhang if (brows != PETSC_DEFAULT) matfd->brows = brows; 218*f86b9fbaSHong Zhang if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 219*f86b9fbaSHong Zhang PetscFunctionReturn(0); 220*f86b9fbaSHong Zhang } 221*f86b9fbaSHong Zhang 222*f86b9fbaSHong Zhang #undef __FUNCT__ 223*f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp" 224*f86b9fbaSHong Zhang /*@ 225*f86b9fbaSHong Zhang MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 226*f86b9fbaSHong Zhang 227*f86b9fbaSHong Zhang Collective on Mat 228*f86b9fbaSHong Zhang 229*f86b9fbaSHong Zhang Input Parameters: 230*f86b9fbaSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 231*f86b9fbaSHong Zhang . iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 232*f86b9fbaSHong Zhang - color - the matrix coloring context 233*f86b9fbaSHong Zhang 234*f86b9fbaSHong Zhang Level: beginner 235*f86b9fbaSHong Zhang 236*f86b9fbaSHong Zhang .keywords: MatFDColoring, setup 237*f86b9fbaSHong Zhang 238*f86b9fbaSHong Zhang .seealso: MatFDColoringCreate(), MatFDColoringDestroy() 239*f86b9fbaSHong Zhang @*/ 240*f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color) 241*f86b9fbaSHong Zhang { 242*f86b9fbaSHong Zhang PetscErrorCode ierr; 243*f86b9fbaSHong Zhang 244*f86b9fbaSHong Zhang PetscFunctionBegin; 245*f86b9fbaSHong Zhang if (mat->ops->fdcoloringsetup) { 246*f86b9fbaSHong Zhang ierr = (*mat->ops->fdcoloringsetup)(mat,iscoloring,color);CHKERRQ(ierr); 247*f86b9fbaSHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name); 248*f86b9fbaSHong Zhang PetscFunctionReturn(0); 249*f86b9fbaSHong Zhang } 250ff0cfa39SBarry Smith 2514a2ae208SSatish Balay #undef __FUNCT__ 2524a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction" 2534a9d489dSBarry Smith /*@C 2544a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2554a9d489dSBarry Smith 2563f9fe445SBarry Smith Not Collective 2574a9d489dSBarry Smith 2584a9d489dSBarry Smith Input Parameters: 2594a9d489dSBarry Smith . coloring - the coloring context 2604a9d489dSBarry Smith 2614a9d489dSBarry Smith Output Parameters: 2624a9d489dSBarry Smith + f - the function 2634a9d489dSBarry Smith - fctx - the optional user-defined function context 2644a9d489dSBarry Smith 2654a9d489dSBarry Smith Level: intermediate 2664a9d489dSBarry Smith 2674a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 26895b89fc3SBarry Smith 26995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 27095b89fc3SBarry Smith 2714a9d489dSBarry Smith @*/ 2727087cfbeSBarry Smith PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2734a9d489dSBarry Smith { 2744a9d489dSBarry Smith PetscFunctionBegin; 2750700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 2764a9d489dSBarry Smith if (f) *f = matfd->f; 2774a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2784a9d489dSBarry Smith PetscFunctionReturn(0); 2794a9d489dSBarry Smith } 2804a9d489dSBarry Smith 2814a9d489dSBarry Smith #undef __FUNCT__ 2824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 283d64ed03dSBarry Smith /*@C 284005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 285005c665bSBarry Smith 2863f9fe445SBarry Smith Logically Collective on MatFDColoring 287fee21e36SBarry Smith 288ef5ee4d1SLois Curfman McInnes Input Parameters: 289ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 290ef5ee4d1SLois Curfman McInnes . f - the function 291ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 292ef5ee4d1SLois Curfman McInnes 2937850c7c0SBarry Smith Calling sequence of (*f) function: 2947850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 295ab637aeaSJed Brown If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 29615091d37SBarry Smith 2977850c7c0SBarry Smith Level: advanced 2987850c7c0SBarry Smith 299ab637aeaSJed Brown Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 3008d359177SBarry Smith SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by 301ab637aeaSJed Brown calling MatFDColoringApply() 3027850c7c0SBarry Smith 3037850c7c0SBarry Smith Fortran Notes: 3047850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 305ab637aeaSJed Brown be used without SNES or within the SNES solvers. 306f881d145SBarry Smith 307b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 30895b89fc3SBarry Smith 30995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 31095b89fc3SBarry Smith 311005c665bSBarry Smith @*/ 3127087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 313005c665bSBarry Smith { 3143a40ed3dSBarry Smith PetscFunctionBegin; 3150700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 316005c665bSBarry Smith matfd->f = f; 317005c665bSBarry Smith matfd->fctx = fctx; 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 319005c665bSBarry Smith } 320005c665bSBarry Smith 3214a2ae208SSatish Balay #undef __FUNCT__ 3224a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 323639f9d9dSBarry Smith /*@ 324b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 325639f9d9dSBarry Smith the options database. 326639f9d9dSBarry Smith 327fee21e36SBarry Smith Collective on MatFDColoring 328fee21e36SBarry Smith 32965f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 330ef5ee4d1SLois Curfman McInnes .vb 33165f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 332f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 333f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 334ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 335ef5ee4d1SLois Curfman McInnes .ve 336ef5ee4d1SLois Curfman McInnes 337ef5ee4d1SLois Curfman McInnes Input Parameter: 338ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 339ef5ee4d1SLois Curfman McInnes 340b4fc646aSLois Curfman McInnes Options Database Keys: 341d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 342ef5ee4d1SLois Curfman McInnes of relative error in the function) 343f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 3443ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 345ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 346e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info 347e350db43SBarry Smith - -mat_fd_coloring_view draw - Activates drawing 348639f9d9dSBarry Smith 34915091d37SBarry Smith Level: intermediate 35015091d37SBarry Smith 351b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 352d1c39f55SBarry Smith 353d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 354d1c39f55SBarry Smith 355639f9d9dSBarry Smith @*/ 3567087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 357639f9d9dSBarry Smith { 358dfbe8321SBarry Smith PetscErrorCode ierr; 359ace3abfcSBarry Smith PetscBool flg; 360efb30889SBarry Smith char value[3]; 3613a40ed3dSBarry Smith 3623a40ed3dSBarry Smith PetscFunctionBegin; 3630700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 364639f9d9dSBarry Smith 3653194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr); 36687828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 36787828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 3681bc75a8dSBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr); 369efb30889SBarry Smith if (flg) { 370efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 371efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 372e32f2f54SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 373efb30889SBarry Smith } 374*f86b9fbaSHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL);CHKERRQ(ierr); 375*f86b9fbaSHong Zhang ierr = PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,NULL);CHKERRQ(ierr); 376*f86b9fbaSHong Zhang 3775d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 3785d973c19SBarry Smith ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr); 379b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3803a40ed3dSBarry Smith PetscFunctionReturn(0); 381005c665bSBarry Smith } 382005c665bSBarry Smith 3834a2ae208SSatish Balay #undef __FUNCT__ 384e350db43SBarry Smith #define __FUNCT__ "MatFDColoringViewFromOptions" 385146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[]) 386005c665bSBarry Smith { 387dfbe8321SBarry Smith PetscErrorCode ierr; 388e350db43SBarry Smith PetscBool flg; 3893050cee2SBarry Smith PetscViewer viewer; 390cffb1e40SBarry Smith PetscViewerFormat format; 391005c665bSBarry Smith 3923a40ed3dSBarry Smith PetscFunctionBegin; 393146574abSBarry Smith if (prefix) { 394146574abSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 395146574abSBarry Smith } else { 396ce94432eSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg);CHKERRQ(ierr); 397146574abSBarry Smith } 398005c665bSBarry Smith if (flg) { 399cffb1e40SBarry Smith ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 4003050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 401cffb1e40SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 402cffb1e40SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 403005c665bSBarry Smith } 4043a40ed3dSBarry Smith PetscFunctionReturn(0); 405bbf0e169SBarry Smith } 406bbf0e169SBarry Smith 4074a2ae208SSatish Balay #undef __FUNCT__ 4084a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 40905869f15SSatish Balay /*@ 410639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 411639f9d9dSBarry Smith computation of Jacobians. 412bbf0e169SBarry Smith 413ef5ee4d1SLois Curfman McInnes Collective on Mat 414ef5ee4d1SLois Curfman McInnes 415639f9d9dSBarry Smith Input Parameters: 416ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 417e727c939SJed Brown - iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 418bbf0e169SBarry Smith 419bbf0e169SBarry Smith Output Parameter: 420639f9d9dSBarry Smith . color - the new coloring context 421bbf0e169SBarry Smith 42215091d37SBarry Smith Level: intermediate 42315091d37SBarry Smith 4248d359177SBarry Smith .seealso: MatFDColoringDestroy(),SNESComputeJacobianDefaultColor(), ISColoringCreate(), 425d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 426e727c939SJed Brown MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring() 427bbf0e169SBarry Smith @*/ 4287087cfbeSBarry Smith PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 429bbf0e169SBarry Smith { 430639f9d9dSBarry Smith MatFDColoring c; 431639f9d9dSBarry Smith MPI_Comm comm; 432dfbe8321SBarry Smith PetscErrorCode ierr; 43338baddfdSBarry Smith PetscInt M,N; 434639f9d9dSBarry Smith 4353a40ed3dSBarry Smith PetscFunctionBegin; 436f141eedfSHong Zhang if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();"); 437d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 438639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 439ce94432eSBarry Smith if (M != N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices"); 440f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 44167c2884eSBarry Smith ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 442639f9d9dSBarry Smith 443b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 444b8f8c88eSHong Zhang 445f77a5242SKarl Rupp ierr = MatGetVecs(mat,NULL,&c->w1);CHKERRQ(ierr); 4463bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1);CHKERRQ(ierr); 447b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 4483bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2);CHKERRQ(ierr); 449b8f8c88eSHong Zhang 45077d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 45177d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 45249b058dcSBarry Smith c->currentcolor = -1; 453efb30889SBarry Smith c->htype = "wp"; 4544e269d77SPeter Brune c->fset = PETSC_FALSE; 455639f9d9dSBarry Smith 456639f9d9dSBarry Smith *color = c; 4574e269d77SPeter Brune ierr = PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c);CHKERRQ(ierr); 458d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4593a40ed3dSBarry Smith PetscFunctionReturn(0); 460bbf0e169SBarry Smith } 461bbf0e169SBarry Smith 4624a2ae208SSatish Balay #undef __FUNCT__ 4634a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 46405869f15SSatish Balay /*@ 465639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 466639f9d9dSBarry Smith via MatFDColoringCreate(). 467bbf0e169SBarry Smith 468ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 469ef5ee4d1SLois Curfman McInnes 470b4fc646aSLois Curfman McInnes Input Parameter: 471639f9d9dSBarry Smith . c - coloring context 472bbf0e169SBarry Smith 47315091d37SBarry Smith Level: intermediate 47415091d37SBarry Smith 475639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 476bbf0e169SBarry Smith @*/ 4776bf464f9SBarry Smith PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 478bbf0e169SBarry Smith { 4796849ba73SBarry Smith PetscErrorCode ierr; 48038baddfdSBarry Smith PetscInt i; 481bbf0e169SBarry Smith 4823a40ed3dSBarry Smith PetscFunctionBegin; 4836bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 4846bf464f9SBarry Smith if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);} 485d4bb536fSBarry Smith 4866bf464f9SBarry Smith for (i=0; i<(*c)->ncolors; i++) { 4876bf464f9SBarry Smith ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr); 488bbf0e169SBarry Smith } 4896bf464f9SBarry Smith ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr); 4906bf464f9SBarry Smith ierr = PetscFree((*c)->columns);CHKERRQ(ierr); 4916bf464f9SBarry Smith ierr = PetscFree((*c)->nrows);CHKERRQ(ierr); 492573f477fSHong Zhang ierr = PetscFree((*c)->matentry);CHKERRQ(ierr); 493b7799381SHong Zhang ierr = PetscFree((*c)->dy);CHKERRQ(ierr); 4949e917edbSHong Zhang if ((*c)->vscale) {ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);} 4956bf464f9SBarry Smith ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr); 4966bf464f9SBarry Smith ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr); 4976bf464f9SBarry Smith ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr); 498d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 4993a40ed3dSBarry Smith PetscFunctionReturn(0); 500bbf0e169SBarry Smith } 50143a90d84SBarry Smith 5024a2ae208SSatish Balay #undef __FUNCT__ 50349b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 50449b058dcSBarry Smith /*@C 50549b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 50649b058dcSBarry Smith that are currently being perturbed. 50749b058dcSBarry Smith 50849b058dcSBarry Smith Not Collective 50949b058dcSBarry Smith 51049b058dcSBarry Smith Input Parameters: 51149b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 51249b058dcSBarry Smith 51349b058dcSBarry Smith Output Parameters: 51449b058dcSBarry Smith + n - the number of local columns being perturbed 51549b058dcSBarry Smith - cols - the column indices, in global numbering 51649b058dcSBarry Smith 51749b058dcSBarry Smith Level: intermediate 51849b058dcSBarry Smith 51949b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 52049b058dcSBarry Smith 52149b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 52249b058dcSBarry Smith @*/ 5237087cfbeSBarry Smith PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 52449b058dcSBarry Smith { 52549b058dcSBarry Smith PetscFunctionBegin; 52649b058dcSBarry Smith if (coloring->currentcolor >= 0) { 52749b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 52849b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 52949b058dcSBarry Smith } else { 53049b058dcSBarry Smith *n = 0; 53149b058dcSBarry Smith } 53249b058dcSBarry Smith PetscFunctionReturn(0); 53349b058dcSBarry Smith } 53449b058dcSBarry Smith 53549b058dcSBarry Smith #undef __FUNCT__ 5364a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 53743a90d84SBarry Smith /*@ 538e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 539e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 54043a90d84SBarry Smith 541fee21e36SBarry Smith Collective on MatFDColoring 542fee21e36SBarry Smith 543ef5ee4d1SLois Curfman McInnes Input Parameters: 544ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 545ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 546ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 5477850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 548ef5ee4d1SLois Curfman McInnes 5498bba8e72SBarry Smith Options Database Keys: 550ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 551b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 552e350db43SBarry Smith . -mat_fd_coloring_view draw - Activates drawing of coloring 553e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 5548bba8e72SBarry Smith 55515091d37SBarry Smith Level: intermediate 55698d79826SSatish Balay 5577850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 55843a90d84SBarry Smith 55943a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 56043a90d84SBarry Smith @*/ 5617087cfbeSBarry Smith PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 56243a90d84SBarry Smith { 5633acb8795SBarry Smith PetscErrorCode ierr; 564684f2004SHong Zhang PetscBool flg = PETSC_FALSE; 5653acb8795SBarry Smith 5663acb8795SBarry Smith PetscFunctionBegin; 5670700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 5680700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 5690700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 570ce94432eSBarry Smith if (!coloring->f) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 571e32f2f54SBarry Smith if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 572684f2004SHong Zhang 573684f2004SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 574684f2004SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 575684f2004SHong Zhang if (flg) { 576684f2004SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 577684f2004SHong Zhang } else { 578684f2004SHong Zhang PetscBool assembled; 579684f2004SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 580684f2004SHong Zhang if (assembled) { 581684f2004SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 582684f2004SHong Zhang } 583684f2004SHong Zhang } 584684f2004SHong Zhang 5855922145eSHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 5863acb8795SBarry Smith ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 5875922145eSHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 5883acb8795SBarry Smith PetscFunctionReturn(0); 5893acb8795SBarry Smith } 590