1*8cd53115SBarry Smith 24036925cSBarry Smith #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 39804daf3SBarry Smith #include <petscdraw.h> 44036925cSBarry Smith 54036925cSBarry Smith /* ------------------------------------------------------------------------*/ 64036925cSBarry Smith struct _n_TSMonitorSPEigCtx { 74036925cSBarry Smith PetscDrawSP drawsp; 84036925cSBarry Smith KSP ksp; 94036925cSBarry Smith PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 10b51f60e1SBarry Smith PetscBool computeexplicitly; 11b51f60e1SBarry Smith MPI_Comm comm; 12b51f60e1SBarry Smith PetscRandom rand; 130d18c744SBarry Smith PetscReal xmin,xmax,ymin,ymax; 144036925cSBarry Smith }; 154036925cSBarry Smith 164036925cSBarry Smith 174036925cSBarry Smith #undef __FUNCT__ 184036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxCreate" 194036925cSBarry Smith /*@C 204036925cSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator 214036925cSBarry Smith 224036925cSBarry Smith Collective on TS 234036925cSBarry Smith 244036925cSBarry Smith Input Parameters: 254036925cSBarry Smith + host - the X display to open, or null for the local machine 264036925cSBarry Smith . label - the title to put in the title bar 274036925cSBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 284036925cSBarry Smith . m, n - the screen width and height in pixels 294036925cSBarry Smith - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 304036925cSBarry Smith 314036925cSBarry Smith Output Parameter: 324036925cSBarry Smith . ctx - the context 334036925cSBarry Smith 344036925cSBarry Smith Options Database Key: 35b51f60e1SBarry Smith . -ts_monitor_sp_eig - plot egienvalues of linearized right hand side 364036925cSBarry Smith 374036925cSBarry Smith Notes: 384036925cSBarry Smith Use TSMonitorSPEigCtxDestroy() to destroy. 394036925cSBarry Smith 40b51f60e1SBarry Smith Currently only works if the Jacobian is provided explicitly. 41b51f60e1SBarry Smith 42b51f60e1SBarry Smith Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix. 43b51f60e1SBarry Smith 444036925cSBarry Smith Level: intermediate 454036925cSBarry Smith 464036925cSBarry Smith .keywords: TS, monitor, line graph, residual, seealso 474036925cSBarry Smith 484036925cSBarry Smith .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError() 494036925cSBarry Smith 504036925cSBarry Smith @*/ 514036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx) 524036925cSBarry Smith { 534036925cSBarry Smith PetscDraw win; 544036925cSBarry Smith PetscErrorCode ierr; 554036925cSBarry Smith PC pc; 564036925cSBarry Smith 574036925cSBarry Smith PetscFunctionBegin; 58b00a9115SJed Brown ierr = PetscNew(ctx);CHKERRQ(ierr); 59b51f60e1SBarry Smith ierr = PetscRandomCreate(comm,&(*ctx)->rand);CHKERRQ(ierr); 60b51f60e1SBarry Smith ierr = PetscRandomSetFromOptions((*ctx)->rand);CHKERRQ(ierr); 614036925cSBarry Smith ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr); 624036925cSBarry Smith ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr); 634036925cSBarry Smith ierr = PetscDrawSPCreate(win,1,&(*ctx)->drawsp);CHKERRQ(ierr); 644036925cSBarry Smith ierr = KSPCreate(comm,&(*ctx)->ksp);CHKERRQ(ierr); 65b51f60e1SBarry Smith ierr = KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_sp_eig_");CHKERRQ(ierr); /* this is wrong, used use also prefix from the TS */ 664036925cSBarry Smith ierr = KSPSetType((*ctx)->ksp,KSPGMRES);CHKERRQ(ierr); 674036925cSBarry Smith ierr = KSPGMRESSetRestart((*ctx)->ksp,200);CHKERRQ(ierr); 684036925cSBarry Smith ierr = KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);CHKERRQ(ierr); 694036925cSBarry Smith ierr = KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);CHKERRQ(ierr); 704036925cSBarry Smith ierr = KSPSetFromOptions((*ctx)->ksp);CHKERRQ(ierr); 714036925cSBarry Smith ierr = KSPGetPC((*ctx)->ksp,&pc);CHKERRQ(ierr); 724036925cSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 73bbd56ea5SKarl Rupp 744036925cSBarry Smith (*ctx)->howoften = howoften; 75b51f60e1SBarry Smith (*ctx)->computeexplicitly = PETSC_FALSE; 76bbd56ea5SKarl Rupp 770298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,NULL);CHKERRQ(ierr); 78bbd56ea5SKarl Rupp 79b51f60e1SBarry Smith (*ctx)->comm = comm; 800d18c744SBarry Smith (*ctx)->xmin = -2.1; 81f9c1d6abSBarry Smith (*ctx)->xmax = 1.1; 820d18c744SBarry Smith (*ctx)->ymin = -1.1; 830d18c744SBarry Smith (*ctx)->ymax = 1.1; 844036925cSBarry Smith PetscFunctionReturn(0); 854036925cSBarry Smith } 864036925cSBarry Smith 87f9c1d6abSBarry Smith #undef __FUNCT__ 88f9c1d6abSBarry Smith #define __FUNCT__ "TSLinearStabilityIndicator" 89f9c1d6abSBarry Smith static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr,PetscReal xi,PetscBool *flg) 900d18c744SBarry Smith { 91f9c1d6abSBarry Smith PetscErrorCode ierr; 92f9c1d6abSBarry Smith PetscReal yr,yi; 930d18c744SBarry Smith 94f9c1d6abSBarry Smith PetscFunctionBegin; 95f9c1d6abSBarry Smith ierr = TSComputeLinearStability(ts,xr,xi,&yr,&yi);CHKERRQ(ierr); 96f9c1d6abSBarry Smith if ((yr*yr + yi*yi) <= 1.0) *flg = PETSC_TRUE; 97f9c1d6abSBarry Smith else *flg = PETSC_FALSE; 98f9c1d6abSBarry Smith PetscFunctionReturn(0); 990d18c744SBarry Smith } 1000d18c744SBarry Smith 1014036925cSBarry Smith #undef __FUNCT__ 1024036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEig" 1034036925cSBarry Smith PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx) 1044036925cSBarry Smith { 1054036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx; 1064036925cSBarry Smith PetscErrorCode ierr; 1074036925cSBarry Smith KSP ksp = ctx->ksp; 108b51f60e1SBarry Smith PetscInt n,N,nits,neig,i,its = 200; 109b296d7d5SJed Brown PetscReal *r,*c,time_step_save; 1104036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 1114036925cSBarry Smith Mat A,B; 1124036925cSBarry Smith Vec xdot; 1134036925cSBarry Smith SNES snes; 1144036925cSBarry Smith 1154036925cSBarry Smith PetscFunctionBegin; 1164036925cSBarry Smith if (!step) PetscFunctionReturn(0); 117b06615a5SLisandro Dalcin if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1184036925cSBarry Smith ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr); 1194036925cSBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1200298fd71SBarry Smith ierr = SNESGetJacobian(snes,&A,&B,NULL,NULL);CHKERRQ(ierr); 121a174af7bSBarry Smith ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);CHKERRQ(ierr); 122b51f60e1SBarry Smith /* 123b51f60e1SBarry Smith This doesn't work because methods keep and use internal information about the shift so it 124b51f60e1SBarry Smith seems we would need code for each method to trick the correct Jacobian in being computed. 125b296d7d5SJed Brown */ 126b296d7d5SJed Brown time_step_save = ts->time_step; 127b51f60e1SBarry Smith ts->time_step = PETSC_MAX_REAL; 128bbd56ea5SKarl Rupp 129d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,v,A,B);CHKERRQ(ierr); 130bbd56ea5SKarl Rupp 131b296d7d5SJed Brown ts->time_step = time_step_save; 132a174af7bSBarry Smith 13323ee1639SBarry Smith ierr = KSPSetOperators(ksp,B,B);CHKERRQ(ierr); 134b51f60e1SBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 135b51f60e1SBarry Smith if (n < 200) its = n; 136b51f60e1SBarry Smith ierr = KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);CHKERRQ(ierr); 1374036925cSBarry Smith ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr); 1384036925cSBarry Smith ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr); 1394036925cSBarry Smith ierr = VecDestroy(&xdot);CHKERRQ(ierr); 1404036925cSBarry Smith ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr); 141b51f60e1SBarry Smith N = nits+2; 1424036925cSBarry Smith 1434036925cSBarry Smith if (nits) { 144a174af7bSBarry Smith PetscDraw draw; 145d98bdc7eSBarry Smith PetscReal pause; 146f9c1d6abSBarry Smith PetscDrawAxis axis; 147f9c1d6abSBarry Smith PetscReal xmin,xmax,ymin,ymax; 148d98bdc7eSBarry Smith 149b51f60e1SBarry Smith ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr); 1500d18c744SBarry Smith ierr = PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);CHKERRQ(ierr); 151dcca6d9dSJed Brown ierr = PetscMalloc2(PetscMax(n,N),&r,PetscMax(n,N),&c);CHKERRQ(ierr); 152b51f60e1SBarry Smith if (ctx->computeexplicitly) { 153b51f60e1SBarry Smith ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr); 154b51f60e1SBarry Smith neig = n; 155b51f60e1SBarry Smith } else { 156b51f60e1SBarry Smith ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr); 157b51f60e1SBarry Smith } 158b296d7d5SJed 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 */ 159b296d7d5SJed Brown for (i=0; i<neig; i++) r[i] = -r[i]; 1604036925cSBarry Smith for (i=0; i<neig; i++) { 1615f737190SBarry Smith if (ts->ops->linearstability) { 1625f737190SBarry Smith PetscReal fr,fi; 1635f737190SBarry Smith ierr = TSComputeLinearStability(ts,r[i],c[i],&fr,&fi);CHKERRQ(ierr); 1645f737190SBarry Smith if ((fr*fr + fi*fi) > 1.0) { 1655f737190SBarry 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); 1665f737190SBarry Smith } 1675d726895SBarry Smith } 1684036925cSBarry Smith ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr); 1694036925cSBarry Smith } 1700d18c744SBarry Smith ierr = PetscFree2(r,c);CHKERRQ(ierr); 171a174af7bSBarry Smith ierr = PetscDrawSPGetDraw(drawsp,&draw);CHKERRQ(ierr); 172d98bdc7eSBarry Smith ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr); 173d98bdc7eSBarry Smith ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr); 174d98bdc7eSBarry Smith ierr = PetscDrawSPDraw(drawsp,PETSC_TRUE);CHKERRQ(ierr); 175d98bdc7eSBarry Smith ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr); 1760d18c744SBarry Smith 177f9c1d6abSBarry Smith if (ts->ops->linearstability) { 1780d18c744SBarry Smith ierr = PetscDrawSPGetAxis(drawsp,&axis);CHKERRQ(ierr); 1790d18c744SBarry Smith ierr = PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);CHKERRQ(ierr); 180f9c1d6abSBarry Smith ierr = PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,(PetscErrorCode (*)(void*,PetscReal,PetscReal,PetscBool*))TSLinearStabilityIndicator,ts);CHKERRQ(ierr); 181a174af7bSBarry Smith ierr = PetscDrawSPDraw(drawsp,PETSC_FALSE);CHKERRQ(ierr); 1820d18c744SBarry Smith } 1834036925cSBarry Smith } 184a174af7bSBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 1854036925cSBarry Smith } 1864036925cSBarry Smith PetscFunctionReturn(0); 1874036925cSBarry Smith } 1884036925cSBarry Smith 1894036925cSBarry Smith #undef __FUNCT__ 1904036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy" 1914036925cSBarry Smith /*@C 1924036925cSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate(). 1934036925cSBarry Smith 1944036925cSBarry Smith Collective on TSMonitorSPEigCtx 1954036925cSBarry Smith 1964036925cSBarry Smith Input Parameter: 1974036925cSBarry Smith . ctx - the monitor context 1984036925cSBarry Smith 1994036925cSBarry Smith Level: intermediate 2004036925cSBarry Smith 2014036925cSBarry Smith .keywords: TS, monitor, line graph, destroy 2024036925cSBarry Smith 2034036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(), TSMonitorSet(), TSMonitorSPEig(); 2044036925cSBarry Smith @*/ 2054036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) 2064036925cSBarry Smith { 2074036925cSBarry Smith PetscDraw draw; 2084036925cSBarry Smith PetscErrorCode ierr; 2094036925cSBarry Smith 2104036925cSBarry Smith PetscFunctionBegin; 2114036925cSBarry Smith ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr); 2124036925cSBarry Smith ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 2134036925cSBarry Smith ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr); 2144036925cSBarry Smith ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr); 2154036925cSBarry Smith ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr); 2164036925cSBarry Smith ierr = PetscFree(*ctx);CHKERRQ(ierr); 2174036925cSBarry Smith PetscFunctionReturn(0); 2184036925cSBarry Smith } 2194036925cSBarry Smith 2200d18c744SBarry Smith 2210d18c744SBarry Smith 222