Actual source code: tseig.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
  3: #include <petscdraw.h>

  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: };


 19: /*@C
 20:    TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator

 22:    Collective on TS

 24:    Input Parameters:
 25: +  host - the X display to open, or null for the local machine
 26: .  label - the title to put in the title bar
 27: .  x, y - the screen coordinates of the upper left coordinate of the window
 28: .  m, n - the screen width and height in pixels
 29: -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time

 31:    Output Parameter:
 32: .  ctx - the context

 34:    Options Database Key:
 35: .  -ts_monitor_sp_eig - plot egienvalues of linearized right hand side

 37:    Notes:
 38:    Use TSMonitorSPEigCtxDestroy() to destroy.

 40:    Currently only works if the Jacobian is provided explicitly.

 42:    Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix.

 44:    Level: intermediate

 46: .keywords: TS, monitor, line graph, residual, seealso

 48: .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()

 50: @*/
 51: PetscErrorCode  TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx)
 52: {
 53:   PetscDraw      win;
 55:   PC             pc;

 58:   PetscNew(ctx);
 59:   PetscRandomCreate(comm,&(*ctx)->rand);
 60:   PetscRandomSetFromOptions((*ctx)->rand);
 61:   PetscDrawCreate(comm,host,label,x,y,m,n,&win);
 62:   PetscDrawSetFromOptions(win);
 63:   PetscDrawSPCreate(win,1,&(*ctx)->drawsp);
 64:   KSPCreate(comm,&(*ctx)->ksp);
 65:   KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_sp_eig_"); /* this is wrong, used use also prefix from the TS */
 66:   KSPSetType((*ctx)->ksp,KSPGMRES);
 67:   KSPGMRESSetRestart((*ctx)->ksp,200);
 68:   KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);
 69:   KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);
 70:   KSPSetFromOptions((*ctx)->ksp);
 71:   KSPGetPC((*ctx)->ksp,&pc);
 72:   PCSetType(pc,PCNONE);

 74:   (*ctx)->howoften          = howoften;
 75:   (*ctx)->computeexplicitly = PETSC_FALSE;

 77:   PetscOptionsGetBool(NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,NULL);

 79:   (*ctx)->comm = comm;
 80:   (*ctx)->xmin = -2.1;
 81:   (*ctx)->xmax = 1.1;
 82:   (*ctx)->ymin = -1.1;
 83:   (*ctx)->ymax = 1.1;
 84:   return(0);
 85: }

 89: static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr,PetscReal xi,PetscBool *flg)
 90: {
 92:   PetscReal      yr,yi;

 95:   TSComputeLinearStability(ts,xr,xi,&yr,&yi);
 96:   if ((yr*yr + yi*yi) <= 1.0) *flg = PETSC_TRUE;
 97:   else *flg = PETSC_FALSE;
 98:   return(0);
 99: }

103: PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
104: {
105:   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx;
106:   PetscErrorCode    ierr;
107:   KSP               ksp = ctx->ksp;
108:   PetscInt          n,N,nits,neig,i,its = 200;
109:   PetscReal         *r,*c,time_step_save;
110:   PetscDrawSP       drawsp = ctx->drawsp;
111:   Mat               A,B;
112:   Vec               xdot;
113:   SNES              snes;

116:   if (!step) return(0);
117:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
118:     VecDuplicate(v,&xdot);
119:     TSGetSNES(ts,&snes);
120:     SNESGetJacobian(snes,&A,&B,NULL,NULL);
121:     MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);
122:     /*
123:        This doesn't work because methods keep and use internal information about the shift so it
124:        seems we would need code for each method to trick the correct Jacobian in being computed.
125:      */
126:     time_step_save = ts->time_step;
127:     ts->time_step  = PETSC_MAX_REAL;

129:     SNESComputeJacobian(snes,v,A,B);

131:     ts->time_step  = time_step_save;

133:     KSPSetOperators(ksp,B,B);
134:     VecGetSize(v,&n);
135:     if (n < 200) its = n;
136:     KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);
137:     VecSetRandom(xdot,ctx->rand);
138:     KSPSolve(ksp,xdot,xdot);
139:     VecDestroy(&xdot);
140:     KSPGetIterationNumber(ksp,&nits);
141:     N    = nits+2;

143:     if (nits) {
144:       PetscDraw     draw;
145:       PetscReal     pause;
146:       PetscDrawAxis axis;
147:       PetscReal     xmin,xmax,ymin,ymax;

149:       PetscDrawSPReset(drawsp);
150:       PetscDrawSPSetLimits(drawsp,ctx->xmin,ctx->xmax,ctx->ymin,ctx->ymax);
151:       PetscMalloc2(PetscMax(n,N),&r,PetscMax(n,N),&c);
152:       if (ctx->computeexplicitly) {
153:         KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
154:         neig = n;
155:       } else {
156:         KSPComputeEigenvalues(ksp,N,r,c,&neig);
157:       }
158:       /* We used the positive operator to be able to reuse KSPs that require positive definiteness, now flip the spectrum as is conventional for ODEs */
159:       for (i=0; i<neig; i++) r[i] = -r[i];
160:       for (i=0; i<neig; i++) {
161:         if (ts->ops->linearstability) {
162:           PetscReal fr,fi;
163:           TSComputeLinearStability(ts,r[i],c[i],&fr,&fi);
164:           if ((fr*fr + fi*fi) > 1.0) {
165:             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));
166:           }
167:         }
168:         PetscDrawSPAddPoint(drawsp,r+i,c+i);
169:       }
170:       PetscFree2(r,c);
171:       PetscDrawSPGetDraw(drawsp,&draw);
172:       PetscDrawGetPause(draw,&pause);
173:       PetscDrawSetPause(draw,0.0);
174:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
175:       PetscDrawSetPause(draw,pause);

177:       if (ts->ops->linearstability) {
178:         PetscDrawSPGetAxis(drawsp,&axis);
179:         PetscDrawAxisGetLimits(axis,&xmin,&xmax,&ymin,&ymax);
180:         PetscDrawIndicatorFunction(draw,xmin,xmax,ymin,ymax,PETSC_DRAW_CYAN,(PetscErrorCode (*)(void*,PetscReal,PetscReal,PetscBool*))TSLinearStabilityIndicator,ts);
181:         PetscDrawSPDraw(drawsp,PETSC_FALSE);
182:       }
183:     }
184:     MatDestroy(&B);
185:   }
186:   return(0);
187: }

191: /*@C
192:    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate().

194:    Collective on TSMonitorSPEigCtx

196:    Input Parameter:
197: .  ctx - the monitor context

199:    Level: intermediate

201: .keywords: TS, monitor, line graph, destroy

203: .seealso: TSMonitorSPEigCtxCreate(),  TSMonitorSet(), TSMonitorSPEig();
204: @*/
205: PetscErrorCode  TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
206: {
207:   PetscDraw      draw;

211:   PetscDrawSPGetDraw((*ctx)->drawsp,&draw);
212:   PetscDrawDestroy(&draw);
213:   PetscDrawSPDestroy(&(*ctx)->drawsp);
214:   KSPDestroy(&(*ctx)->ksp);
215:   PetscRandomDestroy(&(*ctx)->rand);
216:   PetscFree(*ctx);
217:   return(0);
218: }