xref: /petsc/src/ts/interface/tseig.c (revision d98bdc7e1efd5abbe89aa02f7a72e93d2ab575d3)
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;
12b51f60e1SBarry Smith   PetscScalar shift;
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;
584036925cSBarry Smith   ierr = PetscNew(struct _n_TSMonitorSPEigCtx,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);
734036925cSBarry Smith   (*ctx)->howoften          = howoften;
74b51f60e1SBarry Smith   (*ctx)->computeexplicitly = PETSC_FALSE;
75b51f60e1SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,PETSC_NULL);CHKERRQ(ierr);
76b51f60e1SBarry Smith   (*ctx)->shift             = 0.0;
77b51f60e1SBarry Smith   ierr = PetscOptionsGetScalar(PETSC_NULL,"-ts_monitor_sp_eig_shift",&(*ctx)->shift,PETSC_NULL);CHKERRQ(ierr);
78b51f60e1SBarry Smith   (*ctx)->comm              = comm;
790d18c744SBarry Smith   (*ctx)->xmin = -2.1;
800d18c744SBarry Smith   (*ctx)->xmax = .1;
810d18c744SBarry Smith   (*ctx)->ymin = -1.1;
820d18c744SBarry Smith   (*ctx)->ymax = 1.1;
834036925cSBarry Smith   PetscFunctionReturn(0);
844036925cSBarry Smith }
854036925cSBarry Smith 
860d18c744SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
870d18c744SBarry Smith PetscComplex Rtheta(PetscReal theta,PetscComplex z)
880d18c744SBarry Smith {
89*d98bdc7eSBarry Smith   return( (1.0 + (1.0 - theta)*z)/(1.0 - theta*z) );
900d18c744SBarry Smith }
910d18c744SBarry Smith 
920d18c744SBarry Smith PetscBool TSIndicator_Theta(PetscReal x,PetscReal y)
930d18c744SBarry Smith {
94*d98bdc7eSBarry Smith   PetscComplex r = Rtheta(.75,x + PETSC_i*y);
950d18c744SBarry Smith   if (PetscAbsComplex(r) <= 1.0) return PETSC_TRUE;
960d18c744SBarry Smith   else return PETSC_FALSE;
970d18c744SBarry Smith }
980d18c744SBarry Smith #endif
990d18c744SBarry Smith 
1004036925cSBarry Smith #undef __FUNCT__
1014036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEig"
1024036925cSBarry Smith PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
1034036925cSBarry Smith {
1044036925cSBarry Smith   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx;
1054036925cSBarry Smith   PetscErrorCode    ierr;
1064036925cSBarry Smith   KSP               ksp = ctx->ksp;
107b51f60e1SBarry Smith   PetscInt          n,N,nits,neig,i,its = 200;
1084036925cSBarry Smith   PetscReal         *r,*c;
1094036925cSBarry Smith   PetscDrawSP       drawsp = ctx->drawsp;
1104036925cSBarry Smith   MatStructure      structure;
1114036925cSBarry Smith   Mat               A,B;
1124036925cSBarry Smith   Vec               xdot;
1134036925cSBarry Smith   SNES              snes;
1144036925cSBarry Smith 
1154036925cSBarry Smith   PetscFunctionBegin;
1164036925cSBarry Smith   if (!step) PetscFunctionReturn(0);
1174036925cSBarry Smith   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && (step == -1))){
1184036925cSBarry Smith     ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr);
1194036925cSBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1204036925cSBarry Smith     ierr = SNESGetJacobian(snes,&A,&B,PETSC_NULL,PETSC_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.
125b51f60e1SBarry Smith     ts->time_step = PETSC_MAX_REAL;
1263d39420eSBarry Smith     ierr = SNESComputeJacobian(snes,v,&A,&B,&structure);CHKERRQ(ierr);
127b51f60e1SBarry Smith     */
128a174af7bSBarry Smith     ierr = TSComputeIJacobian(ts,ptime,v,xdot,0.0,&B,&B,&structure,PETSC_FALSE);CHKERRQ(ierr);
129a174af7bSBarry Smith     ierr = MatScale(B,-1.0*ts->time_step);CHKERRQ(ierr);
130a174af7bSBarry Smith     ierr = MatShift(B,ctx->shift);CHKERRQ(ierr);
131a174af7bSBarry Smith 
132a174af7bSBarry Smith     ierr = KSPSetOperators(ksp,B,B,structure);CHKERRQ(ierr);
133b51f60e1SBarry Smith     ierr = VecGetSize(v,&n);CHKERRQ(ierr);
134b51f60e1SBarry Smith     if (n < 200) its = n;
135b51f60e1SBarry Smith     ierr = KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);CHKERRQ(ierr);
1364036925cSBarry Smith     ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr);
1374036925cSBarry Smith     ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr);
1384036925cSBarry Smith     ierr = VecDestroy(&xdot);CHKERRQ(ierr);
1394036925cSBarry Smith     ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr);
140b51f60e1SBarry Smith     N = nits+2;
1414036925cSBarry Smith 
1424036925cSBarry Smith     if (nits) {
143a174af7bSBarry Smith       PetscDraw draw;
144*d98bdc7eSBarry Smith       PetscReal pause;
145*d98bdc7eSBarry Smith 
146b51f60e1SBarry Smith       ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr);
1470d18c744SBarry Smith       ierr = PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);CHKERRQ(ierr);
148a174af7bSBarry Smith       ierr = PetscMalloc2(PetscMax(n,N),PetscReal,&r,PetscMax(n,N),PetscReal,&c);CHKERRQ(ierr);
149b51f60e1SBarry Smith       if (ctx->computeexplicitly) {
150b51f60e1SBarry Smith         ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr);
151b51f60e1SBarry Smith         neig = n;
152b51f60e1SBarry Smith       } else {
153b51f60e1SBarry Smith         ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr);
154b51f60e1SBarry Smith       }
1554036925cSBarry Smith       for (i=0; i<neig; i++) {
1567d0d5eacSSatish Balay         r[i] -= PetscAbsScalar(ctx->shift);
157b51f60e1SBarry Smith         ierr = PetscPrintf(ctx->comm,"%g + %g i\n",r[i],c[i]);CHKERRQ(ierr);
1584036925cSBarry Smith         ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr);
1594036925cSBarry Smith       }
1600d18c744SBarry Smith       ierr = PetscFree2(r,c);CHKERRQ(ierr);
161a174af7bSBarry Smith       ierr = PetscDrawSPGetDraw(drawsp,&draw);CHKERRQ(ierr);
162*d98bdc7eSBarry Smith       ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr);
163*d98bdc7eSBarry Smith       ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr);
164*d98bdc7eSBarry Smith       ierr = PetscDrawSPDraw(drawsp,PETSC_TRUE);CHKERRQ(ierr);
165*d98bdc7eSBarry Smith       ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr);
1660d18c744SBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1670d18c744SBarry Smith       {
1680d18c744SBarry Smith       PetscDrawAxis axis;
1690d18c744SBarry Smith       PetscReal     xmin,xmax,ymin,ymax;
1700d18c744SBarry Smith 
1710d18c744SBarry Smith       ierr = PetscDrawSPGetAxis(drawsp,&axis);CHKERRQ(ierr);
1720d18c744SBarry Smith       ierr = PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);CHKERRQ(ierr);
1730d18c744SBarry Smith       ierr = PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,TSIndicator_Theta);CHKERRQ(ierr);
174a174af7bSBarry Smith       ierr = PetscDrawSPDraw(drawsp,PETSC_FALSE);CHKERRQ(ierr);
1750d18c744SBarry Smith       }
1760d18c744SBarry Smith #endif
1774036925cSBarry Smith     }
178a174af7bSBarry Smith     ierr = MatDestroy(&B);CHKERRQ(ierr);
1794036925cSBarry Smith   }
1804036925cSBarry Smith   PetscFunctionReturn(0);
1814036925cSBarry Smith }
1824036925cSBarry Smith 
1834036925cSBarry Smith #undef __FUNCT__
1844036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy"
1854036925cSBarry Smith /*@C
1864036925cSBarry Smith    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate().
1874036925cSBarry Smith 
1884036925cSBarry Smith    Collective on TSMonitorSPEigCtx
1894036925cSBarry Smith 
1904036925cSBarry Smith    Input Parameter:
1914036925cSBarry Smith .  ctx - the monitor context
1924036925cSBarry Smith 
1934036925cSBarry Smith    Level: intermediate
1944036925cSBarry Smith 
1954036925cSBarry Smith .keywords: TS, monitor, line graph, destroy
1964036925cSBarry Smith 
1974036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(),  TSMonitorSet(), TSMonitorSPEig();
1984036925cSBarry Smith @*/
1994036925cSBarry Smith PetscErrorCode  TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
2004036925cSBarry Smith {
2014036925cSBarry Smith   PetscDraw      draw;
2024036925cSBarry Smith   PetscErrorCode ierr;
2034036925cSBarry Smith 
2044036925cSBarry Smith   PetscFunctionBegin;
2054036925cSBarry Smith   ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr);
2064036925cSBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2074036925cSBarry Smith   ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr);
2084036925cSBarry Smith   ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr);
2094036925cSBarry Smith   ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr);
2104036925cSBarry Smith   ierr = PetscFree(*ctx);CHKERRQ(ierr);
2114036925cSBarry Smith   PetscFunctionReturn(0);
2124036925cSBarry Smith }
2134036925cSBarry Smith 
2140d18c744SBarry Smith 
2150d18c744SBarry Smith 
216