1*76d8deaeSBarry Smith /*$Id: fdmatrix.c,v 1.85 2001/04/10 19:35:59 bsmith Exp bsmith $*/ 2bbf0e169SBarry Smith 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 8bbf0e169SBarry Smith #include "petsc.h" 9e090d566SSatish Balay #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 10bbf0e169SBarry Smith 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 13b0a32e0cSBarry Smith static int MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 149194eea9SBarry Smith { 159194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 169194eea9SBarry Smith int ierr,i,j; 179194eea9SBarry Smith PetscReal x,y; 189194eea9SBarry Smith 199194eea9SBarry Smith PetscFunctionBegin; 209194eea9SBarry Smith 219194eea9SBarry Smith /* loop over colors */ 229194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 239194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 249194eea9SBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 259194eea9SBarry Smith x = fd->columnsforrow[i][j]; 26b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 279194eea9SBarry Smith } 289194eea9SBarry Smith } 299194eea9SBarry Smith PetscFunctionReturn(0); 309194eea9SBarry Smith } 319194eea9SBarry Smith 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 34b0a32e0cSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 35005c665bSBarry Smith { 369194eea9SBarry Smith int ierr; 37005c665bSBarry Smith PetscTruth isnull; 38b0a32e0cSBarry Smith PetscDraw draw; 3954d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 40005c665bSBarry Smith 413a40ed3dSBarry Smith PetscFunctionBegin; 42b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 43b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 449194eea9SBarry Smith 459194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 46005c665bSBarry Smith 47005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 48005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 49b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 50b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 519194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 523a40ed3dSBarry Smith PetscFunctionReturn(0); 53005c665bSBarry Smith } 54005c665bSBarry Smith 554a2ae208SSatish Balay #undef __FUNCT__ 564a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 57bbf0e169SBarry Smith /*@C 58639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 59bbf0e169SBarry Smith 604c49b128SBarry Smith Collective on MatFDColoring 61fee21e36SBarry Smith 62ef5ee4d1SLois Curfman McInnes Input Parameters: 63ef5ee4d1SLois Curfman McInnes + c - the coloring context 64ef5ee4d1SLois Curfman McInnes - viewer - visualization context 65ef5ee4d1SLois Curfman McInnes 6615091d37SBarry Smith Level: intermediate 6715091d37SBarry Smith 68b4fc646aSLois Curfman McInnes Notes: 69b4fc646aSLois Curfman McInnes The available visualization contexts include 70b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 71b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 72ef5ee4d1SLois Curfman McInnes output where only the first processor opens 73ef5ee4d1SLois Curfman McInnes the file. All other processors send their 74ef5ee4d1SLois Curfman McInnes data to the first processor to print. 75b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 76bbf0e169SBarry Smith 779194eea9SBarry Smith Notes: 789194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 79b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 809194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 819194eea9SBarry Smith 82639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 83005c665bSBarry Smith 84b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 85bbf0e169SBarry Smith @*/ 86b0a32e0cSBarry Smith int MatFDColoringView(MatFDColoring c,PetscViewer viewer) 87bbf0e169SBarry Smith { 88fb9695e5SSatish Balay int i,j,ierr; 896831982aSBarry Smith PetscTruth isdraw,isascii; 90f3ef73ceSBarry Smith PetscViewerFormat format; 91bbf0e169SBarry Smith 923a40ed3dSBarry Smith PetscFunctionBegin; 93b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 94b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 95b0a32e0cSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 966831982aSBarry Smith PetscCheckSameComm(c,viewer); 97bbf0e169SBarry Smith 98fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 99b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 1000f5bd95cSBarry Smith if (isdraw) { 101b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 1020f5bd95cSBarry Smith } else if (isascii) { 103b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 104b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr); 105b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",c->umin);CHKERRQ(ierr); 106b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors);CHKERRQ(ierr); 107ae09f205SBarry Smith 108b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 109fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 110b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 111b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %d\n",i);CHKERRQ(ierr); 112b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr); 113b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 114b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]);CHKERRQ(ierr); 115639f9d9dSBarry Smith } 116b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr); 117b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 118b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 119b4fc646aSLois Curfman McInnes } 120bbf0e169SBarry Smith } 121bbf0e169SBarry Smith } 122b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1235cd90555SBarry Smith } else { 12429bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 125bbf0e169SBarry Smith } 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 127639f9d9dSBarry Smith } 128639f9d9dSBarry Smith 1294a2ae208SSatish Balay #undef __FUNCT__ 1304a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 131639f9d9dSBarry Smith /*@ 132b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 133b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 134639f9d9dSBarry Smith 135ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 136ef5ee4d1SLois Curfman McInnes 137ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 138ef5ee4d1SLois Curfman McInnes .vb 13965f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 140f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 141f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 142ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 143ef5ee4d1SLois Curfman McInnes .ve 144639f9d9dSBarry Smith 145639f9d9dSBarry Smith Input Parameters: 146ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 147639f9d9dSBarry Smith . error_rel - relative error 148f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 149fee21e36SBarry Smith 15015091d37SBarry Smith Level: advanced 15115091d37SBarry Smith 152b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 153b4fc646aSLois Curfman McInnes 154b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 155639f9d9dSBarry Smith @*/ 156329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 157639f9d9dSBarry Smith { 1583a40ed3dSBarry Smith PetscFunctionBegin; 159639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 160639f9d9dSBarry Smith 161639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 162639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 164639f9d9dSBarry Smith } 165639f9d9dSBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFrequency" 168005c665bSBarry Smith /*@ 169e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 170e0907662SLois Curfman McInnes matrices. 171005c665bSBarry Smith 172fee21e36SBarry Smith Collective on MatFDColoring 173fee21e36SBarry Smith 174ef5ee4d1SLois Curfman McInnes Input Parameters: 175ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 176ef5ee4d1SLois Curfman McInnes - freq - frequency (default is 1) 177ef5ee4d1SLois Curfman McInnes 17815091d37SBarry Smith Options Database Keys: 17915091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 18015091d37SBarry Smith 18115091d37SBarry Smith Level: advanced 18215091d37SBarry Smith 183e0907662SLois Curfman McInnes Notes: 184e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 185e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 186e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 187e0907662SLois Curfman McInnes <freq> nonlinear iterations. 188e0907662SLois Curfman McInnes 189b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 190ef5ee4d1SLois Curfman McInnes 19140964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute() 192005c665bSBarry Smith @*/ 193005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 194005c665bSBarry Smith { 1953a40ed3dSBarry Smith PetscFunctionBegin; 196005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 197005c665bSBarry Smith 198005c665bSBarry Smith matfd->freq = freq; 1993a40ed3dSBarry Smith PetscFunctionReturn(0); 200005c665bSBarry Smith } 201005c665bSBarry Smith 2024a2ae208SSatish Balay #undef __FUNCT__ 2034a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringGetFrequency" 204ff0cfa39SBarry Smith /*@ 205ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 206ff0cfa39SBarry Smith matrices. 207ff0cfa39SBarry Smith 208ef5ee4d1SLois Curfman McInnes Not Collective 209ef5ee4d1SLois Curfman McInnes 210ff0cfa39SBarry Smith Input Parameters: 211ff0cfa39SBarry Smith . coloring - the coloring context 212ff0cfa39SBarry Smith 213ff0cfa39SBarry Smith Output Parameters: 214ff0cfa39SBarry Smith . freq - frequency (default is 1) 215ff0cfa39SBarry Smith 21615091d37SBarry Smith Options Database Keys: 21715091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 21815091d37SBarry Smith 21915091d37SBarry Smith Level: advanced 22015091d37SBarry Smith 221ff0cfa39SBarry Smith Notes: 222ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 223ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 224ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 225ff0cfa39SBarry Smith <freq> nonlinear iterations. 226ff0cfa39SBarry Smith 227ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 228ef5ee4d1SLois Curfman McInnes 229ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency() 230ff0cfa39SBarry Smith @*/ 231ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 232ff0cfa39SBarry Smith { 2333a40ed3dSBarry Smith PetscFunctionBegin; 234ff0cfa39SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 235ff0cfa39SBarry Smith 236ff0cfa39SBarry Smith *freq = matfd->freq; 2373a40ed3dSBarry Smith PetscFunctionReturn(0); 238ff0cfa39SBarry Smith } 239ff0cfa39SBarry Smith 2404a2ae208SSatish Balay #undef __FUNCT__ 2414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 242d64ed03dSBarry Smith /*@C 243005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 244005c665bSBarry Smith 245fee21e36SBarry Smith Collective on MatFDColoring 246fee21e36SBarry Smith 247ef5ee4d1SLois Curfman McInnes Input Parameters: 248ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 249ef5ee4d1SLois Curfman McInnes . f - the function 250ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 251ef5ee4d1SLois Curfman McInnes 25215091d37SBarry Smith Level: intermediate 25315091d37SBarry Smith 254f881d145SBarry Smith Notes: 255f881d145SBarry Smith In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 256f881d145SBarry Smith be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 257f881d145SBarry Smith with the TS solvers. 258f881d145SBarry Smith 259b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 260005c665bSBarry Smith @*/ 261840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 262005c665bSBarry Smith { 2633a40ed3dSBarry Smith PetscFunctionBegin; 264005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 265005c665bSBarry Smith 266005c665bSBarry Smith matfd->f = f; 267005c665bSBarry Smith matfd->fctx = fctx; 268005c665bSBarry Smith 2693a40ed3dSBarry Smith PetscFunctionReturn(0); 270005c665bSBarry Smith } 271005c665bSBarry Smith 2724a2ae208SSatish Balay #undef __FUNCT__ 2734a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 274639f9d9dSBarry Smith /*@ 275b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 276639f9d9dSBarry Smith the options database. 277639f9d9dSBarry Smith 278fee21e36SBarry Smith Collective on MatFDColoring 279fee21e36SBarry Smith 28065f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 281ef5ee4d1SLois Curfman McInnes .vb 28265f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 283f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 284f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 285ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 286ef5ee4d1SLois Curfman McInnes .ve 287ef5ee4d1SLois Curfman McInnes 288ef5ee4d1SLois Curfman McInnes Input Parameter: 289ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 290ef5ee4d1SLois Curfman McInnes 291b4fc646aSLois Curfman McInnes Options Database Keys: 292d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 293ef5ee4d1SLois Curfman McInnes of relative error in the function) 294f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 295ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 296ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 297ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 298ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 299639f9d9dSBarry Smith 30015091d37SBarry Smith Level: intermediate 30115091d37SBarry Smith 302b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 303639f9d9dSBarry Smith @*/ 304639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 305639f9d9dSBarry Smith { 306186905e3SBarry Smith int ierr; 3073a40ed3dSBarry Smith 3083a40ed3dSBarry Smith PetscFunctionBegin; 309639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 310639f9d9dSBarry Smith 311b0a32e0cSBarry Smith ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 312b0a32e0cSBarry Smith ierr = PetscOptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 313b0a32e0cSBarry Smith ierr = PetscOptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 314b0a32e0cSBarry Smith ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 315186905e3SBarry Smith /* not used here; but so they are presented in the GUI */ 316b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 317b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 318b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 319b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 321005c665bSBarry Smith } 322005c665bSBarry Smith 3234a2ae208SSatish Balay #undef __FUNCT__ 3244a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private" 325005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 326005c665bSBarry Smith { 327f1af5d2fSBarry Smith int ierr; 328f1af5d2fSBarry Smith PetscTruth flg; 329005c665bSBarry Smith 3303a40ed3dSBarry Smith PetscFunctionBegin; 331b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 332005c665bSBarry Smith if (flg) { 333b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 334005c665bSBarry Smith } 335b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 336ae09f205SBarry Smith if (flg) { 337fb9695e5SSatish Balay ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 338b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339b0a32e0cSBarry Smith ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 340ae09f205SBarry Smith } 341b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 342005c665bSBarry Smith if (flg) { 343b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 344b0a32e0cSBarry Smith ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 345005c665bSBarry Smith } 3463a40ed3dSBarry Smith PetscFunctionReturn(0); 347bbf0e169SBarry Smith } 348bbf0e169SBarry Smith 3494a2ae208SSatish Balay #undef __FUNCT__ 3504a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 351bbf0e169SBarry Smith /*@C 352639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 353639f9d9dSBarry Smith computation of Jacobians. 354bbf0e169SBarry Smith 355ef5ee4d1SLois Curfman McInnes Collective on Mat 356ef5ee4d1SLois Curfman McInnes 357639f9d9dSBarry Smith Input Parameters: 358ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 359ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 360bbf0e169SBarry Smith 361bbf0e169SBarry Smith Output Parameter: 362639f9d9dSBarry Smith . color - the new coloring context 363bbf0e169SBarry Smith 364b4fc646aSLois Curfman McInnes Options Database Keys: 365ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_view - Activates basic viewing or coloring 366ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 367ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 368639f9d9dSBarry Smith 36915091d37SBarry Smith Level: intermediate 37015091d37SBarry Smith 3716831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 3726831982aSBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 373bbf0e169SBarry Smith @*/ 374639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 375bbf0e169SBarry Smith { 376639f9d9dSBarry Smith MatFDColoring c; 377639f9d9dSBarry Smith MPI_Comm comm; 378639f9d9dSBarry Smith int ierr,M,N; 379639f9d9dSBarry Smith 3803a40ed3dSBarry Smith PetscFunctionBegin; 381b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 382639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 38329bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 384639f9d9dSBarry Smith 385f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 3869194eea9SBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 387b0a32e0cSBarry Smith PetscLogObjectCreate(c); 388639f9d9dSBarry Smith 389f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 390f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 391639f9d9dSBarry Smith } else { 39229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 393639f9d9dSBarry Smith } 394639f9d9dSBarry Smith 395639f9d9dSBarry Smith c->error_rel = 1.e-8; 396ae09f205SBarry Smith c->umin = 1.e-6; 397005c665bSBarry Smith c->freq = 1; 39840964ad5SBarry Smith c->usersetsrecompute = PETSC_FALSE; 39940964ad5SBarry Smith c->recompute = PETSC_FALSE; 400005c665bSBarry Smith 401005c665bSBarry Smith ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 402639f9d9dSBarry Smith 403639f9d9dSBarry Smith *color = c; 404b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4053a40ed3dSBarry Smith PetscFunctionReturn(0); 406bbf0e169SBarry Smith } 407bbf0e169SBarry Smith 4084a2ae208SSatish Balay #undef __FUNCT__ 4094a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 410bbf0e169SBarry Smith /*@C 411639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 412639f9d9dSBarry Smith via MatFDColoringCreate(). 413bbf0e169SBarry Smith 414ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 415ef5ee4d1SLois Curfman McInnes 416b4fc646aSLois Curfman McInnes Input Parameter: 417639f9d9dSBarry Smith . c - coloring context 418bbf0e169SBarry Smith 41915091d37SBarry Smith Level: intermediate 42015091d37SBarry Smith 421639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 422bbf0e169SBarry Smith @*/ 423639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 424bbf0e169SBarry Smith { 425263760aaSBarry Smith int i,ierr; 426bbf0e169SBarry Smith 4273a40ed3dSBarry Smith PetscFunctionBegin; 4283a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 429d4bb536fSBarry Smith 430639f9d9dSBarry Smith for (i=0; i<c->ncolors; i++) { 431606d414cSSatish Balay if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 432606d414cSSatish Balay if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 433606d414cSSatish Balay if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 43430b34957SBarry Smith if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 435bbf0e169SBarry Smith } 436606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 437606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 438606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 439606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 440606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 44130b34957SBarry Smith if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 4424f113abfSBarry Smith if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 443005c665bSBarry Smith if (c->w1) { 444005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 445005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 446005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 447005c665bSBarry Smith } 448b0a32e0cSBarry Smith PetscLogObjectDestroy(c); 449639f9d9dSBarry Smith PetscHeaderDestroy(c); 4503a40ed3dSBarry Smith PetscFunctionReturn(0); 451bbf0e169SBarry Smith } 45243a90d84SBarry Smith 4534a2ae208SSatish Balay #undef __FUNCT__ 4544a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 45543a90d84SBarry Smith /*@ 456e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 457e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 45843a90d84SBarry Smith 459fee21e36SBarry Smith Collective on MatFDColoring 460fee21e36SBarry Smith 461ef5ee4d1SLois Curfman McInnes Input Parameters: 462ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 463ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 464ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 465ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 466ef5ee4d1SLois Curfman McInnes 4678bba8e72SBarry Smith Options Database Keys: 468ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 4698bba8e72SBarry Smith 47015091d37SBarry Smith Level: intermediate 47115091d37SBarry Smith 47243a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 47343a90d84SBarry Smith 47443a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 47543a90d84SBarry Smith @*/ 476005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 47743a90d84SBarry Smith { 4785904e1b1SBarry Smith int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f; 4795904e1b1SBarry Smith int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 4805904e1b1SBarry Smith Scalar dx,mone = -1.0,*y,*xx,*w3_array; 4814f113abfSBarry Smith Scalar *vscale_array; 482329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 483005c665bSBarry Smith Vec w1,w2,w3; 484005c665bSBarry Smith void *fctx = coloring->fctx; 485f1af5d2fSBarry Smith PetscTruth flg; 486005c665bSBarry Smith 487a2e34c3dSBarry Smith 4883a40ed3dSBarry Smith PetscFunctionBegin; 489e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 490e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 491e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 492e0907662SLois Curfman McInnes 49340964ad5SBarry Smith if (coloring->usersetsrecompute) { 49440964ad5SBarry Smith if (!coloring->recompute) { 49540964ad5SBarry Smith *flag = SAME_PRECONDITIONER; 496*76d8deaeSBarry Smith PetscLogInfo(J,"MatFDColoringApply:Skipping Jacobian, since user called MatFDColorSetRecompute()\n"); 49740964ad5SBarry Smith PetscFunctionReturn(0); 49840964ad5SBarry Smith } else { 49940964ad5SBarry Smith coloring->recompute = PETSC_FALSE; 50040964ad5SBarry Smith } 50140964ad5SBarry Smith } 50240964ad5SBarry Smith 503b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 504*76d8deaeSBarry Smith if (J->ops->fdcoloringapply) { 505*76d8deaeSBarry Smith ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 506*76d8deaeSBarry Smith } else { 507*76d8deaeSBarry Smith 508005c665bSBarry Smith if (!coloring->w1) { 509005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 510b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w1); 511005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 512b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w2); 513005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 514b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w3); 515005c665bSBarry Smith } 516005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 51743a90d84SBarry Smith 5184c49b128SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 519b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 520f1af5d2fSBarry Smith if (flg) { 521b0a32e0cSBarry Smith PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 522e0907662SLois Curfman McInnes } else { 52343a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 524e0907662SLois Curfman McInnes } 52543a90d84SBarry Smith 52643a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 52743a90d84SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 52843a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 52943a90d84SBarry Smith 53043a90d84SBarry Smith /* 5319782cf98SBarry Smith Compute all the scale factors and share with other processors 5329782cf98SBarry Smith */ 5339782cf98SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 5344f113abfSBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 5359782cf98SBarry Smith for (k=0; k<coloring->ncolors; k++) { 5369782cf98SBarry Smith /* 5379782cf98SBarry Smith Loop over each column associated with color adding the 5389782cf98SBarry Smith perturbation to the vector w3. 5399782cf98SBarry Smith */ 5409782cf98SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 5419782cf98SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 5429782cf98SBarry Smith dx = xx[col]; 5439782cf98SBarry Smith if (dx == 0.0) dx = 1.0; 5449782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX) 5459782cf98SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 5469782cf98SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 5479782cf98SBarry Smith #else 5489782cf98SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 5499782cf98SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 5509782cf98SBarry Smith #endif 5519782cf98SBarry Smith dx *= epsilon; 55230b34957SBarry Smith vscale_array[col] = 1.0/dx; 5539782cf98SBarry Smith } 5549782cf98SBarry Smith } 555a2e34c3dSBarry Smith vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 55630b34957SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 55730b34957SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5589782cf98SBarry Smith 559ce1411ecSBarry Smith /* ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); 560ce1411ecSBarry Smith ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/ 561b0a32e0cSBarry Smith 56230b34957SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 56330b34957SBarry Smith else vscaleforrow = coloring->columnsforrow; 56430b34957SBarry Smith 56530b34957SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 5669782cf98SBarry Smith /* 56743a90d84SBarry Smith Loop over each color 56843a90d84SBarry Smith */ 56943a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 57043a90d84SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 57142460c72SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 57243a90d84SBarry Smith /* 57343a90d84SBarry Smith Loop over each column associated with color adding the 57443a90d84SBarry Smith perturbation to the vector w3. 57543a90d84SBarry Smith */ 57643a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 57743a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 57842460c72SBarry Smith dx = xx[col]; 579ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 580aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 581ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 582ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 58343a90d84SBarry Smith #else 584329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 585329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 58643a90d84SBarry Smith #endif 58743a90d84SBarry Smith dx *= epsilon; 58829bbc08cSBarry Smith if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter"); 58942460c72SBarry Smith w3_array[col] += dx; 59043a90d84SBarry Smith } 59142460c72SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 5923b28642cSBarry Smith 59343a90d84SBarry Smith /* 594e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 59543a90d84SBarry Smith */ 596a2e34c3dSBarry Smith 59743a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 59843a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 5999782cf98SBarry Smith 60043a90d84SBarry Smith /* 601e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 60243a90d84SBarry Smith */ 6033b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 60443a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 60543a90d84SBarry Smith row = coloring->rows[k][l]; 60643a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 6075904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 60843a90d84SBarry Smith srow = row + start; 60943a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 61043a90d84SBarry Smith } 6113b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 61243a90d84SBarry Smith } 61330b34957SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 61442460c72SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 61543a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61643a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 617*76d8deaeSBarry Smith } 618b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 619a2e34c3dSBarry Smith 620b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 621a2e34c3dSBarry Smith if (flg) { 622a2e34c3dSBarry Smith ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 623a2e34c3dSBarry Smith } 6243a40ed3dSBarry Smith PetscFunctionReturn(0); 62543a90d84SBarry Smith } 626840b8ebdSBarry Smith 6274a2ae208SSatish Balay #undef __FUNCT__ 6284a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS" 629840b8ebdSBarry Smith /*@ 630840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 631840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 632840b8ebdSBarry Smith 633fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 634fee21e36SBarry Smith 635ef5ee4d1SLois Curfman McInnes Input Parameters: 6363b28642cSBarry Smith + mat - location to store Jacobian 637ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 638ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 639ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 640ef5ee4d1SLois Curfman McInnes 641840b8ebdSBarry Smith Options Database Keys: 642ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 643840b8ebdSBarry Smith 64415091d37SBarry Smith Level: intermediate 64515091d37SBarry Smith 646840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 647840b8ebdSBarry Smith 648840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 649840b8ebdSBarry Smith @*/ 650329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 651840b8ebdSBarry Smith { 652329f5518SBarry Smith int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 6535904e1b1SBarry Smith int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 6545904e1b1SBarry Smith Scalar dx,mone = -1.0,*y,*xx,*w3_array; 6555904e1b1SBarry Smith Scalar *vscale_array; 656329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 657840b8ebdSBarry Smith Vec w1,w2,w3; 658840b8ebdSBarry Smith void *fctx = coloring->fctx; 659f1af5d2fSBarry Smith PetscTruth flg; 660840b8ebdSBarry Smith 6613a40ed3dSBarry Smith PetscFunctionBegin; 662840b8ebdSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE); 663840b8ebdSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 664840b8ebdSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE); 665840b8ebdSBarry Smith 666b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 667840b8ebdSBarry Smith if (!coloring->w1) { 668840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 669b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w1); 670840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 671b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w2); 672840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 673b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w3); 674840b8ebdSBarry Smith } 675840b8ebdSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 676840b8ebdSBarry Smith 6775904e1b1SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 678b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 679f1af5d2fSBarry Smith if (flg) { 680b0a32e0cSBarry Smith PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 681840b8ebdSBarry Smith } else { 682840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 683840b8ebdSBarry Smith } 684840b8ebdSBarry Smith 685840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 686840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 687840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 688840b8ebdSBarry Smith 689840b8ebdSBarry Smith /* 6905904e1b1SBarry Smith Compute all the scale factors and share with other processors 691840b8ebdSBarry Smith */ 6925904e1b1SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 6935904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 694840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 695840b8ebdSBarry Smith /* 696840b8ebdSBarry Smith Loop over each column associated with color adding the 697840b8ebdSBarry Smith perturbation to the vector w3. 698840b8ebdSBarry Smith */ 699840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 700840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 7015904e1b1SBarry Smith dx = xx[col]; 702840b8ebdSBarry Smith if (dx == 0.0) dx = 1.0; 703aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 704840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 705840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 706840b8ebdSBarry Smith #else 707329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 708329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 709840b8ebdSBarry Smith #endif 710840b8ebdSBarry Smith dx *= epsilon; 7115904e1b1SBarry Smith vscale_array[col] = 1.0/dx; 712840b8ebdSBarry Smith } 7135904e1b1SBarry Smith } 7145904e1b1SBarry Smith vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7155904e1b1SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7165904e1b1SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7175904e1b1SBarry Smith 7185904e1b1SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 7195904e1b1SBarry Smith else vscaleforrow = coloring->columnsforrow; 7205904e1b1SBarry Smith 7215904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7225904e1b1SBarry Smith /* 7235904e1b1SBarry Smith Loop over each color 7245904e1b1SBarry Smith */ 7255904e1b1SBarry Smith for (k=0; k<coloring->ncolors; k++) { 7265904e1b1SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 7275904e1b1SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 7285904e1b1SBarry Smith /* 7295904e1b1SBarry Smith Loop over each column associated with color adding the 7305904e1b1SBarry Smith perturbation to the vector w3. 7315904e1b1SBarry Smith */ 7325904e1b1SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 7335904e1b1SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 7345904e1b1SBarry Smith dx = xx[col]; 7355904e1b1SBarry Smith if (dx == 0.0) dx = 1.0; 7365904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 7375904e1b1SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 7385904e1b1SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 7395904e1b1SBarry Smith #else 7405904e1b1SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 7415904e1b1SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 7425904e1b1SBarry Smith #endif 7435904e1b1SBarry Smith dx *= epsilon; 7445904e1b1SBarry Smith w3_array[col] += dx; 7455904e1b1SBarry Smith } 7465904e1b1SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 7475904e1b1SBarry Smith 748840b8ebdSBarry Smith /* 749840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 750840b8ebdSBarry Smith */ 751840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 752840b8ebdSBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 7535904e1b1SBarry Smith 754840b8ebdSBarry Smith /* 755840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 756840b8ebdSBarry Smith */ 7573b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 758840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 759840b8ebdSBarry Smith row = coloring->rows[k][l]; 760840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 7615904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 762840b8ebdSBarry Smith srow = row + start; 763840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 764840b8ebdSBarry Smith } 7653b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 766840b8ebdSBarry Smith } 7675904e1b1SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7685904e1b1SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 769840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 770840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 771b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 7723a40ed3dSBarry Smith PetscFunctionReturn(0); 773840b8ebdSBarry Smith } 7743b28642cSBarry Smith 7753b28642cSBarry Smith 7764a2ae208SSatish Balay #undef __FUNCT__ 7774a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetRecompute()" 77840964ad5SBarry Smith /*@C 77940964ad5SBarry Smith MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 78040964ad5SBarry Smith is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 78140964ad5SBarry Smith no additional Jacobian's will be computed (the same one will be used) until this is 78240964ad5SBarry Smith called again. 78340964ad5SBarry Smith 78440964ad5SBarry Smith Collective on MatFDColoring 78540964ad5SBarry Smith 78640964ad5SBarry Smith Input Parameters: 78740964ad5SBarry Smith . c - the coloring context 78840964ad5SBarry Smith 78940964ad5SBarry Smith Level: intermediate 79040964ad5SBarry Smith 79140964ad5SBarry Smith Notes: The MatFDColoringSetFrequency() is ignored once this is called 79240964ad5SBarry Smith 79340964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 79440964ad5SBarry Smith 79540964ad5SBarry Smith .keywords: Mat, finite differences, coloring 79640964ad5SBarry Smith @*/ 79740964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c) 79840964ad5SBarry Smith { 79940964ad5SBarry Smith PetscFunctionBegin; 80040964ad5SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 80140964ad5SBarry Smith c->usersetsrecompute = PETSC_TRUE; 80240964ad5SBarry Smith c->recompute = PETSC_TRUE; 80340964ad5SBarry Smith PetscFunctionReturn(0); 80440964ad5SBarry Smith } 80540964ad5SBarry Smith 8063b28642cSBarry Smith 807