xref: /petsc/src/ts/interface/tseig.c (revision bcf0153e883cfed9568ef4557dcc209048fb58f7)
18cd53115SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
39804daf3SBarry Smith #include <petscdraw.h>
44036925cSBarry Smith 
54036925cSBarry Smith /* ------------------------------------------------------------------------*/
64036925cSBarry Smith struct _n_TSMonitorSPEigCtx {
74036925cSBarry Smith   PetscDrawSP drawsp;
84036925cSBarry Smith   KSP         ksp;
94036925cSBarry Smith   PetscInt    howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
10b51f60e1SBarry Smith   PetscBool   computeexplicitly;
11b51f60e1SBarry Smith   MPI_Comm    comm;
12b51f60e1SBarry Smith   PetscRandom rand;
130d18c744SBarry Smith   PetscReal   xmin, xmax, ymin, ymax;
144036925cSBarry Smith };
154036925cSBarry Smith 
164036925cSBarry Smith /*@C
17*bcf0153eSBarry Smith    TSMonitorSPEigCtxCreate - Creates a context for use with `TS` to monitor the eigenvalues of the linearized operator
184036925cSBarry Smith 
19*bcf0153eSBarry Smith    Collective
204036925cSBarry Smith 
214036925cSBarry Smith    Input Parameters:
224036925cSBarry Smith +  host - the X display to open, or null for the local machine
234036925cSBarry Smith .  label - the title to put in the title bar
244036925cSBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
254036925cSBarry Smith .  m, n - the screen width and height in pixels
264036925cSBarry Smith -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
274036925cSBarry Smith 
284036925cSBarry Smith    Output Parameter:
294036925cSBarry Smith .  ctx - the context
304036925cSBarry Smith 
314036925cSBarry Smith    Options Database Key:
32b51f60e1SBarry Smith .  -ts_monitor_sp_eig - plot egienvalues of linearized right hand side
334036925cSBarry Smith 
34*bcf0153eSBarry Smith    Level: intermediate
35*bcf0153eSBarry Smith 
364036925cSBarry Smith    Notes:
37*bcf0153eSBarry Smith    Use `TSMonitorSPEigCtxDestroy()` to destroy the context
384036925cSBarry Smith 
39b51f60e1SBarry Smith    Currently only works if the Jacobian is provided explicitly.
40b51f60e1SBarry Smith 
41b51f60e1SBarry Smith    Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix.
42b51f60e1SBarry Smith 
43*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSPEigTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()`
444036925cSBarry Smith @*/
45d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorSPEigCtx *ctx)
46d71ae5a4SJacob Faibussowitsch {
474036925cSBarry Smith   PetscDraw win;
484036925cSBarry Smith   PC        pc;
494036925cSBarry Smith 
504036925cSBarry Smith   PetscFunctionBegin;
519566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
529566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &(*ctx)->rand));
539566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions((*ctx)->rand));
549566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &win));
559566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(win));
569566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(win, 1, &(*ctx)->drawsp));
579566063dSJacob Faibussowitsch   PetscCall(KSPCreate(comm, &(*ctx)->ksp));
589566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix((*ctx)->ksp, "ts_monitor_sp_eig_")); /* this is wrong, used use also prefix from the TS */
599566063dSJacob Faibussowitsch   PetscCall(KSPSetType((*ctx)->ksp, KSPGMRES));
609566063dSJacob Faibussowitsch   PetscCall(KSPGMRESSetRestart((*ctx)->ksp, 200));
619566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances((*ctx)->ksp, 1.e-10, PETSC_DEFAULT, PETSC_DEFAULT, 200));
629566063dSJacob Faibussowitsch   PetscCall(KSPSetComputeSingularValues((*ctx)->ksp, PETSC_TRUE));
639566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions((*ctx)->ksp));
649566063dSJacob Faibussowitsch   PetscCall(KSPGetPC((*ctx)->ksp, &pc));
659566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
66bbd56ea5SKarl Rupp 
674036925cSBarry Smith   (*ctx)->howoften          = howoften;
68b51f60e1SBarry Smith   (*ctx)->computeexplicitly = PETSC_FALSE;
69bbd56ea5SKarl Rupp 
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_sp_eig_explicitly", &(*ctx)->computeexplicitly, NULL));
71bbd56ea5SKarl Rupp 
72b51f60e1SBarry Smith   (*ctx)->comm = comm;
730d18c744SBarry Smith   (*ctx)->xmin = -2.1;
74f9c1d6abSBarry Smith   (*ctx)->xmax = 1.1;
750d18c744SBarry Smith   (*ctx)->ymin = -1.1;
760d18c744SBarry Smith   (*ctx)->ymax = 1.1;
774036925cSBarry Smith   PetscFunctionReturn(0);
784036925cSBarry Smith }
794036925cSBarry Smith 
80d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr, PetscReal xi, PetscBool *flg)
81d71ae5a4SJacob Faibussowitsch {
82f9c1d6abSBarry Smith   PetscReal yr, yi;
830d18c744SBarry Smith 
84f9c1d6abSBarry Smith   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(TSComputeLinearStability(ts, xr, xi, &yr, &yi));
86f9c1d6abSBarry Smith   if ((yr * yr + yi * yi) <= 1.0) *flg = PETSC_TRUE;
87f9c1d6abSBarry Smith   else *flg = PETSC_FALSE;
88f9c1d6abSBarry Smith   PetscFunctionReturn(0);
890d18c744SBarry Smith }
900d18c744SBarry Smith 
91d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEig(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx)
92d71ae5a4SJacob Faibussowitsch {
934036925cSBarry Smith   TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx)monctx;
944036925cSBarry Smith   KSP               ksp = ctx->ksp;
95b51f60e1SBarry Smith   PetscInt          n, N, nits, neig, i, its = 200;
96b296d7d5SJed Brown   PetscReal        *r, *c, time_step_save;
974036925cSBarry Smith   PetscDrawSP       drawsp = ctx->drawsp;
984036925cSBarry Smith   Mat               A, B;
994036925cSBarry Smith   Vec               xdot;
1004036925cSBarry Smith   SNES              snes;
1014036925cSBarry Smith 
1024036925cSBarry Smith   PetscFunctionBegin;
10363e21af5SBarry Smith   if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */
1044036925cSBarry Smith   if (!step) PetscFunctionReturn(0);
105b06615a5SLisandro Dalcin   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
1069566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v, &xdot));
1079566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
1089566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &A, &B, NULL, NULL));
1099566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B));
110b51f60e1SBarry Smith     /*
111b51f60e1SBarry Smith        This doesn't work because methods keep and use internal information about the shift so it
112b51f60e1SBarry Smith        seems we would need code for each method to trick the correct Jacobian in being computed.
113b296d7d5SJed Brown      */
114b296d7d5SJed Brown     time_step_save = ts->time_step;
115b51f60e1SBarry Smith     ts->time_step  = PETSC_MAX_REAL;
116bbd56ea5SKarl Rupp 
1179566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, v, A, B));
118bbd56ea5SKarl Rupp 
119b296d7d5SJed Brown     ts->time_step = time_step_save;
120a174af7bSBarry Smith 
1219566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, B, B));
1229566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &n));
123b51f60e1SBarry Smith     if (n < 200) its = n;
1249566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(ksp, 1.e-10, PETSC_DEFAULT, PETSC_DEFAULT, its));
1259566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xdot, ctx->rand));
1269566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ksp, xdot, xdot));
1279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xdot));
1289566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(ksp, &nits));
129b51f60e1SBarry Smith     N = nits + 2;
1304036925cSBarry Smith 
1314036925cSBarry Smith     if (nits) {
132a174af7bSBarry Smith       PetscDraw     draw;
133d98bdc7eSBarry Smith       PetscReal     pause;
134f9c1d6abSBarry Smith       PetscDrawAxis axis;
135f9c1d6abSBarry Smith       PetscReal     xmin, xmax, ymin, ymax;
136d98bdc7eSBarry Smith 
1379566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPReset(drawsp));
1389566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPSetLimits(drawsp, ctx->xmin, ctx->xmax, ctx->ymin, ctx->ymax));
1399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(PetscMax(n, N), &r, PetscMax(n, N), &c));
140b51f60e1SBarry Smith       if (ctx->computeexplicitly) {
1419566063dSJacob Faibussowitsch         PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
142b51f60e1SBarry Smith         neig = n;
143b51f60e1SBarry Smith       } else {
1449566063dSJacob Faibussowitsch         PetscCall(KSPComputeEigenvalues(ksp, N, r, c, &neig));
145b51f60e1SBarry Smith       }
146b296d7d5SJed 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 */
147b296d7d5SJed Brown       for (i = 0; i < neig; i++) r[i] = -r[i];
1484036925cSBarry Smith       for (i = 0; i < neig; i++) {
1495f737190SBarry Smith         if (ts->ops->linearstability) {
1505f737190SBarry Smith           PetscReal fr, fi;
1519566063dSJacob Faibussowitsch           PetscCall(TSComputeLinearStability(ts, r[i], c[i], &fr, &fi));
15248a46eb9SPierre 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)));
1535d726895SBarry Smith         }
1549566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
1554036925cSBarry Smith       }
1569566063dSJacob Faibussowitsch       PetscCall(PetscFree2(r, c));
1579566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPGetDraw(drawsp, &draw));
1589566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetPause(draw, &pause));
1599566063dSJacob Faibussowitsch       PetscCall(PetscDrawSetPause(draw, 0.0));
1609566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
1619566063dSJacob Faibussowitsch       PetscCall(PetscDrawSetPause(draw, pause));
162f9c1d6abSBarry Smith       if (ts->ops->linearstability) {
1639566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPGetAxis(drawsp, &axis));
1649566063dSJacob Faibussowitsch         PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, &ymin, &ymax));
1659566063dSJacob Faibussowitsch         PetscCall(PetscDrawIndicatorFunction(draw, xmin, xmax, ymin, ymax, PETSC_DRAW_CYAN, (PetscErrorCode(*)(void *, PetscReal, PetscReal, PetscBool *))TSLinearStabilityIndicator, ts));
1669566063dSJacob Faibussowitsch         PetscCall(PetscDrawSPDraw(drawsp, PETSC_FALSE));
1670d18c744SBarry Smith       }
1689566063dSJacob Faibussowitsch       PetscCall(PetscDrawSPSave(drawsp));
1694036925cSBarry Smith     }
1709566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
1714036925cSBarry Smith   }
1724036925cSBarry Smith   PetscFunctionReturn(0);
1734036925cSBarry Smith }
1744036925cSBarry Smith 
1754036925cSBarry Smith /*@C
176*bcf0153eSBarry Smith    TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with `TSMonitorSPEigCtxCreate()`.
1774036925cSBarry Smith 
178*bcf0153eSBarry Smith    Collective on ctx
1794036925cSBarry Smith 
1804036925cSBarry Smith    Input Parameter:
1814036925cSBarry Smith .  ctx - the monitor context
1824036925cSBarry Smith 
1834036925cSBarry Smith    Level: intermediate
1844036925cSBarry Smith 
185*bcf0153eSBarry Smith    Note:
186*bcf0153eSBarry Smith    Should be passed to `TSMonitorSet()` along with `TSMonitorSPEig()` an the context created with `TSMonitorSPEigCtxCreate()`
187*bcf0153eSBarry Smith 
188*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSMonitorSPEigCtxCreate()`, `TSMonitorSet()`, `TSMonitorSPEig();`
1894036925cSBarry Smith @*/
190d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx)
191d71ae5a4SJacob Faibussowitsch {
1924036925cSBarry Smith   PetscDraw draw;
1934036925cSBarry Smith 
1944036925cSBarry Smith   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPGetDraw((*ctx)->drawsp, &draw));
1969566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
1979566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&(*ctx)->drawsp));
1989566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&(*ctx)->ksp));
1999566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&(*ctx)->rand));
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
2014036925cSBarry Smith   PetscFunctionReturn(0);
2024036925cSBarry Smith }
203