xref: /petsc/src/ts/interface/tseig.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 
2 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
3 #include <petscdraw.h>
4 
5 /* ------------------------------------------------------------------------*/
6 struct _n_TSMonitorSPEigCtx {
7   PetscDrawSP drawsp;
8   KSP         ksp;
9   PetscInt    howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
10   PetscBool   computeexplicitly;
11   MPI_Comm    comm;
12   PetscRandom rand;
13   PetscReal   xmin,xmax,ymin,ymax;
14 };
15 
16 
17 /*@C
18    TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator
19 
20    Collective on TS
21 
22    Input Parameters:
23 +  host - the X display to open, or null for the local machine
24 .  label - the title to put in the title bar
25 .  x, y - the screen coordinates of the upper left coordinate of the window
26 .  m, n - the screen width and height in pixels
27 -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
28 
29    Output Parameter:
30 .  ctx - the context
31 
32    Options Database Key:
33 .  -ts_monitor_sp_eig - plot egienvalues of linearized right hand side
34 
35    Notes:
36    Use TSMonitorSPEigCtxDestroy() to destroy.
37 
38    Currently only works if the Jacobian is provided explicitly.
39 
40    Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix.
41 
42    Level: intermediate
43 
44 .keywords: TS, monitor, line graph, residual, seealso
45 
46 .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
47 
48 @*/
49 PetscErrorCode  TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx)
50 {
51   PetscDraw      win;
52   PetscErrorCode ierr;
53   PC             pc;
54 
55   PetscFunctionBegin;
56   ierr = PetscNew(ctx);CHKERRQ(ierr);
57   ierr = PetscRandomCreate(comm,&(*ctx)->rand);CHKERRQ(ierr);
58   ierr = PetscRandomSetFromOptions((*ctx)->rand);CHKERRQ(ierr);
59   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
60   ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr);
61   ierr = PetscDrawSPCreate(win,1,&(*ctx)->drawsp);CHKERRQ(ierr);
62   ierr = KSPCreate(comm,&(*ctx)->ksp);CHKERRQ(ierr);
63   ierr = KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_sp_eig_");CHKERRQ(ierr); /* this is wrong, used use also prefix from the TS */
64   ierr = KSPSetType((*ctx)->ksp,KSPGMRES);CHKERRQ(ierr);
65   ierr = KSPGMRESSetRestart((*ctx)->ksp,200);CHKERRQ(ierr);
66   ierr = KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);CHKERRQ(ierr);
67   ierr = KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);CHKERRQ(ierr);
68   ierr = KSPSetFromOptions((*ctx)->ksp);CHKERRQ(ierr);
69   ierr = KSPGetPC((*ctx)->ksp,&pc);CHKERRQ(ierr);
70   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
71 
72   (*ctx)->howoften          = howoften;
73   (*ctx)->computeexplicitly = PETSC_FALSE;
74 
75   ierr = PetscOptionsGetBool(NULL,NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,NULL);CHKERRQ(ierr);
76 
77   (*ctx)->comm = comm;
78   (*ctx)->xmin = -2.1;
79   (*ctx)->xmax = 1.1;
80   (*ctx)->ymin = -1.1;
81   (*ctx)->ymax = 1.1;
82   PetscFunctionReturn(0);
83 }
84 
85 static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr,PetscReal xi,PetscBool *flg)
86 {
87   PetscErrorCode ierr;
88   PetscReal      yr,yi;
89 
90   PetscFunctionBegin;
91   ierr = TSComputeLinearStability(ts,xr,xi,&yr,&yi);CHKERRQ(ierr);
92   if ((yr*yr + yi*yi) <= 1.0) *flg = PETSC_TRUE;
93   else *flg = PETSC_FALSE;
94   PetscFunctionReturn(0);
95 }
96 
97 PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
98 {
99   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx;
100   PetscErrorCode    ierr;
101   KSP               ksp = ctx->ksp;
102   PetscInt          n,N,nits,neig,i,its = 200;
103   PetscReal         *r,*c,time_step_save;
104   PetscDrawSP       drawsp = ctx->drawsp;
105   Mat               A,B;
106   Vec               xdot;
107   SNES              snes;
108 
109   PetscFunctionBegin;
110   if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */
111   if (!step) PetscFunctionReturn(0);
112   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
113     ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr);
114     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
115     ierr = SNESGetJacobian(snes,&A,&B,NULL,NULL);CHKERRQ(ierr);
116     ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);CHKERRQ(ierr);
117     /*
118        This doesn't work because methods keep and use internal information about the shift so it
119        seems we would need code for each method to trick the correct Jacobian in being computed.
120      */
121     time_step_save = ts->time_step;
122     ts->time_step  = PETSC_MAX_REAL;
123 
124     ierr = SNESComputeJacobian(snes,v,A,B);CHKERRQ(ierr);
125 
126     ts->time_step  = time_step_save;
127 
128     ierr = KSPSetOperators(ksp,B,B);CHKERRQ(ierr);
129     ierr = VecGetSize(v,&n);CHKERRQ(ierr);
130     if (n < 200) its = n;
131     ierr = KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);CHKERRQ(ierr);
132     ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr);
133     ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr);
134     ierr = VecDestroy(&xdot);CHKERRQ(ierr);
135     ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr);
136     N    = nits+2;
137 
138     if (nits) {
139       PetscDraw     draw;
140       PetscReal     pause;
141       PetscDrawAxis axis;
142       PetscReal     xmin,xmax,ymin,ymax;
143 
144       ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr);
145       ierr = PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);CHKERRQ(ierr);
146       ierr = PetscMalloc2(PetscMax(n,N),&r,PetscMax(n,N),&c);CHKERRQ(ierr);
147       if (ctx->computeexplicitly) {
148         ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr);
149         neig = n;
150       } else {
151         ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr);
152       }
153       /* 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       for (i=0; i<neig; i++) r[i] = -r[i];
155       for (i=0; i<neig; i++) {
156         if (ts->ops->linearstability) {
157           PetscReal fr,fi;
158           ierr = TSComputeLinearStability(ts,r[i],c[i],&fr,&fi);CHKERRQ(ierr);
159           if ((fr*fr + fi*fi) > 1.0) {
160             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);
161           }
162         }
163         ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr);
164       }
165       ierr = PetscFree2(r,c);CHKERRQ(ierr);
166       ierr = PetscDrawSPGetDraw(drawsp,&draw);CHKERRQ(ierr);
167       ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr);
168       ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr);
169       ierr = PetscDrawSPDraw(drawsp,PETSC_TRUE);CHKERRQ(ierr);
170       ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr);
171       if (ts->ops->linearstability) {
172         ierr = PetscDrawSPGetAxis(drawsp,&axis);CHKERRQ(ierr);
173         ierr = PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);CHKERRQ(ierr);
174         ierr = PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,(PetscErrorCode (*)(void*,PetscReal,PetscReal,PetscBool*))TSLinearStabilityIndicator,ts);CHKERRQ(ierr);
175         ierr = PetscDrawSPDraw(drawsp,PETSC_FALSE);CHKERRQ(ierr);
176       }
177       ierr = PetscDrawSPSave(drawsp);CHKERRQ(ierr);
178     }
179     ierr = MatDestroy(&B);CHKERRQ(ierr);
180   }
181   PetscFunctionReturn(0);
182 }
183 
184 /*@C
185    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate().
186 
187    Collective on TSMonitorSPEigCtx
188 
189    Input Parameter:
190 .  ctx - the monitor context
191 
192    Level: intermediate
193 
194 .keywords: TS, monitor, line graph, destroy
195 
196 .seealso: TSMonitorSPEigCtxCreate(),  TSMonitorSet(), TSMonitorSPEig();
197 @*/
198 PetscErrorCode  TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
199 {
200   PetscDraw      draw;
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr);
205   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
206   ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr);
207   ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr);
208   ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr);
209   ierr = PetscFree(*ctx);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 
214 
215