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 109371c9d4SSatish Balay PetscErrorCode MatFDColoringSetF(MatFDColoring fd, Vec F) { 113a7fca6bSBarry Smith PetscFunctionBegin; 124e269d77SPeter Brune if (F) { 139566063dSJacob Faibussowitsch PetscCall(VecCopy(F, fd->w1)); 144e269d77SPeter Brune fd->fset = PETSC_TRUE; 154e269d77SPeter Brune } else { 164e269d77SPeter Brune fd->fset = PETSC_FALSE; 174e269d77SPeter Brune } 183a7fca6bSBarry Smith PetscFunctionReturn(0); 193a7fca6bSBarry Smith } 203a7fca6bSBarry Smith 219804daf3SBarry Smith #include <petscdraw.h> 229371c9d4SSatish Balay static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw, void *Aa) { 239194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 24b312708bSHong Zhang PetscInt i, j, nz, row; 259194eea9SBarry Smith PetscReal x, y; 26b312708bSHong Zhang MatEntry *Jentry = fd->matentry; 279194eea9SBarry Smith 289194eea9SBarry Smith PetscFunctionBegin; 299194eea9SBarry Smith /* loop over colors */ 30b312708bSHong Zhang nz = 0; 319194eea9SBarry Smith for (i = 0; i < fd->ncolors; i++) { 329194eea9SBarry Smith for (j = 0; j < fd->nrows[i]; j++) { 33b312708bSHong Zhang row = Jentry[nz].row; 34b312708bSHong Zhang y = fd->M - row - fd->rstart; 35b312708bSHong Zhang x = (PetscReal)Jentry[nz++].col; 369566063dSJacob Faibussowitsch PetscCall(PetscDrawRectangle(draw, x, y, x + 1, y + 1, i + 1, i + 1, i + 1, i + 1)); 379194eea9SBarry Smith } 389194eea9SBarry Smith } 399194eea9SBarry Smith PetscFunctionReturn(0); 409194eea9SBarry Smith } 419194eea9SBarry Smith 429371c9d4SSatish Balay static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd, PetscViewer viewer) { 43ace3abfcSBarry Smith PetscBool isnull; 44b0a32e0cSBarry Smith PetscDraw draw; 4554d96fa3SBarry Smith PetscReal xr, yr, xl, yl, h, w; 46005c665bSBarry Smith 473a40ed3dSBarry Smith PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 499566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 50832b7cebSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 51005c665bSBarry Smith 529371c9d4SSatish Balay xr = fd->N; 539371c9d4SSatish Balay yr = fd->M; 549371c9d4SSatish Balay h = yr / 10.0; 559371c9d4SSatish Balay w = xr / 10.0; 569371c9d4SSatish Balay xr += w; 579371c9d4SSatish Balay yr += h; 589371c9d4SSatish Balay xl = -w; 599371c9d4SSatish Balay yl = -h; 609566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr)); 619566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", (PetscObject)viewer)); 629566063dSJacob Faibussowitsch PetscCall(PetscDrawZoom(draw, MatFDColoringView_Draw_Zoom, fd)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)fd, "Zoomviewer", NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 653a40ed3dSBarry Smith PetscFunctionReturn(0); 66005c665bSBarry Smith } 67005c665bSBarry Smith 68bbf0e169SBarry Smith /*@C 69639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 70bbf0e169SBarry Smith 71*11a5261eSBarry Smith Collective on c 72fee21e36SBarry Smith 73ef5ee4d1SLois Curfman McInnes Input Parameters: 74ef5ee4d1SLois Curfman McInnes + c - the coloring context 75ef5ee4d1SLois Curfman McInnes - viewer - visualization context 76ef5ee4d1SLois Curfman McInnes 7715091d37SBarry Smith Level: intermediate 7815091d37SBarry Smith 79b4fc646aSLois Curfman McInnes Notes: 80b4fc646aSLois Curfman McInnes The available visualization contexts include 81*11a5261eSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 82*11a5261eSBarry Smith . `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 83ef5ee4d1SLois Curfman McInnes output where only the first processor opens 84ef5ee4d1SLois Curfman McInnes the file. All other processors send their 85ef5ee4d1SLois Curfman McInnes data to the first processor to print. 86*11a5261eSBarry Smith - `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure 87bbf0e169SBarry Smith 889194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 89b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 909194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 919194eea9SBarry Smith 92*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()` 93bbf0e169SBarry Smith @*/ 949371c9d4SSatish Balay PetscErrorCode MatFDColoringView(MatFDColoring c, PetscViewer viewer) { 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); 10148a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer)); 1020700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 103c9780b6fSBarry Smith PetscCheckSameComm(c, 1, viewer, 2); 104bbf0e169SBarry Smith 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1070f5bd95cSBarry Smith if (isdraw) { 1089566063dSJacob Faibussowitsch PetscCall(MatFDColoringView_Draw(c, viewer)); 10932077d6dSBarry Smith } else if (iascii) { 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)c, viewer)); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Error tolerance=%g\n", (double)c->error_rel)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Umin=%g\n", (double)c->umin)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of colors=%" PetscInt_FMT "\n", c->ncolors)); 114ae09f205SBarry Smith 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 116fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 117b312708bSHong Zhang PetscInt row, col, nz; 118b312708bSHong Zhang nz = 0; 119b4fc646aSLois Curfman McInnes for (i = 0; i < c->ncolors; i++) { 1209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Information for color %" PetscInt_FMT "\n", i)); 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of columns %" PetscInt_FMT "\n", c->ncolumns[i])); 12248a46eb9SPierre Jolivet for (j = 0; j < c->ncolumns[i]; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", c->columns[i][j])); 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of rows %" PetscInt_FMT "\n", c->nrows[i])); 1245bdb020cSBarry Smith if (c->matentry) { 125b4fc646aSLois Curfman McInnes for (j = 0; j < c->nrows[i]; j++) { 126b312708bSHong Zhang row = c->matentry[nz].row; 127b312708bSHong Zhang col = c->matentry[nz++].col; 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " \n", row, col)); 129b4fc646aSLois Curfman McInnes } 130bbf0e169SBarry Smith } 131bbf0e169SBarry Smith } 1325bdb020cSBarry Smith } 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 134bbf0e169SBarry Smith } 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 136639f9d9dSBarry Smith } 137639f9d9dSBarry Smith 138639f9d9dSBarry Smith /*@ 139b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 140b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 141639f9d9dSBarry Smith 142*11a5261eSBarry Smith Logically Collective on matfd 143ef5ee4d1SLois Curfman McInnes 144ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 145ef5ee4d1SLois Curfman McInnes .vb 14665f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 147d461c19dSHong Zhang htype = 'ds': 148f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 149f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 150ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 151d461c19dSHong Zhang 152d461c19dSHong Zhang htype = 'wp': 153d461c19dSHong Zhang h = error_rel * sqrt(1 + ||u||) 154ef5ee4d1SLois Curfman McInnes .ve 155639f9d9dSBarry Smith 156639f9d9dSBarry Smith Input Parameters: 157ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 158639f9d9dSBarry Smith . error_rel - relative error 159f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 160fee21e36SBarry Smith 16115091d37SBarry Smith Level: advanced 16215091d37SBarry Smith 163*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()` 164639f9d9dSBarry Smith @*/ 1659371c9d4SSatish Balay PetscErrorCode MatFDColoringSetParameters(MatFDColoring matfd, PetscReal error, PetscReal umin) { 1663a40ed3dSBarry Smith PetscFunctionBegin; 1670700a824SBarry Smith PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 168c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd, error, 2); 169c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(matfd, umin, 3); 170639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 171639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1723a40ed3dSBarry Smith PetscFunctionReturn(0); 173639f9d9dSBarry Smith } 174639f9d9dSBarry Smith 175f86b9fbaSHong Zhang /*@ 176c8a9c622SHong Zhang MatFDColoringSetBlockSize - Sets block size for efficient inserting entries of Jacobian matrix. 177005c665bSBarry Smith 178*11a5261eSBarry Smith Logically Collective on matfd 179f86b9fbaSHong Zhang 180f86b9fbaSHong Zhang Input Parameters: 181f86b9fbaSHong Zhang + coloring - the coloring context 182f86b9fbaSHong Zhang . brows - number of rows in the block 183f86b9fbaSHong Zhang - bcols - number of columns in the block 184f86b9fbaSHong Zhang 185f86b9fbaSHong Zhang Level: intermediate 186f86b9fbaSHong Zhang 187*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFromOptions()` 188f86b9fbaSHong Zhang @*/ 1899371c9d4SSatish Balay PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring matfd, PetscInt brows, PetscInt bcols) { 190f86b9fbaSHong Zhang PetscFunctionBegin; 191f86b9fbaSHong Zhang PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 192f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd, brows, 2); 193f86b9fbaSHong Zhang PetscValidLogicalCollectiveInt(matfd, bcols, 3); 194f86b9fbaSHong Zhang if (brows != PETSC_DEFAULT) matfd->brows = brows; 195f86b9fbaSHong Zhang if (bcols != PETSC_DEFAULT) matfd->bcols = bcols; 196f86b9fbaSHong Zhang PetscFunctionReturn(0); 197f86b9fbaSHong Zhang } 198f86b9fbaSHong Zhang 199f86b9fbaSHong Zhang /*@ 200f86b9fbaSHong Zhang MatFDColoringSetUp - Sets up the internal data structures of matrix coloring context for the later use. 201f86b9fbaSHong Zhang 202*11a5261eSBarry Smith Collective on mat 203f86b9fbaSHong Zhang 204f86b9fbaSHong Zhang Input Parameters: 205f86b9fbaSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 206*11a5261eSBarry Smith . iscoloring - the coloring of the matrix; usually obtained with `MatGetColoring()` or `DMCreateColoring()` 207f86b9fbaSHong Zhang - color - the matrix coloring context 208f86b9fbaSHong Zhang 209f86b9fbaSHong Zhang Level: beginner 210f86b9fbaSHong Zhang 211*11a5261eSBarry Smith Notes: 212*11a5261eSBarry Smith When the coloring type is `IS_COLORING_LOCAL` the coloring is in the local ordering of the unknowns. 213bdaf1daeSBarry Smith 214*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()` 215f86b9fbaSHong Zhang @*/ 2169371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUp(Mat mat, ISColoring iscoloring, MatFDColoring color) { 217bdaf1daeSBarry Smith PetscBool eq; 218f86b9fbaSHong Zhang 219f86b9fbaSHong Zhang PetscFunctionBegin; 220c8a9c622SHong Zhang PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 221c8a9c622SHong Zhang PetscValidHeaderSpecific(color, MAT_FDCOLORING_CLASSID, 3); 222c8a9c622SHong Zhang if (color->setupcalled) PetscFunctionReturn(0); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)mat, color->matid, &eq)); 22428b400f6SJacob Faibussowitsch PetscCheck(eq, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringSetUp() must be that used with MatFDColoringCreate()"); 225c8a9c622SHong Zhang 2269566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringSetUp, mat, 0, 0, 0)); 227dbbe0bcdSBarry Smith PetscUseTypeMethod(mat, fdcoloringsetup, iscoloring, color); 228c8a9c622SHong Zhang 229c8a9c622SHong Zhang color->setupcalled = PETSC_TRUE; 2309566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringSetUp, mat, 0, 0, 0)); 231f86b9fbaSHong Zhang PetscFunctionReturn(0); 232f86b9fbaSHong Zhang } 233ff0cfa39SBarry Smith 2344a9d489dSBarry Smith /*@C 2354a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 2364a9d489dSBarry Smith 2373f9fe445SBarry Smith Not Collective 2384a9d489dSBarry Smith 239f899ff85SJose E. Roman Input Parameter: 2404a9d489dSBarry Smith . coloring - the coloring context 2414a9d489dSBarry Smith 2424a9d489dSBarry Smith Output Parameters: 2434a9d489dSBarry Smith + f - the function 2444a9d489dSBarry Smith - fctx - the optional user-defined function context 2454a9d489dSBarry Smith 2464a9d489dSBarry Smith Level: intermediate 2474a9d489dSBarry Smith 248*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()` 2494a9d489dSBarry Smith @*/ 2509371c9d4SSatish Balay PetscErrorCode MatFDColoringGetFunction(MatFDColoring matfd, PetscErrorCode (**f)(void), void **fctx) { 2514a9d489dSBarry Smith PetscFunctionBegin; 2520700a824SBarry Smith PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 2534a9d489dSBarry Smith if (f) *f = matfd->f; 2544a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2554a9d489dSBarry Smith PetscFunctionReturn(0); 2564a9d489dSBarry Smith } 2574a9d489dSBarry Smith 258d64ed03dSBarry Smith /*@C 259005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 260005c665bSBarry Smith 261*11a5261eSBarry Smith Logically Collective on matfd 262fee21e36SBarry Smith 263ef5ee4d1SLois Curfman McInnes Input Parameters: 264ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 265ef5ee4d1SLois Curfman McInnes . f - the function 266ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 267ef5ee4d1SLois Curfman McInnes 2687850c7c0SBarry Smith Calling sequence of (*f) function: 269*11a5261eSBarry Smith For `SNES`: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 270*11a5261eSBarry Smith If not using `SNES`: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 27115091d37SBarry Smith 2727850c7c0SBarry Smith Level: advanced 2737850c7c0SBarry Smith 274*11a5261eSBarry Smith Note: 275*11a5261eSBarry Smith This function is usually used automatically by `SNES` (when one uses `SNESSetJacobian()` with the argument 276*11a5261eSBarry Smith `SNESComputeJacobianDefaultColor()`) and only needs to be used by someone computing a matrix via coloring directly by 277*11a5261eSBarry Smith calling `MatFDColoringApply()` 2787850c7c0SBarry Smith 279*11a5261eSBarry Smith Fortran Note: 280*11a5261eSBarry Smith In Fortran you must call `MatFDColoringSetFunction()` for a coloring object to 281*11a5261eSBarry Smith be used without `SNES` or within the `SNES` solvers. 282f881d145SBarry Smith 283*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringGetFunction()`, `MatFDColoringSetFromOptions()` 284005c665bSBarry Smith @*/ 2859371c9d4SSatish Balay PetscErrorCode MatFDColoringSetFunction(MatFDColoring matfd, PetscErrorCode (*f)(void), void *fctx) { 2863a40ed3dSBarry Smith PetscFunctionBegin; 2870700a824SBarry Smith PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 288005c665bSBarry Smith matfd->f = f; 289005c665bSBarry Smith matfd->fctx = fctx; 2903a40ed3dSBarry Smith PetscFunctionReturn(0); 291005c665bSBarry Smith } 292005c665bSBarry Smith 293639f9d9dSBarry Smith /*@ 294b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 295639f9d9dSBarry Smith the options database. 296639f9d9dSBarry Smith 297*11a5261eSBarry Smith Collective on matfd 298fee21e36SBarry Smith 29965f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 300ef5ee4d1SLois Curfman McInnes .vb 30165f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 302f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 303f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 304ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 305ef5ee4d1SLois Curfman McInnes .ve 306ef5ee4d1SLois Curfman McInnes 307ef5ee4d1SLois Curfman McInnes Input Parameter: 308ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 309ef5ee4d1SLois Curfman McInnes 310b4fc646aSLois Curfman McInnes Options Database Keys: 3115620d6dcSBarry Smith + -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 312f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 3133ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 314ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 315e350db43SBarry Smith . -mat_fd_coloring_view ::ascii_info - Activates viewing info 316e350db43SBarry Smith - -mat_fd_coloring_view draw - Activates drawing 317639f9d9dSBarry Smith 31815091d37SBarry Smith Level: intermediate 31915091d37SBarry Smith 320*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 321639f9d9dSBarry Smith @*/ 3229371c9d4SSatish Balay PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring matfd) { 323ace3abfcSBarry Smith PetscBool flg; 324efb30889SBarry Smith char value[3]; 3253a40ed3dSBarry Smith 3263a40ed3dSBarry Smith PetscFunctionBegin; 3270700a824SBarry Smith PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 328639f9d9dSBarry Smith 329d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)matfd); 3309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_fd_coloring_err", "Square root of relative error in function", "MatFDColoringSetParameters", matfd->error_rel, &matfd->error_rel, NULL)); 3319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_fd_coloring_umin", "Minimum allowable u magnitude", "MatFDColoringSetParameters", matfd->umin, &matfd->umin, NULL)); 3329566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-mat_fd_type", "Algorithm to compute h, wp or ds", "MatFDColoringCreate", matfd->htype, value, sizeof(value), &flg)); 333efb30889SBarry Smith if (flg) { 334efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 335efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 33698921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", value); 337efb30889SBarry Smith } 3389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_fd_coloring_brows", "Number of block rows", "MatFDColoringSetBlockSize", matfd->brows, &matfd->brows, NULL)); 3399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_fd_coloring_bcols", "Number of block columns", "MatFDColoringSetBlockSize", matfd->bcols, &matfd->bcols, &flg)); 34093dfae19SHong Zhang if (flg && matfd->bcols > matfd->ncolors) { 34193dfae19SHong Zhang /* input bcols cannot be > matfd->ncolors, thus set it as ncolors */ 34293dfae19SHong Zhang matfd->bcols = matfd->ncolors; 34393dfae19SHong Zhang } 344f86b9fbaSHong Zhang 3455d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 346dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)matfd, PetscOptionsObject)); 347d0609cedSBarry Smith PetscOptionsEnd(); 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 349005c665bSBarry Smith } 350005c665bSBarry Smith 35161ab5df0SBarry Smith /*@C 35261ab5df0SBarry Smith MatFDColoringSetType - Sets the approach for computing the finite difference parameter 35361ab5df0SBarry Smith 354*11a5261eSBarry Smith Collective on matfd 35561ab5df0SBarry Smith 35661ab5df0SBarry Smith Input Parameters: 35761ab5df0SBarry Smith + coloring - the coloring context 358*11a5261eSBarry Smith - type - either `MATMFFD_WP` or `MATMFFD_DS` 35961ab5df0SBarry Smith 36061ab5df0SBarry Smith Options Database Keys: 36161ab5df0SBarry Smith . -mat_fd_type - "wp" or "ds" 36261ab5df0SBarry Smith 363*11a5261eSBarry Smith Note: 364*11a5261eSBarry Smith It is goofy that the argument type is `MatMFFDType` since the `MatFDColoring` actually computes the matrix entries 365*11a5261eSBarry 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 36661ab5df0SBarry Smith introducing another one. 36761ab5df0SBarry Smith 36861ab5df0SBarry Smith Level: intermediate 36961ab5df0SBarry Smith 370*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringView()`, `MatFDColoringSetParameters()` 37161ab5df0SBarry Smith @*/ 3729371c9d4SSatish Balay PetscErrorCode MatFDColoringSetType(MatFDColoring matfd, MatMFFDType type) { 37361ab5df0SBarry Smith PetscFunctionBegin; 37461ab5df0SBarry Smith PetscValidHeaderSpecific(matfd, MAT_FDCOLORING_CLASSID, 1); 37561ab5df0SBarry Smith /* 37661ab5df0SBarry Smith It is goofy to handle the strings this way but currently there is no code to free a dynamically created matfd->htype 37761ab5df0SBarry Smith and this function is being provided as patch for a release so it shouldn't change the implementaton 37861ab5df0SBarry Smith */ 37961ab5df0SBarry Smith if (type[0] == 'w' && type[1] == 'p') matfd->htype = "wp"; 38061ab5df0SBarry Smith else if (type[0] == 'd' && type[1] == 's') matfd->htype = "ds"; 38198921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown finite differencing type %s", type); 38261ab5df0SBarry Smith PetscFunctionReturn(0); 38361ab5df0SBarry Smith } 38461ab5df0SBarry Smith 3859371c9d4SSatish Balay PetscErrorCode MatFDColoringViewFromOptions(MatFDColoring fd, const char prefix[], const char optionname[]) { 386e350db43SBarry Smith PetscBool flg; 3873050cee2SBarry Smith PetscViewer viewer; 388cffb1e40SBarry Smith PetscViewerFormat format; 389005c665bSBarry Smith 3903a40ed3dSBarry Smith PetscFunctionBegin; 391146574abSBarry Smith if (prefix) { 3929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, prefix, optionname, &viewer, &format, &flg)); 393146574abSBarry Smith } else { 3949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)fd), ((PetscObject)fd)->options, ((PetscObject)fd)->prefix, optionname, &viewer, &format, &flg)); 395146574abSBarry Smith } 396005c665bSBarry Smith if (flg) { 3979566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 3989566063dSJacob Faibussowitsch PetscCall(MatFDColoringView(fd, viewer)); 3999566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 4009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 401005c665bSBarry Smith } 4023a40ed3dSBarry Smith PetscFunctionReturn(0); 403bbf0e169SBarry Smith } 404bbf0e169SBarry Smith 40505869f15SSatish Balay /*@ 406639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 407639f9d9dSBarry Smith computation of Jacobians. 408bbf0e169SBarry Smith 409*11a5261eSBarry Smith Collective on mat 410ef5ee4d1SLois Curfman McInnes 411639f9d9dSBarry Smith Input Parameters: 412ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 413*11a5261eSBarry Smith - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()` 414bbf0e169SBarry Smith 415bbf0e169SBarry Smith Output Parameter: 416639f9d9dSBarry Smith . color - the new coloring context 417bbf0e169SBarry Smith 41815091d37SBarry Smith Level: intermediate 41915091d37SBarry Smith 420*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringDestroy()`, `SNESComputeJacobianDefaultColor()`, `ISColoringCreate()`, 421db781477SPatrick Sanan `MatFDColoringSetFunction()`, `MatFDColoringSetFromOptions()`, `MatFDColoringApply()`, 422db781477SPatrick Sanan `MatFDColoringView()`, `MatFDColoringSetParameters()`, `MatColoringCreate()`, `DMCreateColoring()`, `MatFDColoringSetValues()` 423bbf0e169SBarry Smith @*/ 4249371c9d4SSatish Balay PetscErrorCode MatFDColoringCreate(Mat mat, ISColoring iscoloring, MatFDColoring *color) { 425639f9d9dSBarry Smith MatFDColoring c; 426639f9d9dSBarry Smith MPI_Comm comm; 42738baddfdSBarry Smith PetscInt M, N; 428639f9d9dSBarry Smith 4293a40ed3dSBarry Smith PetscFunctionBegin; 430c8a9c622SHong Zhang PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 43128b400f6SJacob Faibussowitsch PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be assembled by calls to MatAssemblyBegin/End();"); 4329566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringCreate, mat, 0, 0, 0)); 4339566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &M, &N)); 43408401ef6SPierre Jolivet PetscCheck(M == N, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Only for square matrices"); 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 4369566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(c, MAT_FDCOLORING_CLASSID, "MatFDColoring", "Jacobian computation via finite differences with coloring", "Mat", comm, MatFDColoringDestroy, MatFDColoringView)); 437639f9d9dSBarry Smith 438b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)mat, &c->matid)); 440b8f8c88eSHong Zhang 441dbbe0bcdSBarry Smith PetscUseTypeMethod(mat, fdcoloringcreate, iscoloring, c); 44293dfae19SHong Zhang 4439566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, NULL, &c->w1)); 444ce911d59SRichard Tran Mills /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 4459566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->w1, PETSC_TRUE)); 4469566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)c, (PetscObject)c->w1)); 4479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(c->w1, &c->w2)); 448ce911d59SRichard Tran Mills /* Vec is used intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */ 4499566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(c->w2, PETSC_TRUE)); 4509566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)c, (PetscObject)c->w2)); 451b8f8c88eSHong Zhang 45277d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 45377d8c4bbSBarry Smith c->umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 45449b058dcSBarry Smith c->currentcolor = -1; 455efb30889SBarry Smith c->htype = "wp"; 4564e269d77SPeter Brune c->fset = PETSC_FALSE; 457c8a9c622SHong Zhang c->setupcalled = PETSC_FALSE; 458639f9d9dSBarry Smith 459639f9d9dSBarry Smith *color = c; 4609566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)mat, "SNESMatFDColoring", (PetscObject)c)); 4619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringCreate, mat, 0, 0, 0)); 4623a40ed3dSBarry Smith PetscFunctionReturn(0); 463bbf0e169SBarry Smith } 464bbf0e169SBarry Smith 46505869f15SSatish Balay /*@ 466639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 467*11a5261eSBarry Smith via `MatFDColoringCreate()`. 468bbf0e169SBarry Smith 469*11a5261eSBarry Smith Collective on c 470ef5ee4d1SLois Curfman McInnes 471b4fc646aSLois Curfman McInnes Input Parameter: 472639f9d9dSBarry Smith . c - coloring context 473bbf0e169SBarry Smith 47415091d37SBarry Smith Level: intermediate 47515091d37SBarry Smith 476*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()` 477bbf0e169SBarry Smith @*/ 4789371c9d4SSatish Balay PetscErrorCode MatFDColoringDestroy(MatFDColoring *c) { 47938baddfdSBarry Smith PetscInt i; 4800df34763SHong Zhang MatFDColoring color = *c; 481bbf0e169SBarry Smith 4823a40ed3dSBarry Smith PetscFunctionBegin; 4836bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 4849371c9d4SSatish Balay if (--((PetscObject)color)->refct > 0) { 4859371c9d4SSatish Balay *c = NULL; 4869371c9d4SSatish Balay PetscFunctionReturn(0); 4879371c9d4SSatish Balay } 488d4bb536fSBarry Smith 489071fcb05SBarry Smith /* we do not free the column arrays since their entries are owned by the ISs in color->isa */ 49048a46eb9SPierre Jolivet for (i = 0; i < color->ncolors; i++) PetscCall(ISDestroy(&color->isa[i])); 4919566063dSJacob Faibussowitsch PetscCall(PetscFree(color->isa)); 4929566063dSJacob Faibussowitsch PetscCall(PetscFree2(color->ncolumns, color->columns)); 4939566063dSJacob Faibussowitsch PetscCall(PetscFree(color->nrows)); 4940df34763SHong Zhang if (color->htype[0] == 'w') { 4959566063dSJacob Faibussowitsch PetscCall(PetscFree(color->matentry2)); 4960df34763SHong Zhang } else { 4979566063dSJacob Faibussowitsch PetscCall(PetscFree(color->matentry)); 4980df34763SHong Zhang } 4999566063dSJacob Faibussowitsch PetscCall(PetscFree(color->dy)); 5009566063dSJacob Faibussowitsch if (color->vscale) PetscCall(VecDestroy(&color->vscale)); 5019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w1)); 5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w2)); 5039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&color->w3)); 5049566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(c)); 5053a40ed3dSBarry Smith PetscFunctionReturn(0); 506bbf0e169SBarry Smith } 50743a90d84SBarry Smith 50849b058dcSBarry Smith /*@C 50949b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 51049b058dcSBarry Smith that are currently being perturbed. 51149b058dcSBarry Smith 51249b058dcSBarry Smith Not Collective 51349b058dcSBarry Smith 514f899ff85SJose E. Roman Input Parameter: 515*11a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()` 51649b058dcSBarry Smith 51749b058dcSBarry Smith Output Parameters: 51849b058dcSBarry Smith + n - the number of local columns being perturbed 51949b058dcSBarry Smith - cols - the column indices, in global numbering 52049b058dcSBarry Smith 52121fcc2ddSBarry Smith Level: advanced 52249b058dcSBarry Smith 523*11a5261eSBarry Smith Note: 524*11a5261eSBarry Smith IF the matrix type is `MATBAIJ`, then the block column indices are returned 525e66d0d93SMatthew Knepley 526edaa7c33SBarry Smith Fortran Note: 527edaa7c33SBarry Smith This routine has a different interface for Fortran 52821fcc2ddSBarry Smith $ #include <petsc/finclude/petscmat.h> 52921fcc2ddSBarry Smith $ use petscmat 530edaa7c33SBarry Smith $ PetscInt, pointer :: array(:) 531edaa7c33SBarry Smith $ PetscErrorCode ierr 532edaa7c33SBarry Smith $ MatFDColoring i 533edaa7c33SBarry Smith $ call MatFDColoringGetPerturbedColumnsF90(i,array,ierr) 534edaa7c33SBarry Smith $ use the entries of array ... 535edaa7c33SBarry Smith $ call MatFDColoringRestorePerturbedColumnsF90(i,array,ierr) 536edaa7c33SBarry Smith 537*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringApply()` 53849b058dcSBarry Smith @*/ 5399371c9d4SSatish Balay PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring coloring, PetscInt *n, const PetscInt *cols[]) { 54049b058dcSBarry Smith PetscFunctionBegin; 54149b058dcSBarry Smith if (coloring->currentcolor >= 0) { 54249b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 54349b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 54449b058dcSBarry Smith } else { 54549b058dcSBarry Smith *n = 0; 54649b058dcSBarry Smith } 54749b058dcSBarry Smith PetscFunctionReturn(0); 54849b058dcSBarry Smith } 54949b058dcSBarry Smith 55043a90d84SBarry Smith /*@ 551*11a5261eSBarry Smith MatFDColoringApply - Given a matrix for which a `MatFDColoring` context 552e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 55343a90d84SBarry Smith 554*11a5261eSBarry Smith Collective on J 555fee21e36SBarry Smith 556ef5ee4d1SLois Curfman McInnes Input Parameters: 557ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 558*11a5261eSBarry Smith . coloring - coloring context created with `MatFDColoringCreate()` 559ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 560*11a5261eSBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is `SNES` object, otherwise it is null 561ef5ee4d1SLois Curfman McInnes 5628bba8e72SBarry Smith Options Database Keys: 563*11a5261eSBarry Smith + -mat_fd_type - "wp" or "ds" (see `MATMFFD_WP` or `MATMFFD_DS`) 564b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 565e350db43SBarry Smith . -mat_fd_coloring_view draw - Activates drawing of coloring 566e350db43SBarry Smith - -mat_fd_coloring_view ::ascii_info - Activates viewing of coloring info 5678bba8e72SBarry Smith 56815091d37SBarry Smith Level: intermediate 56998d79826SSatish Balay 570*11a5261eSBarry Smith .seealso: `MatFDColoring`, `MatFDColoringCreate()`, `MatFDColoringDestroy()`, `MatFDColoringView()`, `MatFDColoringSetFunction()`, `MatFDColoringSetValues()` 57143a90d84SBarry Smith @*/ 5729371c9d4SSatish Balay PetscErrorCode MatFDColoringApply(Mat J, MatFDColoring coloring, Vec x1, void *sctx) { 573bdaf1daeSBarry Smith PetscBool eq; 5743acb8795SBarry Smith 5753acb8795SBarry Smith PetscFunctionBegin; 5760700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 5770700a824SBarry Smith PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 2); 5780700a824SBarry Smith PetscValidHeaderSpecific(x1, VEC_CLASSID, 3); 5799566063dSJacob Faibussowitsch PetscCall(PetscObjectCompareId((PetscObject)J, coloring->matid, &eq)); 58028b400f6SJacob Faibussowitsch PetscCheck(eq, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONG, "Matrix used with MatFDColoringApply() must be that used with MatFDColoringCreate()"); 58128b400f6SJacob Faibussowitsch PetscCheck(coloring->f, PetscObjectComm((PetscObject)J), PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetFunction()"); 58228b400f6SJacob Faibussowitsch PetscCheck(coloring->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatFDColoringSetUp()"); 583684f2004SHong Zhang 5849566063dSJacob Faibussowitsch PetscCall(MatSetUnfactored(J)); 5859566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_FDColoringApply, coloring, J, x1, 0)); 586dbbe0bcdSBarry Smith PetscUseTypeMethod(J, fdcoloringapply, coloring, x1, sctx); 5879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_FDColoringApply, coloring, J, x1, 0)); 5885bdb020cSBarry Smith if (!coloring->viewed) { 5899566063dSJacob Faibussowitsch PetscCall(MatFDColoringViewFromOptions(coloring, NULL, "-mat_fd_coloring_view")); 5905bdb020cSBarry Smith coloring->viewed = PETSC_TRUE; 5915bdb020cSBarry Smith } 5923acb8795SBarry Smith PetscFunctionReturn(0); 5933acb8795SBarry Smith } 594