1*b0a32e0cSBarry Smith /*$Id: fdmatrix.c,v 1.79 2001/01/04 22:21:25 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 115615d1e5SSatish Balay #undef __FUNC__ 12*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringView_Draw_Zoom" 13*b0a32e0cSBarry 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]; 26*b0a32e0cSBarry 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 329194eea9SBarry Smith #undef __FUNC__ 33*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringView_Draw" 34*b0a32e0cSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 35005c665bSBarry Smith { 369194eea9SBarry Smith int ierr; 37005c665bSBarry Smith PetscTruth isnull; 38*b0a32e0cSBarry Smith PetscDraw draw; 3954d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 40005c665bSBarry Smith 413a40ed3dSBarry Smith PetscFunctionBegin; 42*b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 43*b0a32e0cSBarry 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; 49*b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 50*b0a32e0cSBarry 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 55005c665bSBarry Smith #undef __FUNC__ 56*b0a32e0cSBarry Smith #define __FUNC__ "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 70*b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 71*b0a32e0cSBarry 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. 75*b0a32e0cSBarry 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 79*b0a32e0cSBarry 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 @*/ 86*b0a32e0cSBarry Smith int MatFDColoringView(MatFDColoring c,PetscViewer viewer) 87bbf0e169SBarry Smith { 88639f9d9dSBarry Smith int i,j,format,ierr; 896831982aSBarry Smith PetscTruth isdraw,isascii; 90bbf0e169SBarry Smith 913a40ed3dSBarry Smith PetscFunctionBegin; 92b4fc646aSLois Curfman McInnes PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 93*b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(c->comm); 94*b0a32e0cSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 956831982aSBarry Smith PetscCheckSameComm(c,viewer); 96bbf0e169SBarry Smith 97*b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 98*b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 990f5bd95cSBarry Smith if (isdraw) { 100b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 1010f5bd95cSBarry Smith } else if (isascii) { 102*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 103*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr); 104*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%g\n",c->umin);CHKERRQ(ierr); 105*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%d\n",c->ncolors);CHKERRQ(ierr); 106ae09f205SBarry Smith 107*b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 108*b0a32e0cSBarry Smith if (format != PETSC_VIEWER_FORMAT_ASCII_INFO) { 109b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 110*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %d\n",i);CHKERRQ(ierr); 111*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr); 112b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 113*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d\n",c->columns[i][j]);CHKERRQ(ierr); 114639f9d9dSBarry Smith } 115*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr); 116b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 117*b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 118b4fc646aSLois Curfman McInnes } 119bbf0e169SBarry Smith } 120bbf0e169SBarry Smith } 121*b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1225cd90555SBarry Smith } else { 12329bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 124bbf0e169SBarry Smith } 1253a40ed3dSBarry Smith PetscFunctionReturn(0); 126639f9d9dSBarry Smith } 127639f9d9dSBarry Smith 1285615d1e5SSatish Balay #undef __FUNC__ 129*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringSetParameters" 130639f9d9dSBarry Smith /*@ 131b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 132b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 133639f9d9dSBarry Smith 134ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 135ef5ee4d1SLois Curfman McInnes 136ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 137ef5ee4d1SLois Curfman McInnes .vb 13865f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 139f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 140f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 141ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 142ef5ee4d1SLois Curfman McInnes .ve 143639f9d9dSBarry Smith 144639f9d9dSBarry Smith Input Parameters: 145ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 146639f9d9dSBarry Smith . error_rel - relative error 147f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 148fee21e36SBarry Smith 14915091d37SBarry Smith Level: advanced 15015091d37SBarry Smith 151b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 152b4fc646aSLois Curfman McInnes 153b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate() 154639f9d9dSBarry Smith @*/ 155329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 156639f9d9dSBarry Smith { 1573a40ed3dSBarry Smith PetscFunctionBegin; 158639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 159639f9d9dSBarry Smith 160639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 161639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1623a40ed3dSBarry Smith PetscFunctionReturn(0); 163639f9d9dSBarry Smith } 164639f9d9dSBarry Smith 1655615d1e5SSatish Balay #undef __FUNC__ 166*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringSetFrequency" 167005c665bSBarry Smith /*@ 168e0907662SLois Curfman McInnes MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian 169e0907662SLois Curfman McInnes matrices. 170005c665bSBarry Smith 171fee21e36SBarry Smith Collective on MatFDColoring 172fee21e36SBarry Smith 173ef5ee4d1SLois Curfman McInnes Input Parameters: 174ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 175ef5ee4d1SLois Curfman McInnes - freq - frequency (default is 1) 176ef5ee4d1SLois Curfman McInnes 17715091d37SBarry Smith Options Database Keys: 17815091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 17915091d37SBarry Smith 18015091d37SBarry Smith Level: advanced 18115091d37SBarry Smith 182e0907662SLois Curfman McInnes Notes: 183e0907662SLois Curfman McInnes Using a modified Newton strategy, where the Jacobian remains fixed for several 184e0907662SLois Curfman McInnes iterations, can be cost effective in terms of overall nonlinear solution 185e0907662SLois Curfman McInnes efficiency. This parameter indicates that a new Jacobian will be computed every 186e0907662SLois Curfman McInnes <freq> nonlinear iterations. 187e0907662SLois Curfman McInnes 188b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency 189ef5ee4d1SLois Curfman McInnes 19040964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency(), MatFDColoringSetRecompute() 191005c665bSBarry Smith @*/ 192005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq) 193005c665bSBarry Smith { 1943a40ed3dSBarry Smith PetscFunctionBegin; 195005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 196005c665bSBarry Smith 197005c665bSBarry Smith matfd->freq = freq; 1983a40ed3dSBarry Smith PetscFunctionReturn(0); 199005c665bSBarry Smith } 200005c665bSBarry Smith 201005c665bSBarry Smith #undef __FUNC__ 202*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringGetFrequency" 203ff0cfa39SBarry Smith /*@ 204ff0cfa39SBarry Smith MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian 205ff0cfa39SBarry Smith matrices. 206ff0cfa39SBarry Smith 207ef5ee4d1SLois Curfman McInnes Not Collective 208ef5ee4d1SLois Curfman McInnes 209ff0cfa39SBarry Smith Input Parameters: 210ff0cfa39SBarry Smith . coloring - the coloring context 211ff0cfa39SBarry Smith 212ff0cfa39SBarry Smith Output Parameters: 213ff0cfa39SBarry Smith . freq - frequency (default is 1) 214ff0cfa39SBarry Smith 21515091d37SBarry Smith Options Database Keys: 21615091d37SBarry Smith . -mat_fd_coloring_freq <freq> - Sets coloring frequency 21715091d37SBarry Smith 21815091d37SBarry Smith Level: advanced 21915091d37SBarry Smith 220ff0cfa39SBarry Smith Notes: 221ff0cfa39SBarry Smith Using a modified Newton strategy, where the Jacobian remains fixed for several 222ff0cfa39SBarry Smith iterations, can be cost effective in terms of overall nonlinear solution 223ff0cfa39SBarry Smith efficiency. This parameter indicates that a new Jacobian will be computed every 224ff0cfa39SBarry Smith <freq> nonlinear iterations. 225ff0cfa39SBarry Smith 226ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency 227ef5ee4d1SLois Curfman McInnes 228ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency() 229ff0cfa39SBarry Smith @*/ 230ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq) 231ff0cfa39SBarry Smith { 2323a40ed3dSBarry Smith PetscFunctionBegin; 233ff0cfa39SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 234ff0cfa39SBarry Smith 235ff0cfa39SBarry Smith *freq = matfd->freq; 2363a40ed3dSBarry Smith PetscFunctionReturn(0); 237ff0cfa39SBarry Smith } 238ff0cfa39SBarry Smith 239ff0cfa39SBarry Smith #undef __FUNC__ 240*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringSetFunction" 241d64ed03dSBarry Smith /*@C 242005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 243005c665bSBarry Smith 244fee21e36SBarry Smith Collective on MatFDColoring 245fee21e36SBarry Smith 246ef5ee4d1SLois Curfman McInnes Input Parameters: 247ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 248ef5ee4d1SLois Curfman McInnes . f - the function 249ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 250ef5ee4d1SLois Curfman McInnes 25115091d37SBarry Smith Level: intermediate 25215091d37SBarry Smith 253f881d145SBarry Smith Notes: 254f881d145SBarry Smith In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to 255f881d145SBarry Smith be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 256f881d145SBarry Smith with the TS solvers. 257f881d145SBarry Smith 258b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 259005c665bSBarry Smith @*/ 260840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx) 261005c665bSBarry Smith { 2623a40ed3dSBarry Smith PetscFunctionBegin; 263005c665bSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 264005c665bSBarry Smith 265005c665bSBarry Smith matfd->f = f; 266005c665bSBarry Smith matfd->fctx = fctx; 267005c665bSBarry Smith 2683a40ed3dSBarry Smith PetscFunctionReturn(0); 269005c665bSBarry Smith } 270005c665bSBarry Smith 271005c665bSBarry Smith #undef __FUNC__ 272*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringSetFromOptions" 273639f9d9dSBarry Smith /*@ 274b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 275639f9d9dSBarry Smith the options database. 276639f9d9dSBarry Smith 277fee21e36SBarry Smith Collective on MatFDColoring 278fee21e36SBarry Smith 27965f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 280ef5ee4d1SLois Curfman McInnes .vb 28165f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 282f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 283f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 284ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 285ef5ee4d1SLois Curfman McInnes .ve 286ef5ee4d1SLois Curfman McInnes 287ef5ee4d1SLois Curfman McInnes Input Parameter: 288ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 289ef5ee4d1SLois Curfman McInnes 290b4fc646aSLois Curfman McInnes Options Database Keys: 291d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 292ef5ee4d1SLois Curfman McInnes of relative error in the function) 293f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 294ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian 295ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 296ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 297ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 298639f9d9dSBarry Smith 29915091d37SBarry Smith Level: intermediate 30015091d37SBarry Smith 301b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 302639f9d9dSBarry Smith @*/ 303639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd) 304639f9d9dSBarry Smith { 305186905e3SBarry Smith int ierr; 3063a40ed3dSBarry Smith 3073a40ed3dSBarry Smith PetscFunctionBegin; 308639f9d9dSBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE); 309639f9d9dSBarry Smith 310*b0a32e0cSBarry Smith ierr = PetscOptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 311*b0a32e0cSBarry Smith ierr = PetscOptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 312*b0a32e0cSBarry Smith ierr = PetscOptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 313*b0a32e0cSBarry Smith ierr = PetscOptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr); 314186905e3SBarry Smith /* not used here; but so they are presented in the GUI */ 315*b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 316*b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 317*b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 318*b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3193a40ed3dSBarry Smith PetscFunctionReturn(0); 320005c665bSBarry Smith } 321005c665bSBarry Smith 322433994e6SBarry Smith #undef __FUNC__ 323*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringView_Private" 324005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd) 325005c665bSBarry Smith { 326f1af5d2fSBarry Smith int ierr; 327f1af5d2fSBarry Smith PetscTruth flg; 328005c665bSBarry Smith 3293a40ed3dSBarry Smith PetscFunctionBegin; 330*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr); 331005c665bSBarry Smith if (flg) { 332*b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 333005c665bSBarry Smith } 334*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr); 335ae09f205SBarry Smith if (flg) { 336*b0a32e0cSBarry Smith ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(fd->comm),PETSC_VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr); 337*b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 338*b0a32e0cSBarry Smith ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr); 339ae09f205SBarry Smith } 340*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr); 341005c665bSBarry Smith if (flg) { 342*b0a32e0cSBarry Smith ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 343*b0a32e0cSBarry Smith ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(fd->comm));CHKERRQ(ierr); 344005c665bSBarry Smith } 3453a40ed3dSBarry Smith PetscFunctionReturn(0); 346bbf0e169SBarry Smith } 347bbf0e169SBarry Smith 3485615d1e5SSatish Balay #undef __FUNC__ 349*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringCreate" 350bbf0e169SBarry Smith /*@C 351639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 352639f9d9dSBarry Smith computation of Jacobians. 353bbf0e169SBarry Smith 354ef5ee4d1SLois Curfman McInnes Collective on Mat 355ef5ee4d1SLois Curfman McInnes 356639f9d9dSBarry Smith Input Parameters: 357ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 358ef5ee4d1SLois Curfman McInnes - iscoloring - the coloring of the matrix 359bbf0e169SBarry Smith 360bbf0e169SBarry Smith Output Parameter: 361639f9d9dSBarry Smith . color - the new coloring context 362bbf0e169SBarry Smith 363b4fc646aSLois Curfman McInnes Options Database Keys: 364ef5ee4d1SLois Curfman McInnes + -mat_fd_coloring_view - Activates basic viewing or coloring 365ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 366ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 367639f9d9dSBarry Smith 36815091d37SBarry Smith Level: intermediate 36915091d37SBarry Smith 3706831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 3716831982aSBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 372bbf0e169SBarry Smith @*/ 373639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 374bbf0e169SBarry Smith { 375639f9d9dSBarry Smith MatFDColoring c; 376639f9d9dSBarry Smith MPI_Comm comm; 377639f9d9dSBarry Smith int ierr,M,N; 378639f9d9dSBarry Smith 3793a40ed3dSBarry Smith PetscFunctionBegin; 380*b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 381639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 38229bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices"); 383639f9d9dSBarry Smith 384f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 3859194eea9SBarry Smith PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView); 386*b0a32e0cSBarry Smith PetscLogObjectCreate(c); 387639f9d9dSBarry Smith 388f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 389f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 390639f9d9dSBarry Smith } else { 39129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type"); 392639f9d9dSBarry Smith } 393639f9d9dSBarry Smith 394639f9d9dSBarry Smith c->error_rel = 1.e-8; 395ae09f205SBarry Smith c->umin = 1.e-6; 396005c665bSBarry Smith c->freq = 1; 39740964ad5SBarry Smith c->usersetsrecompute = PETSC_FALSE; 39840964ad5SBarry Smith c->recompute = PETSC_FALSE; 399005c665bSBarry Smith 400005c665bSBarry Smith ierr = MatFDColoringView_Private(c);CHKERRQ(ierr); 401639f9d9dSBarry Smith 402639f9d9dSBarry Smith *color = c; 403*b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4043a40ed3dSBarry Smith PetscFunctionReturn(0); 405bbf0e169SBarry Smith } 406bbf0e169SBarry Smith 4075615d1e5SSatish Balay #undef __FUNC__ 408*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringDestroy" 409bbf0e169SBarry Smith /*@C 410639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 411639f9d9dSBarry Smith via MatFDColoringCreate(). 412bbf0e169SBarry Smith 413ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 414ef5ee4d1SLois Curfman McInnes 415b4fc646aSLois Curfman McInnes Input Parameter: 416639f9d9dSBarry Smith . c - coloring context 417bbf0e169SBarry Smith 41815091d37SBarry Smith Level: intermediate 41915091d37SBarry Smith 420639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 421bbf0e169SBarry Smith @*/ 422639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c) 423bbf0e169SBarry Smith { 424263760aaSBarry Smith int i,ierr; 425bbf0e169SBarry Smith 4263a40ed3dSBarry Smith PetscFunctionBegin; 4273a40ed3dSBarry Smith if (--c->refct > 0) PetscFunctionReturn(0); 428d4bb536fSBarry Smith 429639f9d9dSBarry Smith for (i=0; i<c->ncolors; i++) { 430606d414cSSatish Balay if (c->columns[i]) {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);} 431606d414cSSatish Balay if (c->rows[i]) {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);} 432606d414cSSatish Balay if (c->columnsforrow[i]) {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);} 43330b34957SBarry Smith if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 434bbf0e169SBarry Smith } 435606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 436606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 437606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 438606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 439606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 44030b34957SBarry Smith if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);} 4414f113abfSBarry Smith if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 442005c665bSBarry Smith if (c->w1) { 443005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 444005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 445005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 446005c665bSBarry Smith } 447*b0a32e0cSBarry Smith PetscLogObjectDestroy(c); 448639f9d9dSBarry Smith PetscHeaderDestroy(c); 4493a40ed3dSBarry Smith PetscFunctionReturn(0); 450bbf0e169SBarry Smith } 45143a90d84SBarry Smith 4525615d1e5SSatish Balay #undef __FUNC__ 453*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringApply" 45443a90d84SBarry Smith /*@ 455e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 456e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 45743a90d84SBarry Smith 458fee21e36SBarry Smith Collective on MatFDColoring 459fee21e36SBarry Smith 460ef5ee4d1SLois Curfman McInnes Input Parameters: 461ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 462ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 463ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 464ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 465ef5ee4d1SLois Curfman McInnes 4668bba8e72SBarry Smith Options Database Keys: 467ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 4688bba8e72SBarry Smith 46915091d37SBarry Smith Level: intermediate 47015091d37SBarry Smith 47143a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 47243a90d84SBarry Smith 47343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 47443a90d84SBarry Smith @*/ 475005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 47643a90d84SBarry Smith { 4775904e1b1SBarry Smith int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f; 4785904e1b1SBarry Smith int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 4795904e1b1SBarry Smith Scalar dx,mone = -1.0,*y,*xx,*w3_array; 4804f113abfSBarry Smith Scalar *vscale_array; 481329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 482005c665bSBarry Smith Vec w1,w2,w3; 483005c665bSBarry Smith void *fctx = coloring->fctx; 484f1af5d2fSBarry Smith PetscTruth flg; 485005c665bSBarry Smith 486a2e34c3dSBarry Smith 4873a40ed3dSBarry Smith PetscFunctionBegin; 488e0907662SLois Curfman McInnes PetscValidHeaderSpecific(J,MAT_COOKIE); 489e0907662SLois Curfman McInnes PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 490e0907662SLois Curfman McInnes PetscValidHeaderSpecific(x1,VEC_COOKIE); 491e0907662SLois Curfman McInnes 49240964ad5SBarry Smith if (coloring->usersetsrecompute) { 49340964ad5SBarry Smith if (!coloring->recompute) { 49440964ad5SBarry Smith *flag = SAME_PRECONDITIONER; 49540964ad5SBarry Smith PetscFunctionReturn(0); 49640964ad5SBarry Smith } else { 49740964ad5SBarry Smith coloring->recompute = PETSC_FALSE; 49840964ad5SBarry Smith } 49940964ad5SBarry Smith } 50040964ad5SBarry Smith 501*b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 502005c665bSBarry Smith if (!coloring->w1) { 503005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 504*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w1); 505005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 506*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w2); 507005c665bSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 508*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w3); 509005c665bSBarry Smith } 510005c665bSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 51143a90d84SBarry Smith 5124c49b128SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 513*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 514f1af5d2fSBarry Smith if (flg) { 515*b0a32e0cSBarry Smith PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 516e0907662SLois Curfman McInnes } else { 51743a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 518e0907662SLois Curfman McInnes } 51943a90d84SBarry Smith 52043a90d84SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 52143a90d84SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 52243a90d84SBarry Smith ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 52343a90d84SBarry Smith 52443a90d84SBarry Smith /* 5259782cf98SBarry Smith Compute all the scale factors and share with other processors 5269782cf98SBarry Smith */ 5279782cf98SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 5284f113abfSBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 5299782cf98SBarry Smith for (k=0; k<coloring->ncolors; k++) { 5309782cf98SBarry Smith /* 5319782cf98SBarry Smith Loop over each column associated with color adding the 5329782cf98SBarry Smith perturbation to the vector w3. 5339782cf98SBarry Smith */ 5349782cf98SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 5359782cf98SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 5369782cf98SBarry Smith dx = xx[col]; 5379782cf98SBarry Smith if (dx == 0.0) dx = 1.0; 5389782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX) 5399782cf98SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 5409782cf98SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 5419782cf98SBarry Smith #else 5429782cf98SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 5439782cf98SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 5449782cf98SBarry Smith #endif 5459782cf98SBarry Smith dx *= epsilon; 54630b34957SBarry Smith vscale_array[col] = 1.0/dx; 5479782cf98SBarry Smith } 5489782cf98SBarry Smith } 549a2e34c3dSBarry Smith vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 55030b34957SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 55130b34957SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5529782cf98SBarry Smith 553*b0a32e0cSBarry Smith ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); 554*b0a32e0cSBarry Smith ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD); 555*b0a32e0cSBarry Smith 55630b34957SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 55730b34957SBarry Smith else vscaleforrow = coloring->columnsforrow; 55830b34957SBarry Smith 55930b34957SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 5609782cf98SBarry Smith /* 56143a90d84SBarry Smith Loop over each color 56243a90d84SBarry Smith */ 56343a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 56443a90d84SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 56542460c72SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 56643a90d84SBarry Smith /* 56743a90d84SBarry Smith Loop over each column associated with color adding the 56843a90d84SBarry Smith perturbation to the vector w3. 56943a90d84SBarry Smith */ 57043a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 57143a90d84SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 57242460c72SBarry Smith dx = xx[col]; 573ae09f205SBarry Smith if (dx == 0.0) dx = 1.0; 574aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 575ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 576ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 57743a90d84SBarry Smith #else 578329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 579329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 58043a90d84SBarry Smith #endif 58143a90d84SBarry Smith dx *= epsilon; 58229bbc08cSBarry Smith if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter"); 58342460c72SBarry Smith w3_array[col] += dx; 58443a90d84SBarry Smith } 58542460c72SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 5863b28642cSBarry Smith 58743a90d84SBarry Smith /* 588e0907662SLois Curfman McInnes Evaluate function at x1 + dx (here dx is a vector of perturbations) 58943a90d84SBarry Smith */ 590a2e34c3dSBarry Smith 59143a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 59243a90d84SBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 5939782cf98SBarry Smith 59443a90d84SBarry Smith /* 595e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 59643a90d84SBarry Smith */ 5973b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 59843a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 59943a90d84SBarry Smith row = coloring->rows[k][l]; 60043a90d84SBarry Smith col = coloring->columnsforrow[k][l]; 6015904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 60243a90d84SBarry Smith srow = row + start; 60343a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 60443a90d84SBarry Smith } 6053b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 60643a90d84SBarry Smith } 60730b34957SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 60842460c72SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 60943a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61043a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 611*b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 612a2e34c3dSBarry Smith 613*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr); 614a2e34c3dSBarry Smith if (flg) { 615a2e34c3dSBarry Smith ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr); 616a2e34c3dSBarry Smith } 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 61843a90d84SBarry Smith } 619840b8ebdSBarry Smith 620840b8ebdSBarry Smith #undef __FUNC__ 621*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringApplyTS" 622840b8ebdSBarry Smith /*@ 623840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 624840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 625840b8ebdSBarry Smith 626fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 627fee21e36SBarry Smith 628ef5ee4d1SLois Curfman McInnes Input Parameters: 6293b28642cSBarry Smith + mat - location to store Jacobian 630ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 631ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 632ef5ee4d1SLois Curfman McInnes - sctx - optional context required by function (actually a SNES context) 633ef5ee4d1SLois Curfman McInnes 634840b8ebdSBarry Smith Options Database Keys: 635ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Sets coloring frequency 636840b8ebdSBarry Smith 63715091d37SBarry Smith Level: intermediate 63815091d37SBarry Smith 639840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView() 640840b8ebdSBarry Smith 641840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 642840b8ebdSBarry Smith @*/ 643329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 644840b8ebdSBarry Smith { 645329f5518SBarry Smith int (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f; 6465904e1b1SBarry Smith int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow; 6475904e1b1SBarry Smith Scalar dx,mone = -1.0,*y,*xx,*w3_array; 6485904e1b1SBarry Smith Scalar *vscale_array; 649329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 650840b8ebdSBarry Smith Vec w1,w2,w3; 651840b8ebdSBarry Smith void *fctx = coloring->fctx; 652f1af5d2fSBarry Smith PetscTruth flg; 653840b8ebdSBarry Smith 6543a40ed3dSBarry Smith PetscFunctionBegin; 655840b8ebdSBarry Smith PetscValidHeaderSpecific(J,MAT_COOKIE); 656840b8ebdSBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE); 657840b8ebdSBarry Smith PetscValidHeaderSpecific(x1,VEC_COOKIE); 658840b8ebdSBarry Smith 659*b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 660840b8ebdSBarry Smith if (!coloring->w1) { 661840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); 662*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w1); 663840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); 664*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w2); 665840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 666*b0a32e0cSBarry Smith PetscLogObjectParent(coloring,coloring->w3); 667840b8ebdSBarry Smith } 668840b8ebdSBarry Smith w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; 669840b8ebdSBarry Smith 6705904e1b1SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 671*b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); 672f1af5d2fSBarry Smith if (flg) { 673*b0a32e0cSBarry Smith PetscLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n"); 674840b8ebdSBarry Smith } else { 675840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 676840b8ebdSBarry Smith } 677840b8ebdSBarry Smith 678840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 679840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 680840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 681840b8ebdSBarry Smith 682840b8ebdSBarry Smith /* 6835904e1b1SBarry Smith Compute all the scale factors and share with other processors 684840b8ebdSBarry Smith */ 6855904e1b1SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 6865904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 687840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 688840b8ebdSBarry Smith /* 689840b8ebdSBarry Smith Loop over each column associated with color adding the 690840b8ebdSBarry Smith perturbation to the vector w3. 691840b8ebdSBarry Smith */ 692840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 693840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 6945904e1b1SBarry Smith dx = xx[col]; 695840b8ebdSBarry Smith if (dx == 0.0) dx = 1.0; 696aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 697840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 698840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 699840b8ebdSBarry Smith #else 700329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 701329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 702840b8ebdSBarry Smith #endif 703840b8ebdSBarry Smith dx *= epsilon; 7045904e1b1SBarry Smith vscale_array[col] = 1.0/dx; 705840b8ebdSBarry Smith } 7065904e1b1SBarry Smith } 7075904e1b1SBarry Smith vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7085904e1b1SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7095904e1b1SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7105904e1b1SBarry Smith 7115904e1b1SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 7125904e1b1SBarry Smith else vscaleforrow = coloring->columnsforrow; 7135904e1b1SBarry Smith 7145904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7155904e1b1SBarry Smith /* 7165904e1b1SBarry Smith Loop over each color 7175904e1b1SBarry Smith */ 7185904e1b1SBarry Smith for (k=0; k<coloring->ncolors; k++) { 7195904e1b1SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 7205904e1b1SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 7215904e1b1SBarry Smith /* 7225904e1b1SBarry Smith Loop over each column associated with color adding the 7235904e1b1SBarry Smith perturbation to the vector w3. 7245904e1b1SBarry Smith */ 7255904e1b1SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 7265904e1b1SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 7275904e1b1SBarry Smith dx = xx[col]; 7285904e1b1SBarry Smith if (dx == 0.0) dx = 1.0; 7295904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 7305904e1b1SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 7315904e1b1SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 7325904e1b1SBarry Smith #else 7335904e1b1SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 7345904e1b1SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 7355904e1b1SBarry Smith #endif 7365904e1b1SBarry Smith dx *= epsilon; 7375904e1b1SBarry Smith w3_array[col] += dx; 7385904e1b1SBarry Smith } 7395904e1b1SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 7405904e1b1SBarry Smith 741840b8ebdSBarry Smith /* 742840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 743840b8ebdSBarry Smith */ 744840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 745840b8ebdSBarry Smith ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr); 7465904e1b1SBarry Smith 747840b8ebdSBarry Smith /* 748840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 749840b8ebdSBarry Smith */ 7503b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 751840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 752840b8ebdSBarry Smith row = coloring->rows[k][l]; 753840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 7545904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 755840b8ebdSBarry Smith srow = row + start; 756840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 757840b8ebdSBarry Smith } 7583b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 759840b8ebdSBarry Smith } 7605904e1b1SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7615904e1b1SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 762840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 763840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 764*b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 7653a40ed3dSBarry Smith PetscFunctionReturn(0); 766840b8ebdSBarry Smith } 7673b28642cSBarry Smith 7683b28642cSBarry Smith 76940964ad5SBarry Smith #undef __FUNC__ 770*b0a32e0cSBarry Smith #define __FUNC__ "MatFDColoringSetRecompute()" 77140964ad5SBarry Smith /*@C 77240964ad5SBarry Smith MatFDColoringSetRecompute - Indicates that the next time a Jacobian preconditioner 77340964ad5SBarry Smith is needed it sholuld be recomputed. Once this is called and the new Jacobian is computed 77440964ad5SBarry Smith no additional Jacobian's will be computed (the same one will be used) until this is 77540964ad5SBarry Smith called again. 77640964ad5SBarry Smith 77740964ad5SBarry Smith Collective on MatFDColoring 77840964ad5SBarry Smith 77940964ad5SBarry Smith Input Parameters: 78040964ad5SBarry Smith . c - the coloring context 78140964ad5SBarry Smith 78240964ad5SBarry Smith Level: intermediate 78340964ad5SBarry Smith 78440964ad5SBarry Smith Notes: The MatFDColoringSetFrequency() is ignored once this is called 78540964ad5SBarry Smith 78640964ad5SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFrequency() 78740964ad5SBarry Smith 78840964ad5SBarry Smith .keywords: Mat, finite differences, coloring 78940964ad5SBarry Smith @*/ 79040964ad5SBarry Smith int MatFDColoringSetRecompute(MatFDColoring c) 79140964ad5SBarry Smith { 79240964ad5SBarry Smith PetscFunctionBegin; 79340964ad5SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE); 79440964ad5SBarry Smith c->usersetsrecompute = PETSC_TRUE; 79540964ad5SBarry Smith c->recompute = PETSC_TRUE; 79640964ad5SBarry Smith PetscFunctionReturn(0); 79740964ad5SBarry Smith } 79840964ad5SBarry Smith 7993b28642cSBarry Smith 800