xref: /petsc/src/ts/interface/tseig.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
29804daf3SBarry Smith #include <petscdraw.h>
34036925cSBarry Smith 
44036925cSBarry Smith struct _n_TSMonitorSPEigCtx {
54036925cSBarry Smith   PetscDrawSP drawsp;
64036925cSBarry Smith   KSP         ksp;
74036925cSBarry Smith   PetscInt    howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
8b51f60e1SBarry Smith   PetscBool   computeexplicitly;
9b51f60e1SBarry Smith   MPI_Comm    comm;
10b51f60e1SBarry Smith   PetscRandom rand;
110d18c744SBarry Smith   PetscReal   xmin, xmax, ymin, ymax;
124036925cSBarry Smith };
134036925cSBarry Smith 
144036925cSBarry Smith /*@C
15bcf0153eSBarry Smith   TSMonitorSPEigCtxCreate - Creates a context for use with `TS` to monitor the eigenvalues of the linearized operator
164036925cSBarry Smith 
17bcf0153eSBarry Smith   Collective
184036925cSBarry Smith 
194036925cSBarry Smith   Input Parameters:
202fe279fdSBarry Smith + comm     - the communicator to share the monitor
212fe279fdSBarry Smith . host     - the X display to open, or `NULL` for the local machine
224036925cSBarry Smith . label    - the title to put in the title bar
232fe279fdSBarry Smith . x        - the horizontal screen coordinates of the upper left coordinate of the window
242fe279fdSBarry Smith . y        - the vertical coordinates of the upper left coordinate of the window
252fe279fdSBarry Smith . m        - the screen width in pixels
262fe279fdSBarry Smith . n        - the screen height in pixels
274036925cSBarry Smith - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
284036925cSBarry Smith 
294036925cSBarry Smith   Output Parameter:
304036925cSBarry Smith . ctx - the context
314036925cSBarry Smith 
324036925cSBarry Smith   Options Database Key:
33dd8e379bSPierre Jolivet . -ts_monitor_sp_eig - plot egienvalues of linearized right-hand side
344036925cSBarry Smith 
35bcf0153eSBarry Smith   Level: intermediate
36bcf0153eSBarry Smith 
374036925cSBarry Smith   Notes:
38bcf0153eSBarry Smith   Use `TSMonitorSPEigCtxDestroy()` to destroy the context
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 
441cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
454036925cSBarry Smith @*/
TSMonitorSPEigCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPEigCtx * ctx)46d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorSPEigCtx *ctx)
47d71ae5a4SJacob Faibussowitsch {
484036925cSBarry Smith   PetscDraw win;
494036925cSBarry Smith   PC        pc;
504036925cSBarry Smith 
514036925cSBarry Smith   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
539566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &(*ctx)->rand));
549566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions((*ctx)->rand));
559566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &win));
569566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(win));
579566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(win, 1, &(*ctx)->drawsp));
589566063dSJacob Faibussowitsch   PetscCall(KSPCreate(comm, &(*ctx)->ksp));
599566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix((*ctx)->ksp, "ts_monitor_sp_eig_")); /* this is wrong, used use also prefix from the TS */
609566063dSJacob Faibussowitsch   PetscCall(KSPSetType((*ctx)->ksp, KSPGMRES));
619566063dSJacob Faibussowitsch   PetscCall(KSPGMRESSetRestart((*ctx)->ksp, 200));
62*09cb0f53SBarry Smith   PetscCall(KSPSetTolerances((*ctx)->ksp, 1.e-10, PETSC_CURRENT, PETSC_CURRENT, 200));
639566063dSJacob Faibussowitsch   PetscCall(KSPSetComputeSingularValues((*ctx)->ksp, PETSC_TRUE));
649566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions((*ctx)->ksp));
659566063dSJacob Faibussowitsch   PetscCall(KSPGetPC((*ctx)->ksp, &pc));
669566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
67bbd56ea5SKarl Rupp 
684036925cSBarry Smith   (*ctx)->howoften          = howoften;
69b51f60e1SBarry Smith   (*ctx)->computeexplicitly = PETSC_FALSE;
70bbd56ea5SKarl Rupp 
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_sp_eig_explicitly", &(*ctx)->computeexplicitly, NULL));
72bbd56ea5SKarl Rupp 
73b51f60e1SBarry Smith   (*ctx)->comm = comm;
740d18c744SBarry Smith   (*ctx)->xmin = -2.1;
75f9c1d6abSBarry Smith   (*ctx)->xmax = 1.1;
760d18c744SBarry Smith   (*ctx)->ymin = -1.1;
770d18c744SBarry Smith   (*ctx)->ymax = 1.1;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
794036925cSBarry Smith }
804036925cSBarry Smith 
TSLinearStabilityIndicator(TS ts,PetscReal xr,PetscReal xi,PetscBool * flg)81d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr, PetscReal xi, PetscBool *flg)
82d71ae5a4SJacob Faibussowitsch {
83f9c1d6abSBarry Smith   PetscReal yr, yi;
840d18c744SBarry Smith 
85f9c1d6abSBarry Smith   PetscFunctionBegin;
869566063dSJacob Faibussowitsch   PetscCall(TSComputeLinearStability(ts, xr, xi, &yr, &yi));
87f9c1d6abSBarry Smith   if ((yr * yr + yi * yi) <= 1.0) *flg = PETSC_TRUE;
88f9c1d6abSBarry Smith   else *flg = PETSC_FALSE;
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
900d18c744SBarry Smith }
910d18c744SBarry Smith 
TSMonitorSPEig(TS ts,PetscInt step,PetscReal ptime,Vec v,void * monctx)92d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEig(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
93d71ae5a4SJacob Faibussowitsch {
944036925cSBarry Smith   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx)monctx;
954036925cSBarry Smith   KSP               ksp = ctx->ksp;
96b51f60e1SBarry Smith   PetscInt          n, N, nits, neig, i, its = 200;
97b296d7d5SJed Brown   PetscReal        *r, *c, time_step_save;
984036925cSBarry Smith   PetscDrawSP       drawsp = ctx->drawsp;
994036925cSBarry Smith   Mat               A, B;
1004036925cSBarry Smith   Vec               xdot;
1014036925cSBarry Smith   SNES              snes;
1024036925cSBarry Smith 
1034036925cSBarry Smith   PetscFunctionBegin;
1043ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
1053ba16761SJacob Faibussowitsch   if (!step) PetscFunctionReturn(PETSC_SUCCESS);
106b06615a5SLisandro Dalcin   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1079566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v, &xdot));
1089566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
1099566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &A, &B, NULL, NULL));
1109566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B));
111b51f60e1SBarry Smith     /*
112b51f60e1SBarry Smith        This doesn't work because methods keep and use internal information about the shift so it
113b51f60e1SBarry Smith        seems we would need code for each method to trick the correct Jacobian in being computed.
114b296d7d5SJed Brown      */
115b296d7d5SJed Brown     time_step_save = ts->time_step;
116b51f60e1SBarry Smith     ts->time_step  = PETSC_MAX_REAL;
117bbd56ea5SKarl Rupp 
1189566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, v, A, B));
119bbd56ea5SKarl Rupp 
120b296d7d5SJed Brown     ts->time_step = time_step_save;
121a174af7bSBarry Smith 
1229566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, B, B));
1239566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &n));
124b51f60e1SBarry Smith     if (n < 200) its = n;
125*09cb0f53SBarry Smith     PetscCall(KSPSetTolerances(ksp, 1.e-10, PETSC_CURRENT, PETSC_CURRENT, its));
1269566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xdot, ctx->rand));
1279566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ksp, xdot, xdot));
1289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xdot));
1299566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(ksp, &nits));
130b51f60e1SBarry Smith     N = nits + 2;
1314036925cSBarry Smith 
1324036925cSBarry Smith     if (nits) {
133a174af7bSBarry Smith       PetscDraw     draw;
134d98bdc7eSBarry Smith       PetscReal     pause;
135f9c1d6abSBarry Smith       PetscDrawAxis axis;
136f9c1d6abSBarry Smith       PetscReal     xmin, xmax, ymin, ymax;
137d98bdc7eSBarry Smith 
1389566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPReset(drawsp));
1399566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPSetLimits(drawsp, ctx->xmin, ctx->xmax, ctx->ymin, ctx->ymax));
1409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(PetscMax(n, N), &r, PetscMax(n, N), &c));
141b51f60e1SBarry Smith       if (ctx->computeexplicitly) {
1429566063dSJacob Faibussowitsch         PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
143b51f60e1SBarry Smith         neig = n;
144b51f60e1SBarry Smith       } else {
1459566063dSJacob Faibussowitsch         PetscCall(KSPComputeEigenvalues(ksp, N, r, c, &neig));
146b51f60e1SBarry Smith       }
147b296d7d5SJed Brown       /* We used the positive operator to be able to reuse KSPs that require positive definiteness, now flip the spectrum as is conventional for ODEs */
148b296d7d5SJed Brown       for (i = 0; i < neig; i++) r[i] = -r[i];
1494036925cSBarry Smith       for (i = 0; i < neig; i++) {
1505f737190SBarry Smith         if (ts->ops->linearstability) {
1515f737190SBarry Smith           PetscReal fr, fi;
1529566063dSJacob Faibussowitsch           PetscCall(TSComputeLinearStability(ts, r[i], c[i], &fr, &fi));
15348a46eb9SPierre Jolivet           if ((fr * fr + fi * fi) > 1.0) PetscCall(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)));
1545d726895SBarry Smith         }
1559566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
1564036925cSBarry Smith       }
1579566063dSJacob Faibussowitsch       PetscCall(PetscFree2(r, c));
1589566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPGetDraw(drawsp, &draw));
1599566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetPause(draw, &pause));
1609566063dSJacob Faibussowitsch       PetscCall(PetscDrawSetPause(draw, 0.0));
1619566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
1629566063dSJacob Faibussowitsch       PetscCall(PetscDrawSetPause(draw, pause));
163f9c1d6abSBarry Smith       if (ts->ops->linearstability) {
1649566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPGetAxis(drawsp, &axis));
1659566063dSJacob Faibussowitsch         PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, &ymin, &ymax));
1669566063dSJacob Faibussowitsch         PetscCall(PetscDrawIndicatorFunction(draw, xmin, xmax, ymin, ymax, PETSC_DRAW_CYAN, (PetscErrorCode (*)(void *, PetscReal, PetscReal, PetscBool *))TSLinearStabilityIndicator, ts));
1679566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPDraw(drawsp, PETSC_FALSE));
1680d18c744SBarry Smith       }
1699566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPSave(drawsp));
1704036925cSBarry Smith     }
1719566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
1724036925cSBarry Smith   }
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1744036925cSBarry Smith }
1754036925cSBarry Smith 
1764036925cSBarry Smith /*@C
177bcf0153eSBarry Smith   TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with `TSMonitorSPEigCtxCreate()`.
1784036925cSBarry Smith 
179c3339decSBarry Smith   Collective
1804036925cSBarry Smith 
1814036925cSBarry Smith   Input Parameter:
1824036925cSBarry Smith . ctx - the monitor context
1834036925cSBarry Smith 
1844036925cSBarry Smith   Level: intermediate
1854036925cSBarry Smith 
186bcf0153eSBarry Smith   Note:
187bcf0153eSBarry Smith   Should be passed to `TSMonitorSet()` along with `TSMonitorSPEig()` an the context created with `TSMonitorSPEigCtxCreate()`
188bcf0153eSBarry Smith 
189a54bb2a9SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigCtxCreate()`, `TSMonitorSet()`, `TSMonitorSPEig()`
1904036925cSBarry Smith @*/
TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx * ctx)191d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
192d71ae5a4SJacob Faibussowitsch {
1934036925cSBarry Smith   PetscDraw draw;
1944036925cSBarry Smith 
1954036925cSBarry Smith   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPGetDraw((*ctx)->drawsp, &draw));
1979566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
1989566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->drawsp));
1999566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&(*ctx)->ksp));
2009566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&(*ctx)->rand));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2034036925cSBarry Smith }
204