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 { 123a7fca6bSBarry Smith PetscFunctionBegin; 134e269d77SPeter Brune if (F) { 149566063dSJacob Faibussowitsch PetscCall(VecCopy(F,fd->w1)); 154e269d77SPeter Brune fd->fset = PETSC_TRUE; 164e269d77SPeter Brune } else { 174e269d77SPeter Brune fd->fset = PETSC_FALSE; 184e269d77SPeter Brune } 193a7fca6bSBarry Smith PetscFunctionReturn(0); 203a7fca6bSBarry Smith } 213a7fca6bSBarry Smith 229804daf3SBarry Smith #include <petscdraw.h> 236849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 249194eea9SBarry Smith { 259194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 26b312708bSHong Zhang PetscInt i,j,nz,row; 279194eea9SBarry Smith PetscReal x,y; 28b312708bSHong Zhang MatEntry *Jentry=fd->matentry; 299194eea9SBarry Smith 309194eea9SBarry Smith PetscFunctionBegin; 319194eea9SBarry Smith /* loop over colors */ 32b312708bSHong Zhang nz = 0; 339194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 349194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 35b312708bSHong Zhang row = Jentry[nz].row; 36b312708bSHong Zhang y = fd->M - row - fd->rstart; 37b312708bSHong Zhang x = (PetscReal)Jentry[nz++].col; 389566063dSJacob Faibussowitsch PetscCall(PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1)); 399194eea9SBarry Smith } 409194eea9SBarry Smith } 419194eea9SBarry Smith PetscFunctionReturn(0); 429194eea9SBarry Smith } 439194eea9SBarry Smith 446849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 45005c665bSBarry Smith { 46ace3abfcSBarry Smith PetscBool isnull; 47b0a32e0cSBarry Smith PetscDraw draw; 4854d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 49005c665bSBarry Smith 503a40ed3dSBarry Smith PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 529566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw,&isnull)); 53832b7cebSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 54005c665bSBarry Smith 55005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 56005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 579566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw,xl,yl,xr,yr)); 589566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer)); 599566063dSJacob Faibussowitsch PetscCall(PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd)); 609566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)fd,"Zoomviewer",NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 623a40ed3dSBarry Smith PetscFunctionReturn(0); 63005c665bSBarry Smith } 64005c665bSBarry Smith 65bbf0e169SBarry Smith /*@C 66639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 67bbf0e169SBarry Smith 684c49b128SBarry Smith Collective on MatFDColoring 69fee21e36SBarry Smith 70ef5ee4d1SLois Curfman McInnes Input Parameters: 71ef5ee4d1SLois Curfman McInnes + c - the coloring context 72ef5ee4d1SLois Curfman McInnes - viewer - visualization context 73ef5ee4d1SLois Curfman McInnes 7415091d37SBarry Smith Level: intermediate 7515091d37SBarry Smith 76b4fc646aSLois Curfman McInnes Notes: 77b4fc646aSLois Curfman McInnes The available visualization contexts include 78b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 79b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 80ef5ee4d1SLois Curfman McInnes output where only the first processor opens 81ef5ee4d1SLois Curfman McInnes the file. All other processors send their 82ef5ee4d1SLois Curfman McInnes data to the first processor to print. 83b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 84bbf0e169SBarry Smith 859194eea9SBarry Smith Notes: 869194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 87b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 889194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 899194eea9SBarry Smith 90db781477SPatrick Sanan .seealso: `MatFDColoringCreate()` 91005c665bSBarry Smith 92bbf0e169SBarry Smith @*/ 937087cfbeSBarry Smith PetscErrorCode MatFDColoringView(MatFDColoring c,PetscViewer viewer) 94bbf0e169SBarry Smith { 9538baddfdSBarry Smith PetscInt i,j; 96ace3abfcSBarry Smith PetscBool isdraw,iascii; 97f3ef73ceSBarry Smith PetscViewerFormat format; 98bbf0e169SBarry Smith 993a40ed3dSBarry Smith PetscFunctionBegin; 1000700a824SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 1013050cee2SBarry Smith if (!viewer) { 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer)); 1033050cee2SBarry Smith } 1040700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 105c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 106bbf0e169SBarry Smith 1079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 1089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 1090f5bd95cSBarry Smith if (isdraw) { 1109566063dSJacob Faibussowitsch PetscCall(MatFDColoringView_Draw(c,viewer)); 11132077d6dSBarry Smith } else if (iascii) { 1129566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",(double)c->error_rel)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Umin=%g\n",(double)c->umin)); 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of colors=%" PetscInt_FMT "\n",c->ncolors)); 116ae09f205SBarry Smith 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 118fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 119b312708bSHong Zhang PetscInt row,col,nz; 120b312708bSHong Zhang nz = 0; 121b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Information for color %" PetscInt_FMT "\n",i)); 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of columns %" PetscInt_FMT "\n",c->ncolumns[i])); 124b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT "\n",c->columns[i][j])); 126639f9d9dSBarry Smith } 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of rows %" PetscInt_FMT "\n",c->nrows[i])); 1285bdb020cSBarry Smith if (c->matentry) { 129b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 130b312708bSHong Zhang row = c->matentry[nz].row; 131b312708bSHong Zhang col = c->matentry[nz++].col; 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT " %" PetscInt_FMT " \n",row,col)); 133b4fc646aSLois Curfman McInnes } 134bbf0e169SBarry Smith } 135bbf0e169SBarry Smith } 1365bdb020cSBarry Smith } 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 138bbf0e169SBarry Smith } 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 140639f9d9dSBarry Smith } 141639f9d9dSBarry Smith 142639f9d9dSBarry Smith /*@ 143b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 144b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 145639f9d9dSBarry Smith 1463f9fe445SBarry Smith Logically Collective on MatFDColoring 147ef5ee4d1SLois Curfman McInnes 148ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 149ef5ee4d1SLois Curfman McInnes .vb 15065f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 151d461c19dSHong Zhang htype = 'ds': 152f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 153f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 154ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 155d461c19dSHong Zhang 156d461c19dSHong Zhang htype = 'wp': 157d461c19dSHong Zhang h = error_rel * sqrt(1 + ||u||) 158ef5ee4d1SLois Curfman McInnes .ve 159639f9d9dSBarry Smith 160639f9d9dSBarry Smith Input Parameters: 161ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 162639f9d9dSBarry Smith . error_rel - relative error 163f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 164fee21e36SBarry Smith 16515091d37SBarry Smith Level: advanced 16615091d37SBarry Smith 167db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()` 16895b89fc3SBarry Smith 169639f9d9dSBarry Smith @*/ 1707087cfbeSBarry Smith PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 171639f9d9dSBarry Smith { 1723a40ed3dSBarry Smith PetscFunctionBegin; 1730700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 174c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,error,2); 175c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd,umin,3); 176639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 177639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1783a40ed3dSBarry Smith PetscFunctionReturn(0); 179639f9d9dSBarry Smith } 180639f9d9dSBarry Smith 181f86b9fbaSHong Zhang /*@ 182c8a9c622SHong Zhang MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix. 183005c665bSBarry Smith 184f86b9fbaSHong Zhang Logically Collective on MatFDColoring 185f86b9fbaSHong Zhang 186f86b9fbaSHong Zhang Input Parameters: 187f86b9fbaSHong Zhang + coloring - the coloring context 188f86b9fbaSHong Zhang . brows - number of rows in the block 189f86b9fbaSHong Zhang - bcols - number of columns in the block 190f86b9fbaSHong Zhang 191f86b9fbaSHong Zhang Level: intermediate 192f86b9fbaSHong Zhang 193db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()` 194f86b9fbaSHong Zhang 195f86b9fbaSHong Zhang @*/ 196f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd,PetscInt brows,PetscInt bcols) 197f86b9fbaSHong Zhang { 198f86b9fbaSHong Zhang PetscFunctionBegin; 199f86b9fbaSHong Zhang PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 200f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,brows,2); 201f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd,bcols,3); 202f86b9fbaSHong Zhang if (brows != PETSC_DEFAULT) matfd->brows = brows; 203f86b9fbaSHong Zhang if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 204f86b9fbaSHong Zhang PetscFunctionReturn(0); 205f86b9fbaSHong Zhang } 206f86b9fbaSHong Zhang 207f86b9fbaSHong Zhang /*@ 208f86b9fbaSHong Zhang MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 209f86b9fbaSHong Zhang 210f86b9fbaSHong Zhang Collective on Mat 211f86b9fbaSHong Zhang 212f86b9fbaSHong Zhang Input Parameters: 213f86b9fbaSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 214f86b9fbaSHong Zhang . iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring() 215f86b9fbaSHong Zhang - color - the matrix coloring context 216f86b9fbaSHong Zhang 217f86b9fbaSHong Zhang Level: beginner 218f86b9fbaSHong Zhang 219bdaf1daeSBarry Smith Notes: When the coloring type is IS_COLORING_LOCAL the coloring is in the local ordering of the unknowns. 220bdaf1daeSBarry Smith 221db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()` 222f86b9fbaSHong Zhang @*/ 223f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp(Mat mat,ISColoring iscoloring,MatFDColoring color) 224f86b9fbaSHong Zhang { 225bdaf1daeSBarry Smith PetscBool eq; 226f86b9fbaSHong Zhang 227f86b9fbaSHong Zhang PetscFunctionBegin; 228c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 229c8a9c622SHong Zhang PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,3); 230c8a9c622SHong Zhang if (color->setupcalled) PetscFunctionReturn(0); 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)mat,color->matid,&eq)); 23228b400f6SJacob Faibussowitsch PetscCheck(eq,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()"); 233c8a9c622SHong Zhang 2349566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp,mat,0,0,0)); 235*dbbe0bcdSBarry Smith PetscUseTypeMethod(mat,fdcoloringsetup ,iscoloring,color); 236c8a9c622SHong Zhang 237c8a9c622SHong Zhang color->setupcalled = PETSC_TRUE; 2389566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp,mat,0,0,0)); 239f86b9fbaSHong Zhang PetscFunctionReturn(0); 240f86b9fbaSHong Zhang } 241ff0cfa39SBarry Smith 2424a9d489dSBarry Smith /*@C 2434a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2444a9d489dSBarry Smith 2453f9fe445SBarry Smith Not Collective 2464a9d489dSBarry Smith 247f899ff85SJose E. Roman Input Parameter: 2484a9d489dSBarry Smith . coloring - the coloring context 2494a9d489dSBarry Smith 2504a9d489dSBarry Smith Output Parameters: 2514a9d489dSBarry Smith + f - the function 2524a9d489dSBarry Smith - fctx - the optional user-defined function context 2534a9d489dSBarry Smith 2544a9d489dSBarry Smith Level: intermediate 2554a9d489dSBarry Smith 256db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()` 25795b89fc3SBarry Smith 2584a9d489dSBarry Smith @*/ 2597087cfbeSBarry Smith PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2604a9d489dSBarry Smith { 2614a9d489dSBarry Smith PetscFunctionBegin; 2620700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 2634a9d489dSBarry Smith if (f) *f = matfd->f; 2644a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2654a9d489dSBarry Smith PetscFunctionReturn(0); 2664a9d489dSBarry Smith } 2674a9d489dSBarry Smith 268d64ed03dSBarry Smith /*@C 269005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 270005c665bSBarry Smith 2713f9fe445SBarry Smith Logically Collective on MatFDColoring 272fee21e36SBarry Smith 273ef5ee4d1SLois Curfman McInnes Input Parameters: 274ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 275ef5ee4d1SLois Curfman McInnes . f - the function 276ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 277ef5ee4d1SLois Curfman McInnes 2787850c7c0SBarry Smith Calling sequence of (*f) function: 2797850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 280ab637aeaSJed Brown If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 28115091d37SBarry Smith 2827850c7c0SBarry Smith Level: advanced 2837850c7c0SBarry Smith 28495452b02SPatrick Sanan Notes: 28595452b02SPatrick Sanan This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument 2868d359177SBarry Smith SNESComputeJacobianDefaultColor()) and only needs to be used by someone computing a matrix via coloring directly by 287ab637aeaSJed Brown calling MatFDColoringApply() 2887850c7c0SBarry Smith 2897850c7c0SBarry Smith Fortran Notes: 2907850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 291ab637aeaSJed Brown be used without SNES or within the SNES solvers. 292f881d145SBarry Smith 293db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()` 29495b89fc3SBarry Smith 295005c665bSBarry Smith @*/ 2967087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 297005c665bSBarry Smith { 2983a40ed3dSBarry Smith PetscFunctionBegin; 2990700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 300005c665bSBarry Smith matfd->f = f; 301005c665bSBarry Smith matfd->fctx = fctx; 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 303005c665bSBarry Smith } 304005c665bSBarry Smith 305639f9d9dSBarry Smith /*@ 306b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 307639f9d9dSBarry Smith the options database. 308639f9d9dSBarry Smith 309fee21e36SBarry Smith Collective on MatFDColoring 310fee21e36SBarry Smith 31165f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 312ef5ee4d1SLois Curfman McInnes .vb 31365f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 314f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 315f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 316ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 317ef5ee4d1SLois Curfman McInnes .ve 318ef5ee4d1SLois Curfman McInnes 319ef5ee4d1SLois Curfman McInnes Input Parameter: 320ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 321ef5ee4d1SLois Curfman McInnes 322b4fc646aSLois Curfman McInnes Options Database Keys: 3235620d6dcSBarry Smith + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 324f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 3253ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 326ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 327e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info 328e350db43SBarry Smith - -mat_fd_coloring_view draw - Activates drawing 329639f9d9dSBarry Smith 33015091d37SBarry Smith Level: intermediate 33115091d37SBarry Smith 332db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 333d1c39f55SBarry Smith 334639f9d9dSBarry Smith @*/ 3357087cfbeSBarry Smith PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) 336639f9d9dSBarry Smith { 337ace3abfcSBarry Smith PetscBool flg; 338efb30889SBarry Smith char value[3]; 3393a40ed3dSBarry Smith 3403a40ed3dSBarry Smith PetscFunctionBegin; 3410700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 342639f9d9dSBarry Smith 343d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)matfd); 3449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,NULL)); 3459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,NULL)); 3469566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,sizeof(value),&flg)); 347efb30889SBarry Smith if (flg) { 348efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 349efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 35098921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 351efb30889SBarry Smith } 3529566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_fd_coloring_brows","Number of block rows","MatFDColoringSetBlockSize",matfd->brows,&matfd->brows,NULL)); 3539566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols","Number of block columns","MatFDColoringSetBlockSize",matfd->bcols,&matfd->bcols,&flg)); 35493dfae19SHong Zhang if (flg && matfd->bcols > matfd->ncolors) { 35593dfae19SHong Zhang /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */ 35693dfae19SHong Zhang matfd->bcols = matfd->ncolors; 35793dfae19SHong Zhang } 358f86b9fbaSHong Zhang 3595d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 360*dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd,PetscOptionsObject)); 361d0609cedSBarry Smith PetscOptionsEnd(); 3623a40ed3dSBarry Smith PetscFunctionReturn(0); 363005c665bSBarry Smith } 364005c665bSBarry Smith 36561ab5df0SBarry Smith /*@C 36661ab5df0SBarry Smith MatFDColoringSetType - Sets the approach for computing the finite difference parameter 36761ab5df0SBarry Smith 36861ab5df0SBarry Smith Collective on MatFDColoring 36961ab5df0SBarry Smith 37061ab5df0SBarry Smith Input Parameters: 37161ab5df0SBarry Smith + coloring - the coloring context 37261ab5df0SBarry Smith - type - either MATMFFD_WP or MATMFFD_DS 37361ab5df0SBarry Smith 37461ab5df0SBarry Smith Options Database Keys: 37561ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds" 37661ab5df0SBarry Smith 37761ab5df0SBarry Smith Note: It is goofy that the argument type is MatMFFDType since the MatFDColoring actually computes the matrix entries 37861ab5df0SBarry 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 37961ab5df0SBarry Smith introducing another one. 38061ab5df0SBarry Smith 38161ab5df0SBarry Smith Level: intermediate 38261ab5df0SBarry Smith 383db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 38461ab5df0SBarry Smith 38561ab5df0SBarry Smith @*/ 38661ab5df0SBarry Smith PetscErrorCode MatFDColoringSetType(MatFDColoring matfd,MatMFFDType type) 38761ab5df0SBarry Smith { 38861ab5df0SBarry Smith PetscFunctionBegin; 38961ab5df0SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 39061ab5df0SBarry Smith /* 39161ab5df0SBarry Smith It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype 39261ab5df0SBarry Smith and this function is being provided as patch for a release so it shouldn't change the implementaton 39361ab5df0SBarry Smith */ 39461ab5df0SBarry Smith if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp"; 39561ab5df0SBarry Smith else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds"; 39698921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",type); 39761ab5df0SBarry Smith PetscFunctionReturn(0); 39861ab5df0SBarry Smith } 39961ab5df0SBarry Smith 400146574abSBarry Smith PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd,const char prefix[],const char optionname[]) 401005c665bSBarry Smith { 402e350db43SBarry Smith PetscBool flg; 4033050cee2SBarry Smith PetscViewer viewer; 404cffb1e40SBarry Smith PetscViewerFormat format; 405005c665bSBarry Smith 4063a40ed3dSBarry Smith PetscFunctionBegin; 407146574abSBarry Smith if (prefix) { 4089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,prefix,optionname,&viewer,&format,&flg)); 409146574abSBarry Smith } else { 4109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd),((PetscObject)fd)->options,((PetscObject)fd)->prefix,optionname,&viewer,&format,&flg)); 411146574abSBarry Smith } 412005c665bSBarry Smith if (flg) { 4139566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,format)); 4149566063dSJacob Faibussowitsch PetscCall(MatFDColoringView(fd,viewer)); 4159566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 4169566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 417005c665bSBarry Smith } 4183a40ed3dSBarry Smith PetscFunctionReturn(0); 419bbf0e169SBarry Smith } 420bbf0e169SBarry Smith 42105869f15SSatish Balay /*@ 422639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 423639f9d9dSBarry Smith computation of Jacobians. 424bbf0e169SBarry Smith 425ef5ee4d1SLois Curfman McInnes Collective on Mat 426ef5ee4d1SLois Curfman McInnes 427639f9d9dSBarry Smith Input Parameters: 428ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 4297086a01eSPeter Brune - iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring() 430bbf0e169SBarry Smith 431bbf0e169SBarry Smith Output Parameter: 432639f9d9dSBarry Smith . color - the new coloring context 433bbf0e169SBarry Smith 43415091d37SBarry Smith Level: intermediate 43515091d37SBarry Smith 436c2e3fba1SPatrick Sanan .seealso: `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`, 437db781477SPatrick Sanan `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`, 438db781477SPatrick Sanan `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()` 439bbf0e169SBarry Smith @*/ 4407087cfbeSBarry Smith PetscErrorCode MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 441bbf0e169SBarry Smith { 442639f9d9dSBarry Smith MatFDColoring c; 443639f9d9dSBarry Smith MPI_Comm comm; 44438baddfdSBarry Smith PetscInt M,N; 445639f9d9dSBarry Smith 4463a40ed3dSBarry Smith PetscFunctionBegin; 447c8a9c622SHong Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 44828b400f6SJacob Faibussowitsch PetscCheck(mat->assembled,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();"); 4499566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0)); 4509566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat,&M,&N)); 45108401ef6SPierre Jolivet PetscCheck(M == N,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Only for square matrices"); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mat,&comm)); 4539566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(c,MAT_FDCOLORING_CLASSID,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView)); 454639f9d9dSBarry Smith 455b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)mat,&c->matid)); 457b8f8c88eSHong Zhang 458*dbbe0bcdSBarry Smith PetscUseTypeMethod(mat,fdcoloringcreate ,iscoloring,c); 45993dfae19SHong Zhang 4609566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,NULL,&c->w1)); 461ce911d59SRichard Tran Mills /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 4629566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->w1,PETSC_TRUE)); 4639566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)c,(PetscObject)c->w1)); 4649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(c->w1,&c->w2)); 465ce911d59SRichard Tran Mills /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 4669566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->w2,PETSC_TRUE)); 4679566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)c,(PetscObject)c->w2)); 468b8f8c88eSHong Zhang 46977d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 47077d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 47149b058dcSBarry Smith c->currentcolor = -1; 472efb30889SBarry Smith c->htype = "wp"; 4734e269d77SPeter Brune c->fset = PETSC_FALSE; 474c8a9c622SHong Zhang c->setupcalled = PETSC_FALSE; 475639f9d9dSBarry Smith 476639f9d9dSBarry Smith *color = c; 4779566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)mat,"SNESMatFDColoring",(PetscObject)c)); 4789566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0)); 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 480bbf0e169SBarry Smith } 481bbf0e169SBarry Smith 48205869f15SSatish Balay /*@ 483639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 484639f9d9dSBarry Smith via MatFDColoringCreate(). 485bbf0e169SBarry Smith 486ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 487ef5ee4d1SLois Curfman McInnes 488b4fc646aSLois Curfman McInnes Input Parameter: 489639f9d9dSBarry Smith . c - coloring context 490bbf0e169SBarry Smith 49115091d37SBarry Smith Level: intermediate 49215091d37SBarry Smith 493db781477SPatrick Sanan .seealso: `MatFDColoringCreate()` 494bbf0e169SBarry Smith @*/ 4956bf464f9SBarry Smith PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) 496bbf0e169SBarry Smith { 49738baddfdSBarry Smith PetscInt i; 4980df34763SHong Zhang MatFDColoring color = *c; 499bbf0e169SBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 5016bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 502f4259b30SLisandro Dalcin if (--((PetscObject)color)->refct > 0) {*c = NULL; PetscFunctionReturn(0);} 503d4bb536fSBarry Smith 504071fcb05SBarry Smith /* we do not free the column arrays since their entries are owned by the ISs in color->isa */ 5050df34763SHong Zhang for (i=0; i<color->ncolors; i++) { 5069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&color->isa[i])); 507bbf0e169SBarry Smith } 5089566063dSJacob Faibussowitsch PetscCall(PetscFree(color->isa)); 5099566063dSJacob Faibussowitsch PetscCall(PetscFree2(color->ncolumns,color->columns)); 5109566063dSJacob Faibussowitsch PetscCall(PetscFree(color->nrows)); 5110df34763SHong Zhang if (color->htype[0] == 'w') { 5129566063dSJacob Faibussowitsch PetscCall(PetscFree(color->matentry2)); 5130df34763SHong Zhang } else { 5149566063dSJacob Faibussowitsch PetscCall(PetscFree(color->matentry)); 5150df34763SHong Zhang } 5169566063dSJacob Faibussowitsch PetscCall(PetscFree(color->dy)); 5179566063dSJacob Faibussowitsch if (color->vscale) PetscCall(VecDestroy(&color->vscale)); 5189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w1)); 5199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w2)); 5209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w3)); 5219566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(c)); 5223a40ed3dSBarry Smith PetscFunctionReturn(0); 523bbf0e169SBarry Smith } 52443a90d84SBarry Smith 52549b058dcSBarry Smith /*@C 52649b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 52749b058dcSBarry Smith that are currently being perturbed. 52849b058dcSBarry Smith 52949b058dcSBarry Smith Not Collective 53049b058dcSBarry Smith 531f899ff85SJose E. Roman Input Parameter: 53249b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 53349b058dcSBarry Smith 53449b058dcSBarry Smith Output Parameters: 53549b058dcSBarry Smith + n - the number of local columns being perturbed 53649b058dcSBarry Smith - cols - the column indices, in global numbering 53749b058dcSBarry Smith 53821fcc2ddSBarry Smith Level: advanced 53949b058dcSBarry Smith 540e66d0d93SMatthew Knepley Note: IF the matrix type is BAIJ, then the block column indices are returned 541e66d0d93SMatthew Knepley 542edaa7c33SBarry Smith Fortran Note: 543edaa7c33SBarry Smith This routine has a different interface for Fortran 54421fcc2ddSBarry Smith $ #include <petsc/finclude/petscmat.h> 54521fcc2ddSBarry Smith $ use petscmat 546edaa7c33SBarry Smith $ PetscInt, pointer :: array(:) 547edaa7c33SBarry Smith $ PetscErrorCode ierr 548edaa7c33SBarry Smith $ MatFDColoring i 549edaa7c33SBarry Smith $ call MatFDColoringGetPerturbedColumnsF90(i,array,ierr) 550edaa7c33SBarry Smith $ use the entries of array ... 551edaa7c33SBarry Smith $ call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr) 552edaa7c33SBarry Smith 553db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()` 55449b058dcSBarry Smith 55549b058dcSBarry Smith @*/ 556edaa7c33SBarry Smith PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,const PetscInt *cols[]) 55749b058dcSBarry Smith { 55849b058dcSBarry Smith PetscFunctionBegin; 55949b058dcSBarry Smith if (coloring->currentcolor >= 0) { 56049b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 56149b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 56249b058dcSBarry Smith } else { 56349b058dcSBarry Smith *n = 0; 56449b058dcSBarry Smith } 56549b058dcSBarry Smith PetscFunctionReturn(0); 56649b058dcSBarry Smith } 56749b058dcSBarry Smith 56843a90d84SBarry Smith /*@ 569e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 570e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 57143a90d84SBarry Smith 572fee21e36SBarry Smith Collective on MatFDColoring 573fee21e36SBarry Smith 574ef5ee4d1SLois Curfman McInnes Input Parameters: 575ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 576ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 577ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 5787850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 579ef5ee4d1SLois Curfman McInnes 5808bba8e72SBarry Smith Options Database Keys: 581ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 582b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 583e350db43SBarry Smith . -mat_fd_coloring_view draw - Activates drawing of coloring 584e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 5858bba8e72SBarry Smith 58615091d37SBarry Smith Level: intermediate 58798d79826SSatish Balay 588db781477SPatrick Sanan .seealso: `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()` 58943a90d84SBarry Smith 59043a90d84SBarry Smith @*/ 591d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 59243a90d84SBarry Smith { 593bdaf1daeSBarry Smith PetscBool eq; 5943acb8795SBarry Smith 5953acb8795SBarry Smith PetscFunctionBegin; 5960700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 5970700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 5980700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 5999566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)J,coloring->matid,&eq)); 60028b400f6SJacob Faibussowitsch PetscCheck(eq,PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()"); 60128b400f6SJacob Faibussowitsch PetscCheck(coloring->f,PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 60228b400f6SJacob Faibussowitsch PetscCheck(coloring->setupcalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetUp()"); 603684f2004SHong Zhang 6049566063dSJacob Faibussowitsch PetscCall(MatSetUnfactored(J)); 6059566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0)); 606*dbbe0bcdSBarry Smith PetscUseTypeMethod(J,fdcoloringapply ,coloring,x1,sctx); 6079566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0)); 6085bdb020cSBarry Smith if (!coloring->viewed) { 6099566063dSJacob Faibussowitsch PetscCall(MatFDColoringViewFromOptions(coloring,NULL,"-mat_fd_coloring_view")); 6105bdb020cSBarry Smith coloring->viewed = PETSC_TRUE; 6115bdb020cSBarry Smith } 6123acb8795SBarry Smith PetscFunctionReturn(0); 6133acb8795SBarry Smith } 614