xref: /petsc/src/ts/interface/tseig.c (revision 5f737190a3790417cd725c1776731a3ea6717c8b)
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;
80f9c1d6abSBarry Smith   (*ctx)->xmax = 1.1;
810d18c744SBarry Smith   (*ctx)->ymin = -1.1;
820d18c744SBarry Smith   (*ctx)->ymax = 1.1;
834036925cSBarry Smith   PetscFunctionReturn(0);
844036925cSBarry Smith }
854036925cSBarry Smith 
86f9c1d6abSBarry Smith #undef __FUNCT__
87f9c1d6abSBarry Smith #define __FUNCT__ "TSLinearStabilityIndicator"
88f9c1d6abSBarry Smith static PetscErrorCode TSLinearStabilityIndicator(TS  ts, PetscReal xr,PetscReal xi,PetscBool *flg)
890d18c744SBarry Smith {
90f9c1d6abSBarry Smith   PetscErrorCode ierr;
91f9c1d6abSBarry Smith   PetscReal      yr,yi;
920d18c744SBarry Smith 
93f9c1d6abSBarry Smith   PetscFunctionBegin;
94f9c1d6abSBarry Smith   ierr = TSComputeLinearStability(ts,xr,xi,&yr,&yi);CHKERRQ(ierr);
95f9c1d6abSBarry Smith   if ((yr*yr + yi*yi) <= 1.0) *flg = PETSC_TRUE;
96f9c1d6abSBarry Smith   else *flg = PETSC_FALSE;
97f9c1d6abSBarry Smith   PetscFunctionReturn(0);
980d18c744SBarry Smith }
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;
144d98bdc7eSBarry Smith       PetscReal     pause;
145f9c1d6abSBarry Smith       PetscDrawAxis axis;
146f9c1d6abSBarry Smith       PetscReal     xmin,xmax,ymin,ymax;
147d98bdc7eSBarry Smith 
148b51f60e1SBarry Smith       ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr);
1490d18c744SBarry Smith       ierr = PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);CHKERRQ(ierr);
150a174af7bSBarry Smith       ierr = PetscMalloc2(PetscMax(n,N),PetscReal,&r,PetscMax(n,N),PetscReal,&c);CHKERRQ(ierr);
151b51f60e1SBarry Smith       if (ctx->computeexplicitly) {
152b51f60e1SBarry Smith         ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr);
153b51f60e1SBarry Smith         neig = n;
154b51f60e1SBarry Smith       } else {
155b51f60e1SBarry Smith         ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr);
156b51f60e1SBarry Smith       }
1574036925cSBarry Smith       for (i=0; i<neig; i++) {
1587d0d5eacSSatish Balay         r[i] -= PetscAbsScalar(ctx->shift);
159*5f737190SBarry Smith         if (ts->ops->linearstability) {
160*5f737190SBarry Smith           PetscReal fr,fi;
161*5f737190SBarry Smith           ierr = TSComputeLinearStability(ts,r[i],c[i],&fr,&fi);CHKERRQ(ierr);
162*5f737190SBarry Smith           if ((fr*fr + fi*fi) > 1.0) {
163*5f737190SBarry 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);
164*5f737190SBarry Smith           }
1655d726895SBarry Smith         }
1664036925cSBarry Smith         ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr);
1674036925cSBarry Smith       }
1680d18c744SBarry Smith       ierr = PetscFree2(r,c);CHKERRQ(ierr);
169a174af7bSBarry Smith       ierr = PetscDrawSPGetDraw(drawsp,&draw);CHKERRQ(ierr);
170d98bdc7eSBarry Smith       ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr);
171d98bdc7eSBarry Smith       ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr);
172d98bdc7eSBarry Smith       ierr = PetscDrawSPDraw(drawsp,PETSC_TRUE);CHKERRQ(ierr);
173d98bdc7eSBarry Smith       ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr);
1740d18c744SBarry Smith 
175f9c1d6abSBarry Smith       if (ts->ops->linearstability) {
1760d18c744SBarry Smith         ierr = PetscDrawSPGetAxis(drawsp,&axis);CHKERRQ(ierr);
1770d18c744SBarry Smith         ierr = PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);CHKERRQ(ierr);
178f9c1d6abSBarry Smith         ierr = PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,(PetscErrorCode (*)(void*,PetscReal,PetscReal,PetscBool*))TSLinearStabilityIndicator,ts);CHKERRQ(ierr);
179a174af7bSBarry Smith         ierr = PetscDrawSPDraw(drawsp,PETSC_FALSE);CHKERRQ(ierr);
1800d18c744SBarry Smith       }
1814036925cSBarry Smith     }
182a174af7bSBarry Smith     ierr = MatDestroy(&B);CHKERRQ(ierr);
1834036925cSBarry Smith   }
1844036925cSBarry Smith   PetscFunctionReturn(0);
1854036925cSBarry Smith }
1864036925cSBarry Smith 
1874036925cSBarry Smith #undef __FUNCT__
1884036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy"
1894036925cSBarry Smith /*@C
1904036925cSBarry Smith    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate().
1914036925cSBarry Smith 
1924036925cSBarry Smith    Collective on TSMonitorSPEigCtx
1934036925cSBarry Smith 
1944036925cSBarry Smith    Input Parameter:
1954036925cSBarry Smith .  ctx - the monitor context
1964036925cSBarry Smith 
1974036925cSBarry Smith    Level: intermediate
1984036925cSBarry Smith 
1994036925cSBarry Smith .keywords: TS, monitor, line graph, destroy
2004036925cSBarry Smith 
2014036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(),  TSMonitorSet(), TSMonitorSPEig();
2024036925cSBarry Smith @*/
2034036925cSBarry Smith PetscErrorCode  TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
2044036925cSBarry Smith {
2054036925cSBarry Smith   PetscDraw      draw;
2064036925cSBarry Smith   PetscErrorCode ierr;
2074036925cSBarry Smith 
2084036925cSBarry Smith   PetscFunctionBegin;
2094036925cSBarry Smith   ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr);
2104036925cSBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2114036925cSBarry Smith   ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr);
2124036925cSBarry Smith   ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr);
2134036925cSBarry Smith   ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr);
2144036925cSBarry Smith   ierr = PetscFree(*ctx);CHKERRQ(ierr);
2154036925cSBarry Smith   PetscFunctionReturn(0);
2164036925cSBarry Smith }
2174036925cSBarry Smith 
2180d18c744SBarry Smith 
2190d18c744SBarry Smith 
220