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 174036925cSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with TS to monitor the eigenvalues of the linearized operator 184036925cSBarry Smith 194036925cSBarry Smith Collective on TS 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 344036925cSBarry Smith Notes: 354036925cSBarry Smith Use TSMonitorSPEigCtxDestroy() to destroy. 364036925cSBarry Smith 37b51f60e1SBarry Smith Currently only works if the Jacobian is provided explicitly. 38b51f60e1SBarry Smith 39b51f60e1SBarry Smith Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix. 40b51f60e1SBarry Smith 414036925cSBarry Smith Level: intermediate 424036925cSBarry Smith 43db781477SPatrick Sanan .seealso: `TSMonitorSPEigTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 444036925cSBarry Smith 454036925cSBarry Smith @*/ 469371c9d4SSatish Balay PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorSPEigCtx *ctx) { 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 809371c9d4SSatish Balay static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr, PetscReal xi, PetscBool *flg) { 81f9c1d6abSBarry Smith PetscReal yr, yi; 820d18c744SBarry Smith 83f9c1d6abSBarry Smith PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, xr, xi, &yr, &yi)); 85f9c1d6abSBarry Smith if ((yr * yr + yi * yi) <= 1.0) *flg = PETSC_TRUE; 86f9c1d6abSBarry Smith else *flg = PETSC_FALSE; 87f9c1d6abSBarry Smith PetscFunctionReturn(0); 880d18c744SBarry Smith } 890d18c744SBarry Smith 909371c9d4SSatish Balay PetscErrorCode TSMonitorSPEig(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) { 914036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx)monctx; 924036925cSBarry Smith KSP ksp = ctx->ksp; 93b51f60e1SBarry Smith PetscInt n, N, nits, neig, i, its = 200; 94b296d7d5SJed Brown PetscReal *r, *c, time_step_save; 954036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 964036925cSBarry Smith Mat A, B; 974036925cSBarry Smith Vec xdot; 984036925cSBarry Smith SNES snes; 994036925cSBarry Smith 1004036925cSBarry Smith PetscFunctionBegin; 10163e21af5SBarry Smith if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 1024036925cSBarry Smith if (!step) PetscFunctionReturn(0); 103b06615a5SLisandro Dalcin if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v, &xdot)); 1059566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1069566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &A, &B, NULL, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B)); 108b51f60e1SBarry Smith /* 109b51f60e1SBarry Smith This doesn't work because methods keep and use internal information about the shift so it 110b51f60e1SBarry Smith seems we would need code for each method to trick the correct Jacobian in being computed. 111b296d7d5SJed Brown */ 112b296d7d5SJed Brown time_step_save = ts->time_step; 113b51f60e1SBarry Smith ts->time_step = PETSC_MAX_REAL; 114bbd56ea5SKarl Rupp 1159566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, v, A, B)); 116bbd56ea5SKarl Rupp 117b296d7d5SJed Brown ts->time_step = time_step_save; 118a174af7bSBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, B, B)); 1209566063dSJacob Faibussowitsch PetscCall(VecGetSize(v, &n)); 121b51f60e1SBarry Smith if (n < 200) its = n; 1229566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp, 1.e-10, PETSC_DEFAULT, PETSC_DEFAULT, its)); 1239566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xdot, ctx->rand)); 1249566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, xdot, xdot)); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xdot)); 1269566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &nits)); 127b51f60e1SBarry Smith N = nits + 2; 1284036925cSBarry Smith 1294036925cSBarry Smith if (nits) { 130a174af7bSBarry Smith PetscDraw draw; 131d98bdc7eSBarry Smith PetscReal pause; 132f9c1d6abSBarry Smith PetscDrawAxis axis; 133f9c1d6abSBarry Smith PetscReal xmin, xmax, ymin, ymax; 134d98bdc7eSBarry Smith 1359566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(drawsp)); 1369566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(drawsp, ctx->xmin, ctx->xmax, ctx->ymin, ctx->ymax)); 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(n, N), &r, PetscMax(n, N), &c)); 138b51f60e1SBarry Smith if (ctx->computeexplicitly) { 1399566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c)); 140b51f60e1SBarry Smith neig = n; 141b51f60e1SBarry Smith } else { 1429566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvalues(ksp, N, r, c, &neig)); 143b51f60e1SBarry Smith } 144b296d7d5SJed 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 */ 145b296d7d5SJed Brown for (i = 0; i < neig; i++) r[i] = -r[i]; 1464036925cSBarry Smith for (i = 0; i < neig; i++) { 1475f737190SBarry Smith if (ts->ops->linearstability) { 1485f737190SBarry Smith PetscReal fr, fi; 1499566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, r[i], c[i], &fr, &fi)); 150*48a46eb9SPierre 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))); 1515d726895SBarry Smith } 1529566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i)); 1534036925cSBarry Smith } 1549566063dSJacob Faibussowitsch PetscCall(PetscFree2(r, c)); 1559566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(drawsp, &draw)); 1569566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 1579566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 1589566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE)); 1599566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 160f9c1d6abSBarry Smith if (ts->ops->linearstability) { 1619566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(drawsp, &axis)); 1629566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, &ymin, &ymax)); 1639566063dSJacob Faibussowitsch PetscCall(PetscDrawIndicatorFunction(draw, xmin, xmax, ymin, ymax, PETSC_DRAW_CYAN, (PetscErrorCode(*)(void *, PetscReal, PetscReal, PetscBool *))TSLinearStabilityIndicator, ts)); 1649566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_FALSE)); 1650d18c744SBarry Smith } 1669566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(drawsp)); 1674036925cSBarry Smith } 1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1694036925cSBarry Smith } 1704036925cSBarry Smith PetscFunctionReturn(0); 1714036925cSBarry Smith } 1724036925cSBarry Smith 1734036925cSBarry Smith /*@C 1744036925cSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with TSMonitorSPEigCtxCreate(). 1754036925cSBarry Smith 1764036925cSBarry Smith Collective on TSMonitorSPEigCtx 1774036925cSBarry Smith 1784036925cSBarry Smith Input Parameter: 1794036925cSBarry Smith . ctx - the monitor context 1804036925cSBarry Smith 1814036925cSBarry Smith Level: intermediate 1824036925cSBarry Smith 183db781477SPatrick Sanan .seealso: `TSMonitorSPEigCtxCreate()`, `TSMonitorSet()`, `TSMonitorSPEig();` 1844036925cSBarry Smith @*/ 1859371c9d4SSatish Balay PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) { 1864036925cSBarry Smith PetscDraw draw; 1874036925cSBarry Smith 1884036925cSBarry Smith PetscFunctionBegin; 1899566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw((*ctx)->drawsp, &draw)); 1909566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 1919566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->drawsp)); 1929566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&(*ctx)->ksp)); 1939566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&(*ctx)->rand)); 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 1954036925cSBarry Smith PetscFunctionReturn(0); 1964036925cSBarry Smith } 197