1*4036925cSBarry Smith 2*4036925cSBarry Smith #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3*4036925cSBarry Smith 4*4036925cSBarry Smith /* ------------------------------------------------------------------------*/ 5*4036925cSBarry Smith struct _n_TSMonitorSPEigCtx { 6*4036925cSBarry Smith PetscDrawSP drawsp; 7*4036925cSBarry Smith KSP ksp; 8*4036925cSBarry Smith PetscRandom rand; 9*4036925cSBarry Smith PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 10*4036925cSBarry Smith }; 11*4036925cSBarry Smith 12*4036925cSBarry Smith 13*4036925cSBarry Smith #undef __FUNCT__ 14*4036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxCreate" 15*4036925cSBarry Smith /*@C 16*4036925cSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator 17*4036925cSBarry Smith 18*4036925cSBarry Smith Collective on TS 19*4036925cSBarry Smith 20*4036925cSBarry Smith Input Parameters: 21*4036925cSBarry Smith + host - the X display to open, or null for the local machine 22*4036925cSBarry Smith . label - the title to put in the title bar 23*4036925cSBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 24*4036925cSBarry Smith . m, n - the screen width and height in pixels 25*4036925cSBarry Smith - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 26*4036925cSBarry Smith 27*4036925cSBarry Smith Output Parameter: 28*4036925cSBarry Smith . ctx - the context 29*4036925cSBarry Smith 30*4036925cSBarry Smith Options Database Key: 31*4036925cSBarry Smith . -ts_monitor_eig - plot egienvalues of linearized right hand side 32*4036925cSBarry Smith 33*4036925cSBarry Smith Notes: 34*4036925cSBarry Smith Use TSMonitorSPEigCtxDestroy() to destroy. 35*4036925cSBarry Smith 36*4036925cSBarry Smith Level: intermediate 37*4036925cSBarry Smith 38*4036925cSBarry Smith .keywords: TS, monitor, line graph, residual, seealso 39*4036925cSBarry Smith 40*4036925cSBarry Smith .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError() 41*4036925cSBarry Smith 42*4036925cSBarry Smith @*/ 43*4036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx) 44*4036925cSBarry Smith { 45*4036925cSBarry Smith PetscDraw win; 46*4036925cSBarry Smith PetscErrorCode ierr; 47*4036925cSBarry Smith PC pc; 48*4036925cSBarry Smith 49*4036925cSBarry Smith PetscFunctionBegin; 50*4036925cSBarry Smith ierr = PetscNew(struct _n_TSMonitorSPEigCtx,ctx);CHKERRQ(ierr); 51*4036925cSBarry Smith ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr); 52*4036925cSBarry Smith ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr); 53*4036925cSBarry Smith ierr = PetscDrawSPCreate(win,1,&(*ctx)->drawsp);CHKERRQ(ierr); 54*4036925cSBarry Smith ierr = PetscRandomCreate(comm,&(*ctx)->rand);CHKERRQ(ierr); 55*4036925cSBarry Smith ierr = PetscRandomSetFromOptions((*ctx)->rand);CHKERRQ(ierr); 56*4036925cSBarry Smith ierr = KSPCreate(comm,&(*ctx)->ksp);CHKERRQ(ierr); 57*4036925cSBarry Smith ierr = KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_eig_");CHKERRQ(ierr); /* this is wrong, used use also prefix from the TS */ 58*4036925cSBarry Smith ierr = KSPSetType((*ctx)->ksp,KSPGMRES);CHKERRQ(ierr); 59*4036925cSBarry Smith ierr = KSPGMRESSetRestart((*ctx)->ksp,200);CHKERRQ(ierr); 60*4036925cSBarry Smith ierr = KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);CHKERRQ(ierr); 61*4036925cSBarry Smith ierr = KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);CHKERRQ(ierr); 62*4036925cSBarry Smith ierr = KSPSetFromOptions((*ctx)->ksp);CHKERRQ(ierr); 63*4036925cSBarry Smith ierr = KSPGetPC((*ctx)->ksp,&pc);CHKERRQ(ierr); 64*4036925cSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 65*4036925cSBarry Smith (*ctx)->howoften = howoften; 66*4036925cSBarry Smith PetscFunctionReturn(0); 67*4036925cSBarry Smith } 68*4036925cSBarry Smith 69*4036925cSBarry Smith #undef __FUNCT__ 70*4036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEig" 71*4036925cSBarry Smith PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx) 72*4036925cSBarry Smith { 73*4036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx; 74*4036925cSBarry Smith PetscErrorCode ierr; 75*4036925cSBarry Smith KSP ksp = ctx->ksp; 76*4036925cSBarry Smith PetscInt n,nits,neig,i; 77*4036925cSBarry Smith PetscReal *r,*c; 78*4036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 79*4036925cSBarry Smith MatStructure structure; 80*4036925cSBarry Smith Mat A,B; 81*4036925cSBarry Smith Vec xdot; 82*4036925cSBarry Smith SNES snes; 83*4036925cSBarry Smith 84*4036925cSBarry Smith PetscFunctionBegin; 85*4036925cSBarry Smith if (!step) PetscFunctionReturn(0); 86*4036925cSBarry Smith if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && (step == -1))){ 87*4036925cSBarry Smith ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr); 88*4036925cSBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 89*4036925cSBarry Smith ierr = SNESGetJacobian(snes,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 90*4036925cSBarry Smith ierr = TSComputeIJacobian(ts,ptime,v,xdot,1.0,&A,&B,&structure,PETSC_FALSE);CHKERRQ(ierr); 91*4036925cSBarry Smith ierr = KSPSetOperators(ksp,A,B,structure);CHKERRQ(ierr); 92*4036925cSBarry Smith ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr); 93*4036925cSBarry Smith ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr); 94*4036925cSBarry Smith ierr = VecDestroy(&xdot);CHKERRQ(ierr); 95*4036925cSBarry Smith ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr); 96*4036925cSBarry Smith n = nits+2; 97*4036925cSBarry Smith 98*4036925cSBarry Smith if (nits) { 99*4036925cSBarry Smith ierr = PetscMalloc2(n,PetscReal,&r,n,PetscReal,&c);CHKERRQ(ierr); 100*4036925cSBarry Smith ierr = KSPComputeEigenvalues(ksp,n,r,c,&neig);CHKERRQ(ierr); 101*4036925cSBarry Smith for (i=0; i<neig; i++) { 102*4036925cSBarry Smith ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr); 103*4036925cSBarry Smith } 104*4036925cSBarry Smith ierr = PetscDrawSPDraw(drawsp);CHKERRQ(ierr); 105*4036925cSBarry Smith ierr = PetscFree2(r,c);CHKERRQ(ierr); 106*4036925cSBarry Smith } 107*4036925cSBarry Smith } 108*4036925cSBarry Smith PetscFunctionReturn(0); 109*4036925cSBarry Smith } 110*4036925cSBarry Smith 111*4036925cSBarry Smith #undef __FUNCT__ 112*4036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy" 113*4036925cSBarry Smith /*@C 114*4036925cSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate(). 115*4036925cSBarry Smith 116*4036925cSBarry Smith Collective on TSMonitorSPEigCtx 117*4036925cSBarry Smith 118*4036925cSBarry Smith Input Parameter: 119*4036925cSBarry Smith . ctx - the monitor context 120*4036925cSBarry Smith 121*4036925cSBarry Smith Level: intermediate 122*4036925cSBarry Smith 123*4036925cSBarry Smith .keywords: TS, monitor, line graph, destroy 124*4036925cSBarry Smith 125*4036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(), TSMonitorSet(), TSMonitorSPEig(); 126*4036925cSBarry Smith @*/ 127*4036925cSBarry Smith PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) 128*4036925cSBarry Smith { 129*4036925cSBarry Smith PetscDraw draw; 130*4036925cSBarry Smith PetscErrorCode ierr; 131*4036925cSBarry Smith 132*4036925cSBarry Smith PetscFunctionBegin; 133*4036925cSBarry Smith ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr); 134*4036925cSBarry Smith ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 135*4036925cSBarry Smith ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr); 136*4036925cSBarry Smith ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr); 137*4036925cSBarry Smith ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr); 138*4036925cSBarry Smith ierr = PetscFree(*ctx);CHKERRQ(ierr); 139*4036925cSBarry Smith PetscFunctionReturn(0); 140*4036925cSBarry Smith } 141*4036925cSBarry Smith 142