10d18c744SBarry Smith #define PETSC_DESIRE_COMPLEX 24036925cSBarry Smith #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 34036925cSBarry Smith 44036925cSBarry Smith /* ------------------------------------------------------------------------*/ 54036925cSBarry Smith struct _n_TSMonitorSPEigCtx { 64036925cSBarry Smith PetscDrawSP drawsp; 74036925cSBarry Smith KSP ksp; 84036925cSBarry Smith PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 9b51f60e1SBarry Smith PetscBool computeexplicitly; 10b51f60e1SBarry Smith MPI_Comm comm; 11b51f60e1SBarry Smith PetscRandom rand; 120d18c744SBarry Smith PetscReal xmin,xmax,ymin,ymax; 134036925cSBarry Smith }; 144036925cSBarry Smith 154036925cSBarry Smith 164036925cSBarry Smith #undef __FUNCT__ 174036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxCreate" 184036925cSBarry Smith /*@C 194036925cSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator 204036925cSBarry Smith 214036925cSBarry Smith Collective on TS 224036925cSBarry Smith 234036925cSBarry Smith Input Parameters: 244036925cSBarry Smith + host - the X display to open, or null for the local machine 254036925cSBarry Smith . label - the title to put in the title bar 264036925cSBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 274036925cSBarry Smith . m, n - the screen width and height in pixels 284036925cSBarry Smith - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 294036925cSBarry Smith 304036925cSBarry Smith Output Parameter: 314036925cSBarry Smith . ctx - the context 324036925cSBarry Smith 334036925cSBarry Smith Options Database Key: 34b51f60e1SBarry Smith . -ts_monitor_sp_eig - plot egienvalues of linearized right hand side 354036925cSBarry Smith 364036925cSBarry Smith Notes: 374036925cSBarry Smith Use TSMonitorSPEigCtxDestroy() to destroy. 384036925cSBarry Smith 39b51f60e1SBarry Smith Currently only works if the Jacobian is provided explicitly. 40b51f60e1SBarry Smith 41b51f60e1SBarry Smith Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix. 42b51f60e1SBarry Smith 434036925cSBarry Smith Level: intermediate 444036925cSBarry Smith 454036925cSBarry Smith .keywords: TS, monitor, line graph, residual, seealso 464036925cSBarry Smith 474036925cSBarry Smith .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError() 484036925cSBarry Smith 494036925cSBarry Smith @*/ 504036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx) 514036925cSBarry Smith { 524036925cSBarry Smith PetscDraw win; 534036925cSBarry Smith PetscErrorCode ierr; 544036925cSBarry Smith PC pc; 554036925cSBarry Smith 564036925cSBarry Smith PetscFunctionBegin; 574036925cSBarry Smith ierr = PetscNew(struct _n_TSMonitorSPEigCtx,ctx);CHKERRQ(ierr); 58b51f60e1SBarry Smith ierr = PetscRandomCreate(comm,&(*ctx)->rand);CHKERRQ(ierr); 59b51f60e1SBarry Smith ierr = PetscRandomSetFromOptions((*ctx)->rand);CHKERRQ(ierr); 604036925cSBarry Smith ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr); 614036925cSBarry Smith ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr); 624036925cSBarry Smith ierr = PetscDrawSPCreate(win,1,&(*ctx)->drawsp);CHKERRQ(ierr); 634036925cSBarry Smith ierr = KSPCreate(comm,&(*ctx)->ksp);CHKERRQ(ierr); 64b51f60e1SBarry Smith ierr = KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_sp_eig_");CHKERRQ(ierr); /* this is wrong, used use also prefix from the TS */ 654036925cSBarry Smith ierr = KSPSetType((*ctx)->ksp,KSPGMRES);CHKERRQ(ierr); 664036925cSBarry Smith ierr = KSPGMRESSetRestart((*ctx)->ksp,200);CHKERRQ(ierr); 674036925cSBarry Smith ierr = KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);CHKERRQ(ierr); 684036925cSBarry Smith ierr = KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);CHKERRQ(ierr); 694036925cSBarry Smith ierr = KSPSetFromOptions((*ctx)->ksp);CHKERRQ(ierr); 704036925cSBarry Smith ierr = KSPGetPC((*ctx)->ksp,&pc);CHKERRQ(ierr); 714036925cSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 724036925cSBarry Smith (*ctx)->howoften = howoften; 73b51f60e1SBarry Smith (*ctx)->computeexplicitly = PETSC_FALSE; 74b51f60e1SBarry Smith ierr = PetscOptionsGetBool(PETSC_NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,PETSC_NULL);CHKERRQ(ierr); 75b51f60e1SBarry Smith (*ctx)->comm = comm; 760d18c744SBarry Smith (*ctx)->xmin = -2.1; 77f9c1d6abSBarry Smith (*ctx)->xmax = 1.1; 780d18c744SBarry Smith (*ctx)->ymin = -1.1; 790d18c744SBarry Smith (*ctx)->ymax = 1.1; 804036925cSBarry Smith PetscFunctionReturn(0); 814036925cSBarry Smith } 824036925cSBarry Smith 83f9c1d6abSBarry Smith #undef __FUNCT__ 84f9c1d6abSBarry Smith #define __FUNCT__ "TSLinearStabilityIndicator" 85f9c1d6abSBarry Smith static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr,PetscReal xi,PetscBool *flg) 860d18c744SBarry Smith { 87f9c1d6abSBarry Smith PetscErrorCode ierr; 88f9c1d6abSBarry Smith PetscReal yr,yi; 890d18c744SBarry Smith 90f9c1d6abSBarry Smith PetscFunctionBegin; 91f9c1d6abSBarry Smith ierr = TSComputeLinearStability(ts,xr,xi,&yr,&yi);CHKERRQ(ierr); 92f9c1d6abSBarry Smith if ((yr*yr + yi*yi) <= 1.0) *flg = PETSC_TRUE; 93f9c1d6abSBarry Smith else *flg = PETSC_FALSE; 94f9c1d6abSBarry Smith PetscFunctionReturn(0); 950d18c744SBarry Smith } 960d18c744SBarry Smith 974036925cSBarry Smith #undef __FUNCT__ 984036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEig" 994036925cSBarry Smith PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx) 1004036925cSBarry Smith { 1014036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx; 1024036925cSBarry Smith PetscErrorCode ierr; 1034036925cSBarry Smith KSP ksp = ctx->ksp; 104b51f60e1SBarry Smith PetscInt n,N,nits,neig,i,its = 200; 105*b296d7d5SJed Brown PetscReal *r,*c,time_step_save; 1064036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 1074036925cSBarry Smith MatStructure structure; 1084036925cSBarry Smith Mat A,B; 1094036925cSBarry Smith Vec xdot; 1104036925cSBarry Smith SNES snes; 1114036925cSBarry Smith 1124036925cSBarry Smith PetscFunctionBegin; 1134036925cSBarry Smith if (!step) PetscFunctionReturn(0); 1144036925cSBarry Smith if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && (step == -1))){ 1154036925cSBarry Smith ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr); 1164036925cSBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1174036925cSBarry Smith ierr = SNESGetJacobian(snes,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 118a174af7bSBarry Smith ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);CHKERRQ(ierr); 119b51f60e1SBarry Smith /* 120b51f60e1SBarry Smith This doesn't work because methods keep and use internal information about the shift so it 121b51f60e1SBarry Smith seems we would need code for each method to trick the correct Jacobian in being computed. 122*b296d7d5SJed Brown */ 123*b296d7d5SJed Brown time_step_save = ts->time_step; 124b51f60e1SBarry Smith ts->time_step = PETSC_MAX_REAL; 1253d39420eSBarry Smith ierr = SNESComputeJacobian(snes,v,&A,&B,&structure);CHKERRQ(ierr); 126*b296d7d5SJed Brown ts->time_step = time_step_save; 127a174af7bSBarry Smith 128a174af7bSBarry Smith ierr = KSPSetOperators(ksp,B,B,structure);CHKERRQ(ierr); 129b51f60e1SBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 130b51f60e1SBarry Smith if (n < 200) its = n; 131b51f60e1SBarry Smith ierr = KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);CHKERRQ(ierr); 1324036925cSBarry Smith ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr); 1334036925cSBarry Smith ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr); 1344036925cSBarry Smith ierr = VecDestroy(&xdot);CHKERRQ(ierr); 1354036925cSBarry Smith ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr); 136b51f60e1SBarry Smith N = nits+2; 1374036925cSBarry Smith 1384036925cSBarry Smith if (nits) { 139a174af7bSBarry Smith PetscDraw draw; 140d98bdc7eSBarry Smith PetscReal pause; 141f9c1d6abSBarry Smith PetscDrawAxis axis; 142f9c1d6abSBarry Smith PetscReal xmin,xmax,ymin,ymax; 143d98bdc7eSBarry Smith 144b51f60e1SBarry Smith ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr); 1450d18c744SBarry Smith ierr = PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);CHKERRQ(ierr); 146a174af7bSBarry Smith ierr = PetscMalloc2(PetscMax(n,N),PetscReal,&r,PetscMax(n,N),PetscReal,&c);CHKERRQ(ierr); 147b51f60e1SBarry Smith if (ctx->computeexplicitly) { 148b51f60e1SBarry Smith ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr); 149b51f60e1SBarry Smith neig = n; 150b51f60e1SBarry Smith } else { 151b51f60e1SBarry Smith ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr); 152b51f60e1SBarry Smith } 153*b296d7d5SJed Brown /* We used the positive operator to be able to reuse KSPs that require positive definiteness, now flip the spectrum as is conventional for ODEs */ 154*b296d7d5SJed Brown for (i=0; i<neig; i++) r[i] = -r[i]; 1554036925cSBarry Smith for (i=0; i<neig; i++) { 1565f737190SBarry Smith if (ts->ops->linearstability) { 1575f737190SBarry Smith PetscReal fr,fi; 1585f737190SBarry Smith ierr = TSComputeLinearStability(ts,r[i],c[i],&fr,&fi);CHKERRQ(ierr); 1595f737190SBarry Smith if ((fr*fr + fi*fi) > 1.0) { 1605f737190SBarry Smith ierr = PetscPrintf(ctx->comm,"Linearized Eigenvalue %g + %g i linear stability function %g norm indicates unstable scheme \n",(double)r[i],(double)c[i],(double)(fr*fr + fi*fi));CHKERRQ(ierr); 1615f737190SBarry Smith } 1625d726895SBarry Smith } 1634036925cSBarry Smith ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr); 1644036925cSBarry Smith } 1650d18c744SBarry Smith ierr = PetscFree2(r,c);CHKERRQ(ierr); 166a174af7bSBarry Smith ierr = PetscDrawSPGetDraw(drawsp,&draw);CHKERRQ(ierr); 167d98bdc7eSBarry Smith ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr); 168d98bdc7eSBarry Smith ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr); 169d98bdc7eSBarry Smith ierr = PetscDrawSPDraw(drawsp,PETSC_TRUE);CHKERRQ(ierr); 170d98bdc7eSBarry Smith ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr); 1710d18c744SBarry Smith 172f9c1d6abSBarry Smith if (ts->ops->linearstability) { 1730d18c744SBarry Smith ierr = PetscDrawSPGetAxis(drawsp,&axis);CHKERRQ(ierr); 1740d18c744SBarry Smith ierr = PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);CHKERRQ(ierr); 175f9c1d6abSBarry Smith ierr = PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,(PetscErrorCode (*)(void*,PetscReal,PetscReal,PetscBool*))TSLinearStabilityIndicator,ts);CHKERRQ(ierr); 176a174af7bSBarry Smith ierr = PetscDrawSPDraw(drawsp,PETSC_FALSE);CHKERRQ(ierr); 1770d18c744SBarry Smith } 1784036925cSBarry Smith } 179a174af7bSBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 1804036925cSBarry Smith } 1814036925cSBarry Smith PetscFunctionReturn(0); 1824036925cSBarry Smith } 1834036925cSBarry Smith 1844036925cSBarry Smith #undef __FUNCT__ 1854036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy" 1864036925cSBarry Smith /*@C 1874036925cSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate(). 1884036925cSBarry Smith 1894036925cSBarry Smith Collective on TSMonitorSPEigCtx 1904036925cSBarry Smith 1914036925cSBarry Smith Input Parameter: 1924036925cSBarry Smith . ctx - the monitor context 1934036925cSBarry Smith 1944036925cSBarry Smith Level: intermediate 1954036925cSBarry Smith 1964036925cSBarry Smith .keywords: TS, monitor, line graph, destroy 1974036925cSBarry Smith 1984036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(), TSMonitorSet(), TSMonitorSPEig(); 1994036925cSBarry Smith @*/ 2004036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) 2014036925cSBarry Smith { 2024036925cSBarry Smith PetscDraw draw; 2034036925cSBarry Smith PetscErrorCode ierr; 2044036925cSBarry Smith 2054036925cSBarry Smith PetscFunctionBegin; 2064036925cSBarry Smith ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr); 2074036925cSBarry Smith ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 2084036925cSBarry Smith ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr); 2094036925cSBarry Smith ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr); 2104036925cSBarry Smith ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr); 2114036925cSBarry Smith ierr = PetscFree(*ctx);CHKERRQ(ierr); 2124036925cSBarry Smith PetscFunctionReturn(0); 2134036925cSBarry Smith } 2144036925cSBarry Smith 2150d18c744SBarry Smith 2160d18c744SBarry Smith 217