xref: /petsc/src/ts/interface/tseig.c (revision b51f60e1c60b87fd2d76e6847c03218ad92bfe04)
14036925cSBarry Smith 
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 */
9*b51f60e1SBarry Smith   PetscBool   computeexplicitly;
10*b51f60e1SBarry Smith   MPI_Comm    comm;
11*b51f60e1SBarry Smith   PetscRandom rand;
12*b51f60e1SBarry Smith   PetscScalar shift;
134036925cSBarry Smith };
144036925cSBarry Smith 
154036925cSBarry Smith 
164036925cSBarry Smith #undef __FUNCT__
174036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxCreate"
184036925cSBarry Smith /*@C
194036925cSBarry Smith    TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator
204036925cSBarry Smith 
214036925cSBarry Smith    Collective on TS
224036925cSBarry Smith 
234036925cSBarry Smith    Input Parameters:
244036925cSBarry Smith +  host - the X display to open, or null for the local machine
254036925cSBarry Smith .  label - the title to put in the title bar
264036925cSBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
274036925cSBarry Smith .  m, n - the screen width and height in pixels
284036925cSBarry Smith -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
294036925cSBarry Smith 
304036925cSBarry Smith    Output Parameter:
314036925cSBarry Smith .  ctx - the context
324036925cSBarry Smith 
334036925cSBarry Smith    Options Database Key:
34*b51f60e1SBarry Smith .  -ts_monitor_sp_eig - plot egienvalues of linearized right hand side
354036925cSBarry Smith 
364036925cSBarry Smith    Notes:
374036925cSBarry Smith    Use TSMonitorSPEigCtxDestroy() to destroy.
384036925cSBarry Smith 
39*b51f60e1SBarry Smith    Currently only works if the Jacobian is provided explicitly.
40*b51f60e1SBarry Smith 
41*b51f60e1SBarry Smith    Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix.
42*b51f60e1SBarry Smith 
434036925cSBarry Smith    Level: intermediate
444036925cSBarry Smith 
454036925cSBarry Smith .keywords: TS, monitor, line graph, residual, seealso
464036925cSBarry Smith 
474036925cSBarry Smith .seealso: TSMonitorSPEigTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
484036925cSBarry Smith 
494036925cSBarry Smith @*/
504036925cSBarry Smith PetscErrorCode  TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx *ctx)
514036925cSBarry Smith {
524036925cSBarry Smith   PetscDraw      win;
534036925cSBarry Smith   PetscErrorCode ierr;
544036925cSBarry Smith   PC             pc;
554036925cSBarry Smith 
564036925cSBarry Smith   PetscFunctionBegin;
574036925cSBarry Smith   ierr = PetscNew(struct _n_TSMonitorSPEigCtx,ctx);CHKERRQ(ierr);
58*b51f60e1SBarry Smith   ierr = PetscRandomCreate(comm,&(*ctx)->rand);CHKERRQ(ierr);
59*b51f60e1SBarry Smith   ierr = PetscRandomSetFromOptions((*ctx)->rand);CHKERRQ(ierr);
604036925cSBarry Smith   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
614036925cSBarry Smith   ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr);
624036925cSBarry Smith   ierr = PetscDrawSPCreate(win,1,&(*ctx)->drawsp);CHKERRQ(ierr);
634036925cSBarry Smith   ierr = KSPCreate(comm,&(*ctx)->ksp);CHKERRQ(ierr);
64*b51f60e1SBarry Smith   ierr = KSPSetOptionsPrefix((*ctx)->ksp,"ts_monitor_sp_eig_");CHKERRQ(ierr); /* this is wrong, used use also prefix from the TS */
654036925cSBarry Smith   ierr = KSPSetType((*ctx)->ksp,KSPGMRES);CHKERRQ(ierr);
664036925cSBarry Smith   ierr = KSPGMRESSetRestart((*ctx)->ksp,200);CHKERRQ(ierr);
674036925cSBarry Smith   ierr = KSPSetTolerances((*ctx)->ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,200);CHKERRQ(ierr);
684036925cSBarry Smith   ierr = KSPSetComputeSingularValues((*ctx)->ksp,PETSC_TRUE);CHKERRQ(ierr);
694036925cSBarry Smith   ierr = KSPSetFromOptions((*ctx)->ksp);CHKERRQ(ierr);
704036925cSBarry Smith   ierr = KSPGetPC((*ctx)->ksp,&pc);CHKERRQ(ierr);
714036925cSBarry Smith   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
724036925cSBarry Smith   (*ctx)->howoften          = howoften;
73*b51f60e1SBarry Smith   (*ctx)->computeexplicitly = PETSC_FALSE;
74*b51f60e1SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-ts_monitor_sp_eig_explicitly",&(*ctx)->computeexplicitly,PETSC_NULL);CHKERRQ(ierr);
75*b51f60e1SBarry Smith   (*ctx)->shift             = 0.0;
76*b51f60e1SBarry Smith   ierr = PetscOptionsGetScalar(PETSC_NULL,"-ts_monitor_sp_eig_shift",&(*ctx)->shift,PETSC_NULL);CHKERRQ(ierr);
77*b51f60e1SBarry Smith   (*ctx)->comm              = comm;
784036925cSBarry Smith   PetscFunctionReturn(0);
794036925cSBarry Smith }
804036925cSBarry Smith 
814036925cSBarry Smith #undef __FUNCT__
824036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEig"
834036925cSBarry Smith PetscErrorCode TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
844036925cSBarry Smith {
854036925cSBarry Smith   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx) monctx;
864036925cSBarry Smith   PetscErrorCode    ierr;
874036925cSBarry Smith   KSP               ksp = ctx->ksp;
88*b51f60e1SBarry Smith   PetscInt          n,N,nits,neig,i,its = 200;
894036925cSBarry Smith   PetscReal         *r,*c;
904036925cSBarry Smith   PetscDrawSP       drawsp = ctx->drawsp;
914036925cSBarry Smith   MatStructure      structure;
924036925cSBarry Smith   Mat               A,B;
934036925cSBarry Smith   Vec               xdot;
944036925cSBarry Smith   SNES              snes;
954036925cSBarry Smith 
964036925cSBarry Smith   PetscFunctionBegin;
974036925cSBarry Smith   if (!step) PetscFunctionReturn(0);
984036925cSBarry Smith   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && (step == -1))){
994036925cSBarry Smith     ierr = VecDuplicate(v,&xdot);CHKERRQ(ierr);
1004036925cSBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1014036925cSBarry Smith     ierr = SNESGetJacobian(snes,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
102*b51f60e1SBarry Smith     /*
103*b51f60e1SBarry Smith        This doesn't work because methods keep and use internal information about the shift so it
104*b51f60e1SBarry Smith        seems we would need code for each method to trick the correct Jacobian in being computed.
105*b51f60e1SBarry Smith     ts->time_step = PETSC_MAX_REAL;
1063d39420eSBarry Smith     ierr = SNESComputeJacobian(snes,v,&A,&B,&structure);CHKERRQ(ierr);
107*b51f60e1SBarry Smith     */
108*b51f60e1SBarry Smith     ierr = TSComputeIJacobian(ts,ptime,v,xdot,0.0,&A,&B,&structure,PETSC_FALSE);CHKERRQ(ierr);
1093d39420eSBarry Smith     ierr = MatScale(A,-1.0);CHKERRQ(ierr);
110*b51f60e1SBarry Smith     ierr = MatShift(A,ctx->shift);CHKERRQ(ierr);
111*b51f60e1SBarry Smith     ierr = KSPSetOperators(ksp,A,A,structure);CHKERRQ(ierr);
112*b51f60e1SBarry Smith     ierr = VecGetSize(v,&n);CHKERRQ(ierr);
113*b51f60e1SBarry Smith     if (n < 200) its = n;
114*b51f60e1SBarry Smith     ierr = KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,its);CHKERRQ(ierr);
1154036925cSBarry Smith     ierr = VecSetRandom(xdot,ctx->rand);CHKERRQ(ierr);
1164036925cSBarry Smith     ierr = KSPSolve(ksp,xdot,xdot);CHKERRQ(ierr);
1174036925cSBarry Smith     ierr = VecDestroy(&xdot);CHKERRQ(ierr);
1184036925cSBarry Smith     ierr = KSPGetIterationNumber(ksp,&nits);CHKERRQ(ierr);
119*b51f60e1SBarry Smith     N = nits+2;
1204036925cSBarry Smith 
1214036925cSBarry Smith     if (nits) {
122*b51f60e1SBarry Smith       ierr = PetscDrawSPReset(drawsp);CHKERRQ(ierr);
123*b51f60e1SBarry Smith       ierr = PetscMalloc2(N,PetscReal,&r,N,PetscReal,&c);CHKERRQ(ierr);
124*b51f60e1SBarry Smith       if (ctx->computeexplicitly) {
125*b51f60e1SBarry Smith         ierr = KSPComputeEigenvaluesExplicitly(ksp,n,r,c);CHKERRQ(ierr);
126*b51f60e1SBarry Smith         neig = n;
127*b51f60e1SBarry Smith       } else {
128*b51f60e1SBarry Smith         ierr = KSPComputeEigenvalues(ksp,N,r,c,&neig);CHKERRQ(ierr);
129*b51f60e1SBarry Smith       }
1304036925cSBarry Smith       for (i=0; i<neig; i++) {
131*b51f60e1SBarry Smith         r[i] -= ctx->shift;
132*b51f60e1SBarry Smith         ierr = PetscPrintf(ctx->comm,"%g + %g i\n",r[i],c[i]);CHKERRQ(ierr);
1334036925cSBarry Smith         ierr = PetscDrawSPAddPoint(drawsp,r+i,c+i);CHKERRQ(ierr);
1344036925cSBarry Smith       }
1354036925cSBarry Smith       ierr = PetscDrawSPDraw(drawsp);CHKERRQ(ierr);
1364036925cSBarry Smith       ierr = PetscFree2(r,c);CHKERRQ(ierr);
1374036925cSBarry Smith     }
138*b51f60e1SBarry Smith     ierr = MatShift(A,-ctx->shift);CHKERRQ(ierr);
139*b51f60e1SBarry Smith     ierr = MatScale(A,-1.0);CHKERRQ(ierr);
1404036925cSBarry Smith   }
1414036925cSBarry Smith   PetscFunctionReturn(0);
1424036925cSBarry Smith }
1434036925cSBarry Smith 
1444036925cSBarry Smith #undef __FUNCT__
1454036925cSBarry Smith #define __FUNCT__ "TSMonitorSPEigCtxDestroy"
1464036925cSBarry Smith /*@C
1474036925cSBarry Smith    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate().
1484036925cSBarry Smith 
1494036925cSBarry Smith    Collective on TSMonitorSPEigCtx
1504036925cSBarry Smith 
1514036925cSBarry Smith    Input Parameter:
1524036925cSBarry Smith .  ctx - the monitor context
1534036925cSBarry Smith 
1544036925cSBarry Smith    Level: intermediate
1554036925cSBarry Smith 
1564036925cSBarry Smith .keywords: TS, monitor, line graph, destroy
1574036925cSBarry Smith 
1584036925cSBarry Smith .seealso: TSMonitorSPEigCtxCreate(),  TSMonitorSet(), TSMonitorSPEig();
1594036925cSBarry Smith @*/
1604036925cSBarry Smith PetscErrorCode  TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
1614036925cSBarry Smith {
1624036925cSBarry Smith   PetscDraw      draw;
1634036925cSBarry Smith   PetscErrorCode ierr;
1644036925cSBarry Smith 
1654036925cSBarry Smith   PetscFunctionBegin;
1664036925cSBarry Smith   ierr = PetscDrawSPGetDraw((*ctx)->drawsp,&draw);CHKERRQ(ierr);
1674036925cSBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1684036925cSBarry Smith   ierr = PetscDrawSPDestroy(&(*ctx)->drawsp);CHKERRQ(ierr);
1694036925cSBarry Smith   ierr = KSPDestroy(&(*ctx)->ksp);CHKERRQ(ierr);
1704036925cSBarry Smith   ierr = PetscRandomDestroy(&(*ctx)->rand);CHKERRQ(ierr);
1714036925cSBarry Smith   ierr = PetscFree(*ctx);CHKERRQ(ierr);
1724036925cSBarry Smith   PetscFunctionReturn(0);
1734036925cSBarry Smith }
1744036925cSBarry Smith 
175