1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3bbf0e169SBarry Smith /* 4639f9d9dSBarry Smith This is where the abstract matrix operations are defined that are 5639f9d9dSBarry Smith used for finite difference computations of Jacobians using coloring. 6bbf0e169SBarry Smith */ 7bbf0e169SBarry Smith 87c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 9bbf0e169SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 113a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF" 12be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring fd,Vec F) 133a7fca6bSBarry Smith { 143a7fca6bSBarry Smith PetscFunctionBegin; 153a7fca6bSBarry Smith fd->F = F; 163a7fca6bSBarry Smith PetscFunctionReturn(0); 173a7fca6bSBarry Smith } 183a7fca6bSBarry Smith 193a7fca6bSBarry Smith #undef __FUNCT__ 204a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom" 216849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa) 229194eea9SBarry Smith { 239194eea9SBarry Smith MatFDColoring fd = (MatFDColoring)Aa; 24dfbe8321SBarry Smith PetscErrorCode ierr; 2538baddfdSBarry Smith PetscInt i,j; 269194eea9SBarry Smith PetscReal x,y; 279194eea9SBarry Smith 289194eea9SBarry Smith PetscFunctionBegin; 299194eea9SBarry Smith 309194eea9SBarry Smith /* loop over colors */ 319194eea9SBarry Smith for (i=0; i<fd->ncolors; i++) { 329194eea9SBarry Smith for (j=0; j<fd->nrows[i]; j++) { 339194eea9SBarry Smith y = fd->M - fd->rows[i][j] - fd->rstart; 349194eea9SBarry Smith x = fd->columnsforrow[i][j]; 35b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr); 369194eea9SBarry Smith } 379194eea9SBarry Smith } 389194eea9SBarry Smith PetscFunctionReturn(0); 399194eea9SBarry Smith } 409194eea9SBarry Smith 414a2ae208SSatish Balay #undef __FUNCT__ 424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw" 436849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer) 44005c665bSBarry Smith { 45dfbe8321SBarry Smith PetscErrorCode ierr; 46005c665bSBarry Smith PetscTruth isnull; 47b0a32e0cSBarry Smith PetscDraw draw; 4854d96fa3SBarry Smith PetscReal xr,yr,xl,yl,h,w; 49005c665bSBarry Smith 503a40ed3dSBarry Smith PetscFunctionBegin; 51b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 52b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 539194eea9SBarry Smith 549194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 55005c665bSBarry Smith 56005c665bSBarry Smith xr = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0; 57005c665bSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 58b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 59b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr); 609194eea9SBarry Smith ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 613a40ed3dSBarry Smith PetscFunctionReturn(0); 62005c665bSBarry Smith } 63005c665bSBarry Smith 644a2ae208SSatish Balay #undef __FUNCT__ 654a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView" 66bbf0e169SBarry Smith /*@C 67639f9d9dSBarry Smith MatFDColoringView - Views a finite difference coloring context. 68bbf0e169SBarry Smith 694c49b128SBarry Smith Collective on MatFDColoring 70fee21e36SBarry Smith 71ef5ee4d1SLois Curfman McInnes Input Parameters: 72ef5ee4d1SLois Curfman McInnes + c - the coloring context 73ef5ee4d1SLois Curfman McInnes - viewer - visualization context 74ef5ee4d1SLois Curfman McInnes 7515091d37SBarry Smith Level: intermediate 7615091d37SBarry Smith 77b4fc646aSLois Curfman McInnes Notes: 78b4fc646aSLois Curfman McInnes The available visualization contexts include 79b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 80b0a32e0cSBarry Smith . PETSC_VIEWER_STDOUT_WORLD - synchronized standard 81ef5ee4d1SLois Curfman McInnes output where only the first processor opens 82ef5ee4d1SLois Curfman McInnes the file. All other processors send their 83ef5ee4d1SLois Curfman McInnes data to the first processor to print. 84b0a32e0cSBarry Smith - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure 85bbf0e169SBarry Smith 869194eea9SBarry Smith Notes: 879194eea9SBarry Smith Since PETSc uses only a small number of basic colors (currently 33), if the coloring 88b0a32e0cSBarry Smith involves more than 33 then some seemingly identical colors are displayed making it look 899194eea9SBarry Smith like an illegal coloring. This is just a graphical artifact. 909194eea9SBarry Smith 91639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 92005c665bSBarry Smith 93b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view 94bbf0e169SBarry Smith @*/ 95be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring c,PetscViewer viewer) 96bbf0e169SBarry Smith { 976849ba73SBarry Smith PetscErrorCode ierr; 9838baddfdSBarry Smith PetscInt i,j; 9932077d6dSBarry Smith PetscTruth isdraw,iascii; 100f3ef73ceSBarry Smith PetscViewerFormat format; 101bbf0e169SBarry Smith 1023a40ed3dSBarry Smith PetscFunctionBegin; 1030700a824SBarry Smith PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1); 1043050cee2SBarry Smith if (!viewer) { 1057adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr); 1063050cee2SBarry Smith } 1070700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 108c9780b6fSBarry Smith PetscCheckSameComm(c,1,viewer,2); 109bbf0e169SBarry Smith 1102692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1112692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1120f5bd95cSBarry Smith if (isdraw) { 113b4fc646aSLois Curfman McInnes ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr); 11432077d6dSBarry Smith } else if (iascii) { 115b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr); 116a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr); 117a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Umin=%G\n",c->umin);CHKERRQ(ierr); 11877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of colors=%D\n",c->ncolors);CHKERRQ(ierr); 119ae09f205SBarry Smith 120b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 121fb9695e5SSatish Balay if (format != PETSC_VIEWER_ASCII_INFO) { 122b4fc646aSLois Curfman McInnes for (i=0; i<c->ncolors; i++) { 12377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Information for color %D\n",i);CHKERRQ(ierr); 12477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr); 125b4fc646aSLois Curfman McInnes for (j=0; j<c->ncolumns[i]; j++) { 12677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D\n",c->columns[i][j]);CHKERRQ(ierr); 127639f9d9dSBarry Smith } 12877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr); 129b4fc646aSLois Curfman McInnes for (j=0; j<c->nrows[i]; j++) { 13077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr); 131b4fc646aSLois Curfman McInnes } 132bbf0e169SBarry Smith } 133bbf0e169SBarry Smith } 134b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1355cd90555SBarry Smith } else { 136e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name); 137bbf0e169SBarry Smith } 1383a40ed3dSBarry Smith PetscFunctionReturn(0); 139639f9d9dSBarry Smith } 140639f9d9dSBarry Smith 1414a2ae208SSatish Balay #undef __FUNCT__ 1424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters" 143639f9d9dSBarry Smith /*@ 144b4fc646aSLois Curfman McInnes MatFDColoringSetParameters - Sets the parameters for the sparse approximation of 145b4fc646aSLois Curfman McInnes a Jacobian matrix using finite differences. 146639f9d9dSBarry Smith 147ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 148ef5ee4d1SLois Curfman McInnes 149ef5ee4d1SLois Curfman McInnes The Jacobian is estimated with the differencing approximation 150ef5ee4d1SLois Curfman McInnes .vb 15165f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 152f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 153f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 154ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 155ef5ee4d1SLois Curfman McInnes .ve 156639f9d9dSBarry Smith 157639f9d9dSBarry Smith Input Parameters: 158ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 159639f9d9dSBarry Smith . error_rel - relative error 160f23b5b22SLois Curfman McInnes - umin - minimum allowable u-value magnitude 161fee21e36SBarry Smith 16215091d37SBarry Smith Level: advanced 16315091d37SBarry Smith 164b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters 165b4fc646aSLois Curfman McInnes 16695b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions() 16795b89fc3SBarry Smith 168639f9d9dSBarry Smith @*/ 169be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin) 170639f9d9dSBarry Smith { 1713a40ed3dSBarry Smith PetscFunctionBegin; 1720700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 173639f9d9dSBarry Smith 174639f9d9dSBarry Smith if (error != PETSC_DEFAULT) matfd->error_rel = error; 175639f9d9dSBarry Smith if (umin != PETSC_DEFAULT) matfd->umin = umin; 1763a40ed3dSBarry Smith PetscFunctionReturn(0); 177639f9d9dSBarry Smith } 178639f9d9dSBarry Smith 179005c665bSBarry Smith 180ff0cfa39SBarry Smith 1814a2ae208SSatish Balay #undef __FUNCT__ 1824a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction" 1834a9d489dSBarry Smith /*@C 1844a9d489dSBarry Smith MatFDColoringGetFunction - Gets the function to use for computing the Jacobian. 1854a9d489dSBarry Smith 1864a9d489dSBarry Smith Collective on MatFDColoring 1874a9d489dSBarry Smith 1884a9d489dSBarry Smith Input Parameters: 1894a9d489dSBarry Smith . coloring - the coloring context 1904a9d489dSBarry Smith 1914a9d489dSBarry Smith Output Parameters: 1924a9d489dSBarry Smith + f - the function 1934a9d489dSBarry Smith - fctx - the optional user-defined function context 1944a9d489dSBarry Smith 1954a9d489dSBarry Smith Level: intermediate 1964a9d489dSBarry Smith 1974a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function 19895b89fc3SBarry Smith 19995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions() 20095b89fc3SBarry Smith 2014a9d489dSBarry Smith @*/ 2024a9d489dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx) 2034a9d489dSBarry Smith { 2044a9d489dSBarry Smith PetscFunctionBegin; 2050700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 2064a9d489dSBarry Smith if (f) *f = matfd->f; 2074a9d489dSBarry Smith if (fctx) *fctx = matfd->fctx; 2084a9d489dSBarry Smith PetscFunctionReturn(0); 2094a9d489dSBarry Smith } 2104a9d489dSBarry Smith 2114a9d489dSBarry Smith #undef __FUNCT__ 2124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction" 213d64ed03dSBarry Smith /*@C 214005c665bSBarry Smith MatFDColoringSetFunction - Sets the function to use for computing the Jacobian. 215005c665bSBarry Smith 216fee21e36SBarry Smith Collective on MatFDColoring 217fee21e36SBarry Smith 218ef5ee4d1SLois Curfman McInnes Input Parameters: 219ef5ee4d1SLois Curfman McInnes + coloring - the coloring context 220ef5ee4d1SLois Curfman McInnes . f - the function 221ef5ee4d1SLois Curfman McInnes - fctx - the optional user-defined function context 222ef5ee4d1SLois Curfman McInnes 2237850c7c0SBarry Smith Calling sequence of (*f) function: 2247850c7c0SBarry Smith For SNES: PetscErrorCode (*f)(SNES,Vec,Vec,void*) 2257850c7c0SBarry Smith For TS: PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*) 2267850c7c0SBarry Smith If not using SNES or TS: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored 22715091d37SBarry Smith 2287850c7c0SBarry Smith Level: advanced 2297850c7c0SBarry Smith 2307850c7c0SBarry Smith Notes: This function is usually used automatically by SNES or TS (when one uses SNESSetJacobian() with the argument 2317850c7c0SBarry Smith SNESDefaultComputeJacobianColor() or TSSetRHSJacobian() with the argument TSDefaultComputeJacobianColor()) and only needs to be used 2327850c7c0SBarry Smith by someone computing a matrix via coloring directly by calling MatFDColoringApply() 2337850c7c0SBarry Smith 2347850c7c0SBarry Smith Fortran Notes: 2357850c7c0SBarry Smith In Fortran you must call MatFDColoringSetFunction() for a coloring object to 2367850c7c0SBarry Smith be used without SNES or TS or within the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used 2377850c7c0SBarry Smith within the TS solvers. 238f881d145SBarry Smith 239b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function 24095b89fc3SBarry Smith 24195b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions() 24295b89fc3SBarry Smith 243005c665bSBarry Smith @*/ 244be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx) 245005c665bSBarry Smith { 2463a40ed3dSBarry Smith PetscFunctionBegin; 2470700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 248005c665bSBarry Smith matfd->f = f; 249005c665bSBarry Smith matfd->fctx = fctx; 2503a40ed3dSBarry Smith PetscFunctionReturn(0); 251005c665bSBarry Smith } 252005c665bSBarry Smith 2534a2ae208SSatish Balay #undef __FUNCT__ 2544a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions" 255639f9d9dSBarry Smith /*@ 256b4fc646aSLois Curfman McInnes MatFDColoringSetFromOptions - Sets coloring finite difference parameters from 257639f9d9dSBarry Smith the options database. 258639f9d9dSBarry Smith 259fee21e36SBarry Smith Collective on MatFDColoring 260fee21e36SBarry Smith 26165f2ba5bSLois Curfman McInnes The Jacobian, F'(u), is estimated with the differencing approximation 262ef5ee4d1SLois Curfman McInnes .vb 26365f2ba5bSLois Curfman McInnes F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where 264f23b5b22SLois Curfman McInnes h = error_rel*u[i] if abs(u[i]) > umin 265f23b5b22SLois Curfman McInnes = +/- error_rel*umin otherwise, with +/- determined by the sign of u[i] 266ef5ee4d1SLois Curfman McInnes dx_{i} = (0, ... 1, .... 0) 267ef5ee4d1SLois Curfman McInnes .ve 268ef5ee4d1SLois Curfman McInnes 269ef5ee4d1SLois Curfman McInnes Input Parameter: 270ef5ee4d1SLois Curfman McInnes . coloring - the coloring context 271ef5ee4d1SLois Curfman McInnes 272b4fc646aSLois Curfman McInnes Options Database Keys: 273d15ffeeaSSatish Balay + -mat_fd_coloring_err <err> - Sets <err> (square root 274ef5ee4d1SLois Curfman McInnes of relative error in the function) 275f23b5b22SLois Curfman McInnes . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 2763ec795f1SBarry Smith . -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 277ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing 278ef5ee4d1SLois Curfman McInnes . -mat_fd_coloring_view_info - Activates viewing info 279ef5ee4d1SLois Curfman McInnes - -mat_fd_coloring_view_draw - Activates drawing 280639f9d9dSBarry Smith 28115091d37SBarry Smith Level: intermediate 28215091d37SBarry Smith 283b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters 284d1c39f55SBarry Smith 285d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters() 286d1c39f55SBarry Smith 287639f9d9dSBarry Smith @*/ 288be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring matfd) 289639f9d9dSBarry Smith { 290dfbe8321SBarry Smith PetscErrorCode ierr; 291efb30889SBarry Smith PetscTruth flg; 292efb30889SBarry Smith char value[3]; 2933a40ed3dSBarry Smith 2943a40ed3dSBarry Smith PetscFunctionBegin; 2950700a824SBarry Smith PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1); 296639f9d9dSBarry Smith 2977adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)matfd)->comm,((PetscObject)matfd)->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr); 29887828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr); 29987828ca2SBarry Smith ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr); 300*1bc75a8dSBarry Smith ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr); 301efb30889SBarry Smith if (flg) { 302efb30889SBarry Smith if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp"; 303efb30889SBarry Smith else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds"; 304e32f2f54SBarry Smith else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value); 305efb30889SBarry Smith } 306186905e3SBarry Smith /* not used here; but so they are presented in the GUI */ 307b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr); 308b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr); 309b0a32e0cSBarry Smith ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr); 310b0a32e0cSBarry Smith PetscOptionsEnd();CHKERRQ(ierr); 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 312005c665bSBarry Smith } 313005c665bSBarry Smith 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private" 316dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd) 317005c665bSBarry Smith { 318dfbe8321SBarry Smith PetscErrorCode ierr; 31990d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE; 3203050cee2SBarry Smith PetscViewer viewer; 321005c665bSBarry Smith 3223a40ed3dSBarry Smith PetscFunctionBegin; 3237adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr); 32490d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr); 325005c665bSBarry Smith if (flg) { 3263050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 327005c665bSBarry Smith } 32890d69ab7SBarry Smith flg = PETSC_FALSE; 32990d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr); 330ae09f205SBarry Smith if (flg) { 3313050cee2SBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 3323050cee2SBarry Smith ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr); 3333050cee2SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 334ae09f205SBarry Smith } 33590d69ab7SBarry Smith flg = PETSC_FALSE; 33690d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr); 337005c665bSBarry Smith if (flg) { 3387adad957SLisandro Dalcin ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 3397adad957SLisandro Dalcin ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr); 340005c665bSBarry Smith } 3413a40ed3dSBarry Smith PetscFunctionReturn(0); 342bbf0e169SBarry Smith } 343bbf0e169SBarry Smith 3444a2ae208SSatish Balay #undef __FUNCT__ 3454a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate" 34605869f15SSatish Balay /*@ 347639f9d9dSBarry Smith MatFDColoringCreate - Creates a matrix coloring context for finite difference 348639f9d9dSBarry Smith computation of Jacobians. 349bbf0e169SBarry Smith 350ef5ee4d1SLois Curfman McInnes Collective on Mat 351ef5ee4d1SLois Curfman McInnes 352639f9d9dSBarry Smith Input Parameters: 353ef5ee4d1SLois Curfman McInnes + mat - the matrix containing the nonzero structure of the Jacobian 3541d0fab5eSBarry Smith - iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DAGetColoring() 355bbf0e169SBarry Smith 356bbf0e169SBarry Smith Output Parameter: 357639f9d9dSBarry Smith . color - the new coloring context 358bbf0e169SBarry Smith 35915091d37SBarry Smith Level: intermediate 36015091d37SBarry Smith 3616831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(), 362d1c39f55SBarry Smith MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(), 3631d0fab5eSBarry Smith MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DAGetColoring() 364bbf0e169SBarry Smith @*/ 365be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color) 366bbf0e169SBarry Smith { 367639f9d9dSBarry Smith MatFDColoring c; 368639f9d9dSBarry Smith MPI_Comm comm; 369dfbe8321SBarry Smith PetscErrorCode ierr; 37038baddfdSBarry Smith PetscInt M,N; 371b8f8c88eSHong Zhang PetscMPIInt size; 372639f9d9dSBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 374d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 375639f9d9dSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 37617186662SBarry Smith if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices"); 377639f9d9dSBarry Smith 378f881d145SBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 3790700a824SBarry Smith ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr); 380639f9d9dSBarry Smith 381b8f8c88eSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 382b8f8c88eSHong Zhang c->ctype = iscoloring->ctype; 383b8f8c88eSHong Zhang 384f830108cSBarry Smith if (mat->ops->fdcoloringcreate) { 385f830108cSBarry Smith ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr); 38617186662SBarry Smith } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for this matrix type"); 387639f9d9dSBarry Smith 388b8f8c88eSHong Zhang ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr); 389b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr); 390b8f8c88eSHong Zhang ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr); 391b8f8c88eSHong Zhang ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr); 392b8f8c88eSHong Zhang 39377d8c4bbSBarry Smith c->error_rel = PETSC_SQRT_MACHINE_EPSILON; 39477d8c4bbSBarry Smith c->umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 39549b058dcSBarry Smith c->currentcolor = -1; 396efb30889SBarry Smith c->htype = "wp"; 397639f9d9dSBarry Smith 398639f9d9dSBarry Smith *color = c; 399d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr); 4003a40ed3dSBarry Smith PetscFunctionReturn(0); 401bbf0e169SBarry Smith } 402bbf0e169SBarry Smith 4034a2ae208SSatish Balay #undef __FUNCT__ 4044a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy" 40505869f15SSatish Balay /*@ 406639f9d9dSBarry Smith MatFDColoringDestroy - Destroys a matrix coloring context that was created 407639f9d9dSBarry Smith via MatFDColoringCreate(). 408bbf0e169SBarry Smith 409ef5ee4d1SLois Curfman McInnes Collective on MatFDColoring 410ef5ee4d1SLois Curfman McInnes 411b4fc646aSLois Curfman McInnes Input Parameter: 412639f9d9dSBarry Smith . c - coloring context 413bbf0e169SBarry Smith 41415091d37SBarry Smith Level: intermediate 41515091d37SBarry Smith 416639f9d9dSBarry Smith .seealso: MatFDColoringCreate() 417bbf0e169SBarry Smith @*/ 418be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring c) 419bbf0e169SBarry Smith { 4206849ba73SBarry Smith PetscErrorCode ierr; 42138baddfdSBarry Smith PetscInt i; 422bbf0e169SBarry Smith 4233a40ed3dSBarry Smith PetscFunctionBegin; 4247adad957SLisandro Dalcin if (--((PetscObject)c)->refct > 0) PetscFunctionReturn(0); 425d4bb536fSBarry Smith 426639f9d9dSBarry Smith for (i=0; i<c->ncolors; i++) { 42705b42c5fSBarry Smith ierr = PetscFree(c->columns[i]);CHKERRQ(ierr); 42805b42c5fSBarry Smith ierr = PetscFree(c->rows[i]);CHKERRQ(ierr); 42905b42c5fSBarry Smith ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr); 43005b42c5fSBarry Smith if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);} 431bbf0e169SBarry Smith } 432606d414cSSatish Balay ierr = PetscFree(c->ncolumns);CHKERRQ(ierr); 433606d414cSSatish Balay ierr = PetscFree(c->columns);CHKERRQ(ierr); 434606d414cSSatish Balay ierr = PetscFree(c->nrows);CHKERRQ(ierr); 435606d414cSSatish Balay ierr = PetscFree(c->rows);CHKERRQ(ierr); 436606d414cSSatish Balay ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr); 43705b42c5fSBarry Smith ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr); 4384f113abfSBarry Smith if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);} 439005c665bSBarry Smith if (c->w1) { 440005c665bSBarry Smith ierr = VecDestroy(c->w1);CHKERRQ(ierr); 441005c665bSBarry Smith ierr = VecDestroy(c->w2);CHKERRQ(ierr); 442b8f8c88eSHong Zhang } 443b8f8c88eSHong Zhang if (c->w3){ 444005c665bSBarry Smith ierr = VecDestroy(c->w3);CHKERRQ(ierr); 445005c665bSBarry Smith } 446d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 4473a40ed3dSBarry Smith PetscFunctionReturn(0); 448bbf0e169SBarry Smith } 44943a90d84SBarry Smith 4504a2ae208SSatish Balay #undef __FUNCT__ 45149b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns" 45249b058dcSBarry Smith /*@C 45349b058dcSBarry Smith MatFDColoringGetPerturbedColumns - Returns the indices of the columns that 45449b058dcSBarry Smith that are currently being perturbed. 45549b058dcSBarry Smith 45649b058dcSBarry Smith Not Collective 45749b058dcSBarry Smith 45849b058dcSBarry Smith Input Parameters: 45949b058dcSBarry Smith . coloring - coloring context created with MatFDColoringCreate() 46049b058dcSBarry Smith 46149b058dcSBarry Smith Output Parameters: 46249b058dcSBarry Smith + n - the number of local columns being perturbed 46349b058dcSBarry Smith - cols - the column indices, in global numbering 46449b058dcSBarry Smith 46549b058dcSBarry Smith Level: intermediate 46649b058dcSBarry Smith 46749b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply() 46849b058dcSBarry Smith 46949b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences 47049b058dcSBarry Smith @*/ 471be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[]) 47249b058dcSBarry Smith { 47349b058dcSBarry Smith PetscFunctionBegin; 47449b058dcSBarry Smith if (coloring->currentcolor >= 0) { 47549b058dcSBarry Smith *n = coloring->ncolumns[coloring->currentcolor]; 47649b058dcSBarry Smith *cols = coloring->columns[coloring->currentcolor]; 47749b058dcSBarry Smith } else { 47849b058dcSBarry Smith *n = 0; 47949b058dcSBarry Smith } 48049b058dcSBarry Smith PetscFunctionReturn(0); 48149b058dcSBarry Smith } 48249b058dcSBarry Smith 48349b058dcSBarry Smith #undef __FUNCT__ 4844a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply" 48543a90d84SBarry Smith /*@ 486e0907662SLois Curfman McInnes MatFDColoringApply - Given a matrix for which a MatFDColoring context 487e0907662SLois Curfman McInnes has been created, computes the Jacobian for a function via finite differences. 48843a90d84SBarry Smith 489fee21e36SBarry Smith Collective on MatFDColoring 490fee21e36SBarry Smith 491ef5ee4d1SLois Curfman McInnes Input Parameters: 492ef5ee4d1SLois Curfman McInnes + mat - location to store Jacobian 493ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 494ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 4957850c7c0SBarry Smith - sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null 496ef5ee4d1SLois Curfman McInnes 4978bba8e72SBarry Smith Options Database Keys: 498ebd3b9afSBarry Smith + -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS) 499b9722508SLois Curfman McInnes . -mat_fd_coloring_view - Activates basic viewing or coloring 500b9722508SLois Curfman McInnes . -mat_fd_coloring_view_draw - Activates drawing of coloring 501b9722508SLois Curfman McInnes - -mat_fd_coloring_view_info - Activates viewing of coloring info 5028bba8e72SBarry Smith 50315091d37SBarry Smith Level: intermediate 50415091d37SBarry Smith 5057850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 50643a90d84SBarry Smith 50743a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences 50843a90d84SBarry Smith @*/ 509be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 51043a90d84SBarry Smith { 5113acb8795SBarry Smith PetscErrorCode ierr; 5123acb8795SBarry Smith 5133acb8795SBarry Smith PetscFunctionBegin; 5140700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 5150700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 5160700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 51717186662SBarry Smith if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 518e32f2f54SBarry Smith if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name); 5193acb8795SBarry Smith ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr); 5203acb8795SBarry Smith PetscFunctionReturn(0); 5213acb8795SBarry Smith } 5223acb8795SBarry Smith 5233acb8795SBarry Smith #undef __FUNCT__ 5243acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ" 5253acb8795SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 5263acb8795SBarry Smith { 5276849ba73SBarry Smith PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; 5286849ba73SBarry Smith PetscErrorCode ierr; 529b8f8c88eSHong Zhang PetscInt k,start,end,l,row,col,srow,**vscaleforrow,m1,m2; 530efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 531ea709b57SSatish Balay PetscScalar *vscale_array; 532efb30889SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 533b8f8c88eSHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 534005c665bSBarry Smith void *fctx = coloring->fctx; 53590d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE; 536705d48f7SSatish Balay PetscInt ctype=coloring->ctype,N,col_start=0,col_end=0; 537b8f8c88eSHong Zhang Vec x1_tmp; 538a2e34c3dSBarry Smith 5393a40ed3dSBarry Smith PetscFunctionBegin; 5400700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 5410700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 5420700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,3); 54317186662SBarry Smith if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()"); 544e0907662SLois Curfman McInnes 545d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 5464c49b128SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 54790d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr); 548f1af5d2fSBarry Smith if (flg) { 549ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 550e0907662SLois Curfman McInnes } else { 5510b9b6f31SBarry Smith PetscTruth assembled; 5520b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 5530b9b6f31SBarry Smith if (assembled) { 55443a90d84SBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 555e0907662SLois Curfman McInnes } 5560b9b6f31SBarry Smith } 55743a90d84SBarry Smith 558b8f8c88eSHong Zhang x1_tmp = x1; 559dac7f5fdSBarry Smith if (!coloring->vscale){ 560b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 561b8f8c88eSHong Zhang } 562be2fbe1fSBarry Smith 5633a7fca6bSBarry Smith /* 5643a7fca6bSBarry Smith This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets 5653a7fca6bSBarry Smith coloring->F for the coarser grids from the finest 5663a7fca6bSBarry Smith */ 5673a7fca6bSBarry Smith if (coloring->F) { 5683a7fca6bSBarry Smith ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); 5693a7fca6bSBarry Smith ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); 5703a7fca6bSBarry Smith if (m1 != m2) { 5713a7fca6bSBarry Smith coloring->F = 0; 5723a7fca6bSBarry Smith } 5733a7fca6bSBarry Smith } 5743a7fca6bSBarry Smith 575b8f8c88eSHong Zhang if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 576b8f8c88eSHong Zhang ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 577b8f8c88eSHong Zhang } 578b8f8c88eSHong Zhang ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 579b8f8c88eSHong Zhang 580b8f8c88eSHong Zhang /* Set w1 = F(x1) */ 581be2fbe1fSBarry Smith if (coloring->F) { 582be2fbe1fSBarry Smith w1 = coloring->F; /* use already computed value of function */ 583be2fbe1fSBarry Smith coloring->F = 0; 584be2fbe1fSBarry Smith } else { 58566f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 586b8f8c88eSHong Zhang ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 58766f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 588be2fbe1fSBarry Smith } 58943a90d84SBarry Smith 590b8f8c88eSHong Zhang if (!coloring->w3) { 591b8f8c88eSHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 592b8f8c88eSHong Zhang ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 593efb30889SBarry Smith } 594b8f8c88eSHong Zhang w3 = coloring->w3; 595efb30889SBarry Smith 596b8f8c88eSHong Zhang /* Compute all the local scale factors, including ghost points */ 597b8f8c88eSHong Zhang ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 598b8f8c88eSHong Zhang ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 599b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 600b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED){ 601b8f8c88eSHong Zhang col_start = 0; col_end = N; 6028ee2e534SBarry Smith } else if (ctype == IS_COLORING_GLOBAL){ 603b8f8c88eSHong Zhang xx = xx - start; 604b8f8c88eSHong Zhang vscale_array = vscale_array - start; 605b8f8c88eSHong Zhang col_start = start; col_end = N + start; 606b8f8c88eSHong Zhang } 607b8f8c88eSHong Zhang for (col=col_start; col<col_end; col++){ 608b8f8c88eSHong Zhang /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 609efb30889SBarry Smith if (coloring->htype[0] == 'w') { 610efb30889SBarry Smith dx = 1.0 + unorm; 611efb30889SBarry Smith } else { 6129782cf98SBarry Smith dx = xx[col]; 613efb30889SBarry Smith } 6145b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 6159782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX) 6169782cf98SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 6179782cf98SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 6189782cf98SBarry Smith #else 6199782cf98SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 6209782cf98SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 6219782cf98SBarry Smith #endif 6229782cf98SBarry Smith dx *= epsilon; 62330b34957SBarry Smith vscale_array[col] = 1.0/dx; 6249782cf98SBarry Smith } 6258ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start; 626efb30889SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 6278ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL){ 62830b34957SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 62930b34957SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 630b8f8c88eSHong Zhang } 6319782cf98SBarry Smith 632b8f8c88eSHong Zhang if (coloring->vscaleforrow) { 633b8f8c88eSHong Zhang vscaleforrow = coloring->vscaleforrow; 63417186662SBarry Smith } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 635b0a32e0cSBarry Smith 6369782cf98SBarry Smith /* 63743a90d84SBarry Smith Loop over each color 63843a90d84SBarry Smith */ 639b8f8c88eSHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 64043a90d84SBarry Smith for (k=0; k<coloring->ncolors; k++) { 64149b058dcSBarry Smith coloring->currentcolor = k; 642b8f8c88eSHong Zhang ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 643b8f8c88eSHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 6448ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 64543a90d84SBarry Smith /* 646b8f8c88eSHong Zhang Loop over each column associated with color 647b8f8c88eSHong Zhang adding the perturbation to the vector w3. 64843a90d84SBarry Smith */ 64943a90d84SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 650b8f8c88eSHong Zhang col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 651efb30889SBarry Smith if (coloring->htype[0] == 'w') { 652efb30889SBarry Smith dx = 1.0 + unorm; 653efb30889SBarry Smith } else { 65442460c72SBarry Smith dx = xx[col]; 655efb30889SBarry Smith } 6565b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 657aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 658ae09f205SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 659ae09f205SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 66043a90d84SBarry Smith #else 661329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 662329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 66343a90d84SBarry Smith #endif 66443a90d84SBarry Smith dx *= epsilon; 665e32f2f54SBarry Smith if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 66642460c72SBarry Smith w3_array[col] += dx; 66743a90d84SBarry Smith } 6688ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 669b8f8c88eSHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 6703b28642cSBarry Smith 67143a90d84SBarry Smith /* 672b8f8c88eSHong Zhang Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 673b8f8c88eSHong Zhang w2 = F(x1 + dx) - F(x1) 67443a90d84SBarry Smith */ 67566f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 67643a90d84SBarry Smith ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 67766f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 678efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 6799782cf98SBarry Smith 68043a90d84SBarry Smith /* 681e0907662SLois Curfman McInnes Loop over rows of vector, putting results into Jacobian matrix 68243a90d84SBarry Smith */ 6833b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 68443a90d84SBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 685b8f8c88eSHong Zhang row = coloring->rows[k][l]; /* local row index */ 686b8f8c88eSHong Zhang col = coloring->columnsforrow[k][l]; /* global column index */ 6875904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 68843a90d84SBarry Smith srow = row + start; 68943a90d84SBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 69043a90d84SBarry Smith } 6913b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 692aeb7e98aSMatthew Knepley } /* endof for each color */ 6938ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) xx = xx + start; 69430b34957SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 695b8f8c88eSHong Zhang ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 696b8f8c88eSHong Zhang 697b8f8c88eSHong Zhang coloring->currentcolor = -1; 69843a90d84SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69943a90d84SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 700d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 701a2e34c3dSBarry Smith 70290d69ab7SBarry Smith flg = PETSC_FALSE; 70390d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr); 704a2e34c3dSBarry Smith if (flg) { 70595902228SMatthew Knepley ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr); 706a2e34c3dSBarry Smith } 707b9722508SLois Curfman McInnes ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr); 7083a40ed3dSBarry Smith PetscFunctionReturn(0); 70943a90d84SBarry Smith } 710840b8ebdSBarry Smith 7114a2ae208SSatish Balay #undef __FUNCT__ 7124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApplyTS" 713840b8ebdSBarry Smith /*@ 714840b8ebdSBarry Smith MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context 715840b8ebdSBarry Smith has been created, computes the Jacobian for a function via finite differences. 716840b8ebdSBarry Smith 717fee21e36SBarry Smith Collective on Mat, MatFDColoring, and Vec 718fee21e36SBarry Smith 719ef5ee4d1SLois Curfman McInnes Input Parameters: 7203b28642cSBarry Smith + mat - location to store Jacobian 721ef5ee4d1SLois Curfman McInnes . coloring - coloring context created with MatFDColoringCreate() 722ef5ee4d1SLois Curfman McInnes . x1 - location at which Jacobian is to be computed 7237850c7c0SBarry Smith - sctx - context required by function, if this is being used with the TS solver then it is TS object, otherwise it is null 724ef5ee4d1SLois Curfman McInnes 72515091d37SBarry Smith Level: intermediate 72615091d37SBarry Smith 7277850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction() 728840b8ebdSBarry Smith 729840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences 730840b8ebdSBarry Smith @*/ 731be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx) 732840b8ebdSBarry Smith { 7336849ba73SBarry Smith PetscErrorCode (*f)(void*,PetscReal,Vec,Vec,void*)=(PetscErrorCode (*)(void*,PetscReal,Vec,Vec,void *))coloring->f; 7346849ba73SBarry Smith PetscErrorCode ierr; 73538baddfdSBarry Smith PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow; 736efb30889SBarry Smith PetscScalar dx,*y,*xx,*w3_array; 737ea709b57SSatish Balay PetscScalar *vscale_array; 738329f5518SBarry Smith PetscReal epsilon = coloring->error_rel,umin = coloring->umin; 739ced766e8SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3; 740840b8ebdSBarry Smith void *fctx = coloring->fctx; 741f1af5d2fSBarry Smith PetscTruth flg; 742840b8ebdSBarry Smith 7433a40ed3dSBarry Smith PetscFunctionBegin; 7440700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 7450700a824SBarry Smith PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2); 7460700a824SBarry Smith PetscValidHeaderSpecific(x1,VEC_CLASSID,4); 747840b8ebdSBarry Smith 748d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 749ced766e8SHong Zhang if (!coloring->w3) { 750840b8ebdSBarry Smith ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 75152e6d16bSBarry Smith ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); 752840b8ebdSBarry Smith } 753ced766e8SHong Zhang w3 = coloring->w3; 754840b8ebdSBarry Smith 7555904e1b1SBarry Smith ierr = MatSetUnfactored(J);CHKERRQ(ierr); 75690d69ab7SBarry Smith flg = PETSC_FALSE; 75790d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr); 758f1af5d2fSBarry Smith if (flg) { 759ae15b995SBarry Smith ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 760840b8ebdSBarry Smith } else { 7610b9b6f31SBarry Smith PetscTruth assembled; 7620b9b6f31SBarry Smith ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 7630b9b6f31SBarry Smith if (assembled) { 764840b8ebdSBarry Smith ierr = MatZeroEntries(J);CHKERRQ(ierr); 765840b8ebdSBarry Smith } 7660b9b6f31SBarry Smith } 767840b8ebdSBarry Smith 768840b8ebdSBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 769840b8ebdSBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 77066f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 771840b8ebdSBarry Smith ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr); 77266f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 773840b8ebdSBarry Smith 774840b8ebdSBarry Smith /* 7755904e1b1SBarry Smith Compute all the scale factors and share with other processors 776840b8ebdSBarry Smith */ 7775904e1b1SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; 7785904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; 779840b8ebdSBarry Smith for (k=0; k<coloring->ncolors; k++) { 780840b8ebdSBarry Smith /* 781840b8ebdSBarry Smith Loop over each column associated with color adding the 782840b8ebdSBarry Smith perturbation to the vector w3. 783840b8ebdSBarry Smith */ 784840b8ebdSBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 785840b8ebdSBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 7865904e1b1SBarry Smith dx = xx[col]; 7875b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 788aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 789840b8ebdSBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 790840b8ebdSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 791840b8ebdSBarry Smith #else 792329f5518SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 793329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 794840b8ebdSBarry Smith #endif 795840b8ebdSBarry Smith dx *= epsilon; 7965904e1b1SBarry Smith vscale_array[col] = 1.0/dx; 797840b8ebdSBarry Smith } 7985904e1b1SBarry Smith } 7995904e1b1SBarry Smith vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8005904e1b1SBarry Smith ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8015904e1b1SBarry Smith ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8025904e1b1SBarry Smith 8035904e1b1SBarry Smith if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; 8045904e1b1SBarry Smith else vscaleforrow = coloring->columnsforrow; 8055904e1b1SBarry Smith 8065904e1b1SBarry Smith ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8075904e1b1SBarry Smith /* 8085904e1b1SBarry Smith Loop over each color 8095904e1b1SBarry Smith */ 8105904e1b1SBarry Smith for (k=0; k<coloring->ncolors; k++) { 8115904e1b1SBarry Smith ierr = VecCopy(x1,w3);CHKERRQ(ierr); 8125904e1b1SBarry Smith ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; 8135904e1b1SBarry Smith /* 8145904e1b1SBarry Smith Loop over each column associated with color adding the 8155904e1b1SBarry Smith perturbation to the vector w3. 8165904e1b1SBarry Smith */ 8175904e1b1SBarry Smith for (l=0; l<coloring->ncolumns[k]; l++) { 8185904e1b1SBarry Smith col = coloring->columns[k][l]; /* column of the matrix we are probing for */ 8195904e1b1SBarry Smith dx = xx[col]; 8205b8514ebSBarry Smith if (dx == 0.0) dx = 1.0; 8215904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 8225904e1b1SBarry Smith if (dx < umin && dx >= 0.0) dx = umin; 8235904e1b1SBarry Smith else if (dx < 0.0 && dx > -umin) dx = -umin; 8245904e1b1SBarry Smith #else 8255904e1b1SBarry Smith if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 8265904e1b1SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 8275904e1b1SBarry Smith #endif 8285904e1b1SBarry Smith dx *= epsilon; 8295904e1b1SBarry Smith w3_array[col] += dx; 8305904e1b1SBarry Smith } 8315904e1b1SBarry Smith w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 8325904e1b1SBarry Smith 833840b8ebdSBarry Smith /* 834840b8ebdSBarry Smith Evaluate function at x1 + dx (here dx is a vector of perturbations) 835840b8ebdSBarry Smith */ 83666f9b7ceSBarry Smith ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 837840b8ebdSBarry Smith ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr); 83866f9b7ceSBarry Smith ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 839efb30889SBarry Smith ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 8405904e1b1SBarry Smith 841840b8ebdSBarry Smith /* 842840b8ebdSBarry Smith Loop over rows of vector, putting results into Jacobian matrix 843840b8ebdSBarry Smith */ 8443b28642cSBarry Smith ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 845840b8ebdSBarry Smith for (l=0; l<coloring->nrows[k]; l++) { 846840b8ebdSBarry Smith row = coloring->rows[k][l]; 847840b8ebdSBarry Smith col = coloring->columnsforrow[k][l]; 8485904e1b1SBarry Smith y[row] *= vscale_array[vscaleforrow[k][l]]; 849840b8ebdSBarry Smith srow = row + start; 850840b8ebdSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 851840b8ebdSBarry Smith } 8523b28642cSBarry Smith ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 853840b8ebdSBarry Smith } 8545904e1b1SBarry Smith ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8555904e1b1SBarry Smith xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 856840b8ebdSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 857840b8ebdSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 858d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr); 8593a40ed3dSBarry Smith PetscFunctionReturn(0); 860840b8ebdSBarry Smith } 8613b28642cSBarry Smith 8623b28642cSBarry Smith 8633b28642cSBarry Smith 864