1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3bbf0e169SBarry Smith /* 4639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 5639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 6bbf0e169SBarry Smith */ 7bbf0e169SBarry Smith 8b9147fbbSdalcinl #include "include/private/matimpl.h" /*I "petscmat.h" I*/ 9bbf0e169SBarry Smith 108ba1e511SMatthew Knepley /* Logging support */ 11166c7f25SBarry Smith PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 128ba1e511SMatthew Knepley 134a2ae208SSatish Balay #undef __FUNCT__ 143a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF" 15be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F) 163a7fca6bSBarry Smith { 173a7fca6bSBarry Smith PetscFunctionBegin; 183a7fca6bSBarry Smith fd->F = F; 193a7fca6bSBarry Smith PetscFunctionReturn(0); 203a7fca6bSBarry Smith } 213a7fca6bSBarry Smith 223a7fca6bSBarry Smith #undef __FUNCT__ 234a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 246849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 259194eea9SBarry Smith { 269194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 27dfbe8321SBarry Smith PetscErrorCode ierr; 2838baddfdSBarry Smith PetscInt i,j; 299194eea9SBarry Smith PetscReal x,y; 309194eea9SBarry Smith 319194eea9SBarry Smith PetscFunctionBegin; 329194eea9SBarry Smith 339194eea9SBarry Smith /* loop over colors */ 349194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 359194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 369194eea9SBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 379194eea9SBarry Smith x = fd->columnsforrow[i][j]; 38b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 399194eea9SBarry Smith } 409194eea9SBarry Smith } 419194eea9SBarry Smith PetscFunctionReturn(0); 429194eea9SBarry Smith } 439194eea9SBarry Smith 444a2ae208SSatish Balay #undef __FUNCT__ 454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 466849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 47005c665bSBarry Smith { 48dfbe8321SBarry Smith PetscErrorCode ierr; 49005c665bSBarry Smith PetscTruth isnull; 50b0a32e0cSBarry Smith PetscDraw draw; 5154d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 52005c665bSBarry Smith 533a40ed3dSBarry Smith PetscFunctionBegin; 54b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 55b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 569194eea9SBarry Smith 579194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 58005c665bSBarry Smith 59005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 60005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 61b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 62b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 639194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 65005c665bSBarry Smith } 66005c665bSBarry Smith 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 69bbf0e169SBarry Smith /*@C 70639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 71bbf0e169SBarry Smith 724c49b128SBarry Smith Collective on MatFDColoring 73fee21e36SBarry Smith 74ef5ee4d1SLois Curfman McInnes Input Parameters: 75ef5ee4d1SLois Curfman McInnes + c - the coloring context 76ef5ee4d1SLois Curfman McInnes - viewer - visualization context 77ef5ee4d1SLois Curfman McInnes 7815091d37SBarry Smith Level: intermediate 7915091d37SBarry Smith 80b4fc646aSLois Curfman McInnes Notes: 81b4fc646aSLois Curfman McInnes The available visualization contexts include 82b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 83b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 84ef5ee4d1SLois Curfman McInnes output where only the first processor opens 85ef5ee4d1SLois Curfman McInnes the file. All other processors send their 86ef5ee4d1SLois Curfman McInnes data to the first processor to print. 87b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 88bbf0e169SBarry Smith 899194eea9SBarry Smith Notes: 909194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 91b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 929194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 939194eea9SBarry Smith 94639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 95005c665bSBarry Smith 96b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 97bbf0e169SBarry Smith @*/ 98be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer) 99bbf0e169SBarry Smith { 1006849ba73SBarry Smith PetscErrorCode ierr; 10138baddfdSBarry Smith PetscInt i,j; 10232077d6dSBarry Smith PetscTruth isdraw,iascii; 103f3ef73ceSBarry Smith PetscViewerFormat format; 104bbf0e169SBarry Smith 1053a40ed3dSBarry Smith PetscFunctionBegin; 1064482741eSBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE,1); 1073050cee2SBarry Smith if (!viewer) { 1087adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr); 1093050cee2SBarry Smith } 1104482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 111c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 112bbf0e169SBarry Smith 113fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 11432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1150f5bd95cSBarry Smith if (isdraw) { 116b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 11732077d6dSBarry Smith } else if (iascii) { 118b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 119a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 120a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 12177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 122ae09f205SBarry Smith 123b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 124fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 125b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 12677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 12777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 128b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 12977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 130639f9d9dSBarry Smith } 13177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 132b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 13377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 134b4fc646aSLois Curfman McInnes } 135bbf0e169SBarry Smith } 136bbf0e169SBarry Smith } 137b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1385cd90555SBarry Smith } else { 13979a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 140bbf0e169SBarry Smith } 1413a40ed3dSBarry Smith PetscFunctionReturn(0); 142639f9d9dSBarry Smith } 143639f9d9dSBarry Smith 1444a2ae208SSatish Balay #undef __FUNCT__ 1454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 146639f9d9dSBarry Smith /*@ 147b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 148b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 149639f9d9dSBarry Smith 150ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 151ef5ee4d1SLois Curfman McInnes 152ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 153ef5ee4d1SLois Curfman McInnes .vb 15465f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 155f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 156f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 157ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 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 167b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 168b4fc646aSLois Curfman McInnes 169b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 170639f9d9dSBarry Smith @*/ 171be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 172639f9d9dSBarry Smith { 1733a40ed3dSBarry Smith PetscFunctionBegin; 1744482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 175639f9d9dSBarry Smith 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 181005c665bSBarry Smith 182ff0cfa39SBarry Smith 1834a2ae208SSatish Balay #undef __FUNCT__ 1844a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction" 1854a9d489dSBarry Smith /*@C 1864a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 1874a9d489dSBarry Smith 1884a9d489dSBarry Smith Collective on MatFDColoring 1894a9d489dSBarry Smith 1904a9d489dSBarry Smith Input Parameters: 1914a9d489dSBarry Smith . coloring - the coloring context 1924a9d489dSBarry Smith 1934a9d489dSBarry Smith Output Parameters: 1944a9d489dSBarry Smith + f - the function 1954a9d489dSBarry Smith - fctx - the optional user-defined function context 1964a9d489dSBarry Smith 1974a9d489dSBarry Smith Level: intermediate 1984a9d489dSBarry Smith 1994a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 2004a9d489dSBarry Smith @*/ 2014a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2024a9d489dSBarry Smith { 2034a9d489dSBarry Smith PetscFunctionBegin; 2044a9d489dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 2054a9d489dSBarry Smith if (f) *f = matfd->f; 2064a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2074a9d489dSBarry Smith PetscFunctionReturn(0); 2084a9d489dSBarry Smith } 2094a9d489dSBarry Smith 2104a9d489dSBarry Smith #undef __FUNCT__ 2114a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 212d64ed03dSBarry Smith /*@C 213005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 214005c665bSBarry Smith 215fee21e36SBarry Smith Collective on MatFDColoring 216fee21e36SBarry Smith 217ef5ee4d1SLois Curfman McInnes Input Parameters: 218ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 219ef5ee4d1SLois Curfman McInnes . f - the function 220ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 221ef5ee4d1SLois Curfman McInnes 2227850c7c0SBarry Smith Calling sequence of (*f) function: 2237850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 2247850c7c0SBarry Smith For TS: PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*) 2257850c7c0SBarry Smith If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 22615091d37SBarry Smith 2277850c7c0SBarry Smith Level: advanced 2287850c7c0SBarry Smith 2297850c7c0SBarry Smith Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument 2307850c7c0SBarry Smith SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used 2317850c7c0SBarry Smith by someone computing a matrix via coloring directly by calling MatFDColoringApply() 2327850c7c0SBarry Smith 2337850c7c0SBarry Smith Fortran Notes: 2347850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 2357850c7c0SBarry Smith be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 2367850c7c0SBarry Smith within the TS solvers. 237f881d145SBarry Smith 238b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 239005c665bSBarry Smith @*/ 240be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 241005c665bSBarry Smith { 2423a40ed3dSBarry Smith PetscFunctionBegin; 2434482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 244005c665bSBarry Smith matfd->f = f; 245005c665bSBarry Smith matfd->fctx = fctx; 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 247005c665bSBarry Smith } 248005c665bSBarry Smith 2494a2ae208SSatish Balay #undef __FUNCT__ 2504a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 251639f9d9dSBarry Smith /*@ 252b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 253639f9d9dSBarry Smith the options database. 254639f9d9dSBarry Smith 255fee21e36SBarry Smith Collective on MatFDColoring 256fee21e36SBarry Smith 25765f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 258ef5ee4d1SLois Curfman McInnes .vb 25965f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 260f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 261f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 262ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 263ef5ee4d1SLois Curfman McInnes .ve 264ef5ee4d1SLois Curfman McInnes 265ef5ee4d1SLois Curfman McInnes Input Parameter: 266ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 267ef5ee4d1SLois Curfman McInnes 268b4fc646aSLois Curfman McInnes Options Database Keys: 269d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 270ef5ee4d1SLois Curfman McInnes of relative error in the function) 271f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 2723ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 273ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 274ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 275ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 276639f9d9dSBarry Smith 27715091d37SBarry Smith Level: intermediate 27815091d37SBarry Smith 279b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 280d1c39f55SBarry Smith 281d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 282d1c39f55SBarry Smith 283639f9d9dSBarry Smith @*/ 284be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd) 285639f9d9dSBarry Smith { 286dfbe8321SBarry Smith PetscErrorCode ierr; 287efb30889SBarry Smith PetscTruth flg; 288efb30889SBarry Smith char value[3]; 2893a40ed3dSBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 2914482741eSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE,1); 292639f9d9dSBarry Smith 2937adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)matfd)->comm,((PetscObject)matfd)->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 29487828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 29587828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 296efb30889SBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,2,&flg);CHKERRQ(ierr); 297efb30889SBarry Smith if (flg) { 298efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 299efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 300efb30889SBarry Smith else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 301efb30889SBarry Smith } 302186905e3SBarry Smith /* not used here; but so they are presented in the GUI */ 303b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 304b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 305b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 306b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3073a40ed3dSBarry Smith PetscFunctionReturn(0); 308005c665bSBarry Smith } 309005c665bSBarry Smith 3104a2ae208SSatish Balay #undef __FUNCT__ 3114a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private" 312dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd) 313005c665bSBarry Smith { 314dfbe8321SBarry Smith PetscErrorCode ierr; 315f1af5d2fSBarry Smith PetscTruth flg; 3163050cee2SBarry Smith PetscViewer viewer; 317005c665bSBarry Smith 3183a40ed3dSBarry Smith PetscFunctionBegin; 3197adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr); 320b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 321005c665bSBarry Smith if (flg) { 3223050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 323005c665bSBarry Smith } 324b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 325ae09f205SBarry Smith if (flg) { 3263050cee2SBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 3273050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 3283050cee2SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 329ae09f205SBarry Smith } 330b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 331005c665bSBarry Smith if (flg) { 3327adad957SLisandro Dalcin ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 3337adad957SLisandro Dalcin ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 334005c665bSBarry Smith } 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 336bbf0e169SBarry Smith } 337bbf0e169SBarry Smith 3384a2ae208SSatish Balay #undef __FUNCT__ 3394a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 34005869f15SSatish Balay /*@ 341639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 342639f9d9dSBarry Smith computation of Jacobians. 343bbf0e169SBarry Smith 344ef5ee4d1SLois Curfman McInnes Collective on Mat 345ef5ee4d1SLois Curfman McInnes 346639f9d9dSBarry Smith Input Parameters: 347ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 348ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 349bbf0e169SBarry Smith 350bbf0e169SBarry Smith Output Parameter: 351639f9d9dSBarry Smith . color - the new coloring context 352bbf0e169SBarry Smith 35315091d37SBarry Smith Level: intermediate 35415091d37SBarry Smith 3556831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 356d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 357*ebd3b9afSBarry Smith MatFDColoringView(), MatFDColoringSetParameters() 358bbf0e169SBarry Smith @*/ 359be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 360bbf0e169SBarry Smith { 361639f9d9dSBarry Smith MatFDColoring c; 362639f9d9dSBarry Smith MPI_Comm comm; 363dfbe8321SBarry Smith PetscErrorCode ierr; 36438baddfdSBarry Smith PetscInt M,N; 365b8f8c88eSHong Zhang PetscMPIInt size; 366639f9d9dSBarry Smith 3673a40ed3dSBarry Smith PetscFunctionBegin; 368d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 369639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 37029bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 371639f9d9dSBarry Smith 372f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 37352e6d16bSBarry Smith ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 374639f9d9dSBarry Smith 375b8f8c88eSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 376b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 377b8f8c88eSHong Zhang 378f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 379f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 380639f9d9dSBarry Smith } else { 38129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 382639f9d9dSBarry Smith } 383639f9d9dSBarry Smith 384b8f8c88eSHong Zhang ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr); 385b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr); 386b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 387b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr); 388b8f8c88eSHong Zhang 38977d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 39077d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 39149b058dcSBarry Smith c->currentcolor = -1; 392efb30889SBarry Smith c->htype = "wp"; 393639f9d9dSBarry Smith 394639f9d9dSBarry Smith *color = c; 395d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 3963a40ed3dSBarry Smith PetscFunctionReturn(0); 397bbf0e169SBarry Smith } 398bbf0e169SBarry Smith 3994a2ae208SSatish Balay #undef __FUNCT__ 4004a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 40105869f15SSatish Balay /*@ 402639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 403639f9d9dSBarry Smith via MatFDColoringCreate(). 404bbf0e169SBarry Smith 405ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 406ef5ee4d1SLois Curfman McInnes 407b4fc646aSLois Curfman McInnes Input Parameter: 408639f9d9dSBarry Smith . c - coloring context 409bbf0e169SBarry Smith 41015091d37SBarry Smith Level: intermediate 41115091d37SBarry Smith 412639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 413bbf0e169SBarry Smith @*/ 414be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c) 415bbf0e169SBarry Smith { 4166849ba73SBarry Smith PetscErrorCode ierr; 41738baddfdSBarry Smith PetscInt i; 418bbf0e169SBarry Smith 4193a40ed3dSBarry Smith PetscFunctionBegin; 4207adad957SLisandro Dalcin if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0); 421d4bb536fSBarry Smith 422639f9d9dSBarry Smith for (i=0; i<c->ncolors; i++) { 42305b42c5fSBarry Smith ierr = PetscFree(c->columns[i]);CHKERRQ(ierr); 42405b42c5fSBarry Smith ierr = PetscFree(c->rows[i]);CHKERRQ(ierr); 42505b42c5fSBarry Smith ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr); 42605b42c5fSBarry Smith if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 427bbf0e169SBarry Smith } 428606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 429606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 430606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 431606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 432606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 43305b42c5fSBarry Smith ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr); 4344f113abfSBarry Smith if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 435005c665bSBarry Smith if (c->w1) { 436005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 437005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 438b8f8c88eSHong Zhang } 439b8f8c88eSHong Zhang if (c->w3){ 440005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 441005c665bSBarry Smith } 442d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 4433a40ed3dSBarry Smith PetscFunctionReturn(0); 444bbf0e169SBarry Smith } 44543a90d84SBarry Smith 4464a2ae208SSatish Balay #undef __FUNCT__ 44749b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 44849b058dcSBarry Smith /*@C 44949b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 45049b058dcSBarry Smith that are currently being perturbed. 45149b058dcSBarry Smith 45249b058dcSBarry Smith Not Collective 45349b058dcSBarry Smith 45449b058dcSBarry Smith Input Parameters: 45549b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 45649b058dcSBarry Smith 45749b058dcSBarry Smith Output Parameters: 45849b058dcSBarry Smith + n - the number of local columns being perturbed 45949b058dcSBarry Smith - cols - the column indices, in global numbering 46049b058dcSBarry Smith 46149b058dcSBarry Smith Level: intermediate 46249b058dcSBarry Smith 46349b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 46449b058dcSBarry Smith 46549b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 46649b058dcSBarry Smith @*/ 467be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 46849b058dcSBarry Smith { 46949b058dcSBarry Smith PetscFunctionBegin; 47049b058dcSBarry Smith if (coloring->currentcolor >= 0) { 47149b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 47249b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 47349b058dcSBarry Smith } else { 47449b058dcSBarry Smith *n = 0; 47549b058dcSBarry Smith } 47649b058dcSBarry Smith PetscFunctionReturn(0); 47749b058dcSBarry Smith } 47849b058dcSBarry Smith 479b8f8c88eSHong Zhang #include "petscda.h" /*I "petscda.h" I*/ 480b8f8c88eSHong Zhang 48149b058dcSBarry Smith #undef __FUNCT__ 4824a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 48343a90d84SBarry Smith /*@ 484e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 485e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 48643a90d84SBarry Smith 487fee21e36SBarry Smith Collective on MatFDColoring 488fee21e36SBarry Smith 489ef5ee4d1SLois Curfman McInnes Input Parameters: 490ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 491ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 492ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 4937850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 494ef5ee4d1SLois Curfman McInnes 4958bba8e72SBarry Smith Options Database Keys: 496*ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 497b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 498b9722508SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 499b9722508SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 5008bba8e72SBarry Smith 50115091d37SBarry Smith Level: intermediate 50215091d37SBarry Smith 5037850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 50443a90d84SBarry Smith 50543a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 50643a90d84SBarry Smith @*/ 507be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 50843a90d84SBarry Smith { 5096849ba73SBarry Smith PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; 5106849ba73SBarry Smith PetscErrorCode ierr; 511b8f8c88eSHong Zhang PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 512efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 513ea709b57SSatish Balay PetscScalar *vscale_array; 514efb30889SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 515b8f8c88eSHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 516005c665bSBarry Smith void *fctx = coloring->fctx; 517f1af5d2fSBarry Smith PetscTruth flg; 518705d48f7SSatish Balay PetscInt ctype=coloring->ctype,N,col_start=0,col_end=0; 519b8f8c88eSHong Zhang Vec x1_tmp; 520a2e34c3dSBarry Smith 5213a40ed3dSBarry Smith PetscFunctionBegin; 5224482741eSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE,1); 5234482741eSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 5244482741eSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE,3); 525feba9168SBarry Smith if (!f) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 526e0907662SLois Curfman McInnes 527d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 5284c49b128SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 529b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 530f1af5d2fSBarry Smith if (flg) { 531ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 532e0907662SLois Curfman McInnes } else { 5330b9b6f31SBarry Smith PetscTruth assembled; 5340b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 5350b9b6f31SBarry Smith if (assembled) { 53643a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 537e0907662SLois Curfman McInnes } 5380b9b6f31SBarry Smith } 53943a90d84SBarry Smith 540b8f8c88eSHong Zhang x1_tmp = x1; 541dac7f5fdSBarry Smith if (!coloring->vscale){ 542b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 543b8f8c88eSHong Zhang } 544be2fbe1fSBarry Smith 5453a7fca6bSBarry Smith /* 5463a7fca6bSBarry Smith This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets 5473a7fca6bSBarry Smith coloring->F for the coarser grids from the finest 5483a7fca6bSBarry Smith */ 5493a7fca6bSBarry Smith if (coloring->F) { 5503a7fca6bSBarry Smith ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 5513a7fca6bSBarry Smith ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 5523a7fca6bSBarry Smith if (m1 != m2) { 5533a7fca6bSBarry Smith coloring->F = 0; 5543a7fca6bSBarry Smith } 5553a7fca6bSBarry Smith } 5563a7fca6bSBarry Smith 557b8f8c88eSHong Zhang if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 558b8f8c88eSHong Zhang ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 559b8f8c88eSHong Zhang } 560b8f8c88eSHong Zhang ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 561b8f8c88eSHong Zhang 562b8f8c88eSHong Zhang /* Set w1 = F(x1) */ 563be2fbe1fSBarry Smith if (coloring->F) { 564be2fbe1fSBarry Smith w1 = coloring->F; /* use already computed value of function */ 565be2fbe1fSBarry Smith coloring->F = 0; 566be2fbe1fSBarry Smith } else { 56766f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 568b8f8c88eSHong Zhang ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 56966f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 570be2fbe1fSBarry Smith } 57143a90d84SBarry Smith 572b8f8c88eSHong Zhang if (!coloring->w3) { 573b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 574b8f8c88eSHong Zhang ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 575efb30889SBarry Smith } 576b8f8c88eSHong Zhang w3 = coloring->w3; 577efb30889SBarry Smith 578b8f8c88eSHong Zhang /* Compute all the local scale factors, including ghost points */ 579b8f8c88eSHong Zhang ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 580b8f8c88eSHong Zhang ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 581b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 582b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED){ 583b8f8c88eSHong Zhang col_start = 0; col_end = N; 5848ee2e534SBarry Smith } else if (ctype == IS_COLORING_GLOBAL){ 585b8f8c88eSHong Zhang xx = xx - start; 586b8f8c88eSHong Zhang vscale_array = vscale_array - start; 587b8f8c88eSHong Zhang col_start = start; col_end = N + start; 588b8f8c88eSHong Zhang } 589b8f8c88eSHong Zhang for (col=col_start; col<col_end; col++){ 590b8f8c88eSHong Zhang /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 591efb30889SBarry Smith if (coloring->htype[0] == 'w') { 592efb30889SBarry Smith dx = 1.0 + unorm; 593efb30889SBarry Smith } else { 5949782cf98SBarry Smith dx = xx[col]; 595efb30889SBarry Smith } 5965b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 5979782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX) 5989782cf98SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 5999782cf98SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 6009782cf98SBarry Smith #else 6019782cf98SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 6029782cf98SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 6039782cf98SBarry Smith #endif 6049782cf98SBarry Smith dx *= epsilon; 60530b34957SBarry Smith vscale_array[col] = 1.0/dx; 6069782cf98SBarry Smith } 6078ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start; 608efb30889SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 6098ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL){ 61030b34957SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61130b34957SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 612b8f8c88eSHong Zhang } 6139782cf98SBarry Smith 614b8f8c88eSHong Zhang if (coloring->vscaleforrow) { 615b8f8c88eSHong Zhang vscaleforrow = coloring->vscaleforrow; 616b8f8c88eSHong Zhang } else { 617b8f8c88eSHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 618b8f8c88eSHong Zhang } 619b0a32e0cSBarry Smith 6209782cf98SBarry Smith /* 62143a90d84SBarry Smith Loop over each color 62243a90d84SBarry Smith */ 623b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 62443a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 62549b058dcSBarry Smith coloring->currentcolor = k; 626b8f8c88eSHong Zhang ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 627b8f8c88eSHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 6288ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 62943a90d84SBarry Smith /* 630b8f8c88eSHong Zhang Loop over each column associated with color 631b8f8c88eSHong Zhang adding the perturbation to the vector w3. 63243a90d84SBarry Smith */ 63343a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 634b8f8c88eSHong Zhang col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 635efb30889SBarry Smith if (coloring->htype[0] == 'w') { 636efb30889SBarry Smith dx = 1.0 + unorm; 637efb30889SBarry Smith } else { 63842460c72SBarry Smith dx = xx[col]; 639efb30889SBarry Smith } 6405b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 641aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 642ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 643ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 64443a90d84SBarry Smith #else 645329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 646329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 64743a90d84SBarry Smith #endif 64843a90d84SBarry Smith dx *= epsilon; 6491302d50aSBarry Smith if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 65042460c72SBarry Smith w3_array[col] += dx; 65143a90d84SBarry Smith } 6528ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 653b8f8c88eSHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 6543b28642cSBarry Smith 65543a90d84SBarry Smith /* 656b8f8c88eSHong Zhang Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 657b8f8c88eSHong Zhang w2 = F(x1 + dx) - F(x1) 65843a90d84SBarry Smith */ 65966f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 66043a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 66166f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 662efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 6639782cf98SBarry Smith 66443a90d84SBarry Smith /* 665e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 66643a90d84SBarry Smith */ 6673b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 66843a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 669b8f8c88eSHong Zhang row = coloring->rows[k][l]; /* local row index */ 670b8f8c88eSHong Zhang col = coloring->columnsforrow[k][l]; /* global column index */ 6715904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 67243a90d84SBarry Smith srow = row + start; 67343a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 67443a90d84SBarry Smith } 6753b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 676aeb7e98aSMatthew Knepley } /* endof for each color */ 6778ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) xx = xx + start; 67830b34957SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 679b8f8c88eSHong Zhang ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 680b8f8c88eSHong Zhang 681b8f8c88eSHong Zhang coloring->currentcolor = -1; 68243a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68343a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 684d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 685a2e34c3dSBarry Smith 686b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 687a2e34c3dSBarry Smith if (flg) { 688a2e34c3dSBarry Smith ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 689a2e34c3dSBarry Smith } 690b9722508SLois Curfman McInnes ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 6913a40ed3dSBarry Smith PetscFunctionReturn(0); 69243a90d84SBarry Smith } 693840b8ebdSBarry Smith 6944a2ae208SSatish Balay #undef __FUNCT__ 6954a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS" 696840b8ebdSBarry Smith /*@ 697840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 698840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 699840b8ebdSBarry Smith 700fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 701fee21e36SBarry Smith 702ef5ee4d1SLois Curfman McInnes Input Parameters: 7033b28642cSBarry Smith + mat - location to store Jacobian 704ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 705ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 7067850c7c0SBarry Smith - sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null 707ef5ee4d1SLois Curfman McInnes 70815091d37SBarry Smith Level: intermediate 70915091d37SBarry Smith 7107850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 711840b8ebdSBarry Smith 712840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 713840b8ebdSBarry Smith @*/ 714be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 715840b8ebdSBarry Smith { 7166849ba73SBarry Smith PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f; 7176849ba73SBarry Smith PetscErrorCode ierr; 71838baddfdSBarry Smith PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow; 719efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 720ea709b57SSatish Balay PetscScalar *vscale_array; 721329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 722ced766e8SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 723840b8ebdSBarry Smith void *fctx = coloring->fctx; 724f1af5d2fSBarry Smith PetscTruth flg; 725840b8ebdSBarry Smith 7263a40ed3dSBarry Smith PetscFunctionBegin; 7274482741eSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE,1); 7284482741eSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE,2); 7294482741eSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE,4); 730840b8ebdSBarry Smith 731d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 732ced766e8SHong Zhang if (!coloring->w3) { 733840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 73452e6d16bSBarry Smith ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 735840b8ebdSBarry Smith } 736ced766e8SHong Zhang w3 = coloring->w3; 737840b8ebdSBarry Smith 7385904e1b1SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 739b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 740f1af5d2fSBarry Smith if (flg) { 741ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 742840b8ebdSBarry Smith } else { 7430b9b6f31SBarry Smith PetscTruth assembled; 7440b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 7450b9b6f31SBarry Smith if (assembled) { 746840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 747840b8ebdSBarry Smith } 7480b9b6f31SBarry Smith } 749840b8ebdSBarry Smith 750840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 751840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 75266f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 753840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 75466f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 755840b8ebdSBarry Smith 756840b8ebdSBarry Smith /* 7575904e1b1SBarry Smith Compute all the scale factors and share with other processors 758840b8ebdSBarry Smith */ 7595904e1b1SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 7605904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 761840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 762840b8ebdSBarry Smith /* 763840b8ebdSBarry Smith Loop over each column associated with color adding the 764840b8ebdSBarry Smith perturbation to the vector w3. 765840b8ebdSBarry Smith */ 766840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 767840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 7685904e1b1SBarry Smith dx = xx[col]; 7695b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 770aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 771840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 772840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 773840b8ebdSBarry Smith #else 774329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 775329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 776840b8ebdSBarry Smith #endif 777840b8ebdSBarry Smith dx *= epsilon; 7785904e1b1SBarry Smith vscale_array[col] = 1.0/dx; 779840b8ebdSBarry Smith } 7805904e1b1SBarry Smith } 7815904e1b1SBarry Smith vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7825904e1b1SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7835904e1b1SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7845904e1b1SBarry Smith 7855904e1b1SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 7865904e1b1SBarry Smith else vscaleforrow = coloring->columnsforrow; 7875904e1b1SBarry Smith 7885904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7895904e1b1SBarry Smith /* 7905904e1b1SBarry Smith Loop over each color 7915904e1b1SBarry Smith */ 7925904e1b1SBarry Smith for (k=0; k<coloring->ncolors; k++) { 7935904e1b1SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 7945904e1b1SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 7955904e1b1SBarry Smith /* 7965904e1b1SBarry Smith Loop over each column associated with color adding the 7975904e1b1SBarry Smith perturbation to the vector w3. 7985904e1b1SBarry Smith */ 7995904e1b1SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 8005904e1b1SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 8015904e1b1SBarry Smith dx = xx[col]; 8025b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 8035904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 8045904e1b1SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 8055904e1b1SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 8065904e1b1SBarry Smith #else 8075904e1b1SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 8085904e1b1SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 8095904e1b1SBarry Smith #endif 8105904e1b1SBarry Smith dx *= epsilon; 8115904e1b1SBarry Smith w3_array[col] += dx; 8125904e1b1SBarry Smith } 8135904e1b1SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 8145904e1b1SBarry Smith 815840b8ebdSBarry Smith /* 816840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 817840b8ebdSBarry Smith */ 81866f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 819840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 82066f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 821efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 8225904e1b1SBarry Smith 823840b8ebdSBarry Smith /* 824840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 825840b8ebdSBarry Smith */ 8263b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 827840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 828840b8ebdSBarry Smith row = coloring->rows[k][l]; 829840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 8305904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 831840b8ebdSBarry Smith srow = row + start; 832840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 833840b8ebdSBarry Smith } 8343b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 835840b8ebdSBarry Smith } 8365904e1b1SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8375904e1b1SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 838840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 839840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 840d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 8413a40ed3dSBarry Smith PetscFunctionReturn(0); 842840b8ebdSBarry Smith } 8433b28642cSBarry Smith 8443b28642cSBarry Smith 8453b28642cSBarry Smith 846