1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 29804daf3SBarry Smith #include <petscdraw.h> 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 */ 9b51f60e1SBarry Smith PetscBool computeexplicitly; 10b51f60e1SBarry Smith MPI_Comm comm; 11b51f60e1SBarry Smith PetscRandom rand; 120d18c744SBarry Smith PetscReal xmin, xmax, ymin, ymax; 134036925cSBarry Smith }; 144036925cSBarry Smith 154036925cSBarry Smith /*@C 16bcf0153eSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with `TS` to monitor the eigenvalues of the linearized operator 174036925cSBarry Smith 18bcf0153eSBarry Smith Collective 194036925cSBarry Smith 204036925cSBarry Smith Input Parameters: 212fe279fdSBarry Smith + comm - the communicator to share the monitor 222fe279fdSBarry Smith . host - the X display to open, or `NULL` for the local machine 234036925cSBarry Smith . label - the title to put in the title bar 242fe279fdSBarry Smith . x - the horizontal screen coordinates of the upper left coordinate of the window 252fe279fdSBarry Smith . y - the vertical coordinates of the upper left coordinate of the window 262fe279fdSBarry Smith . m - the screen width in pixels 272fe279fdSBarry Smith . n - the screen 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: 34dd8e379bSPierre Jolivet . -ts_monitor_sp_eig - plot egienvalues of linearized right-hand side 354036925cSBarry Smith 36bcf0153eSBarry Smith Level: intermediate 37bcf0153eSBarry Smith 384036925cSBarry Smith Notes: 39bcf0153eSBarry Smith Use `TSMonitorSPEigCtxDestroy()` to destroy the context 404036925cSBarry Smith 41b51f60e1SBarry Smith Currently only works if the Jacobian is provided explicitly. 42b51f60e1SBarry Smith 43b51f60e1SBarry Smith Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix. 44b51f60e1SBarry Smith 451cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 464036925cSBarry Smith @*/ 47d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorSPEigCtx *ctx) 48d71ae5a4SJacob Faibussowitsch { 494036925cSBarry Smith PetscDraw win; 504036925cSBarry Smith PC pc; 514036925cSBarry Smith 524036925cSBarry Smith PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 549566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &(*ctx)->rand)); 559566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions((*ctx)->rand)); 569566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &win)); 579566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(win)); 589566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(win, 1, &(*ctx)->drawsp)); 599566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &(*ctx)->ksp)); 609566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix((*ctx)->ksp, "ts_monitor_sp_eig_")); /* this is wrong, used use also prefix from the TS */ 619566063dSJacob Faibussowitsch PetscCall(KSPSetType((*ctx)->ksp, KSPGMRES)); 629566063dSJacob Faibussowitsch PetscCall(KSPGMRESSetRestart((*ctx)->ksp, 200)); 63*09cb0f53SBarry Smith PetscCall(KSPSetTolerances((*ctx)->ksp, 1.e-10, PETSC_CURRENT, PETSC_CURRENT, 200)); 649566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues((*ctx)->ksp, PETSC_TRUE)); 659566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions((*ctx)->ksp)); 669566063dSJacob Faibussowitsch PetscCall(KSPGetPC((*ctx)->ksp, &pc)); 679566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 68bbd56ea5SKarl Rupp 694036925cSBarry Smith (*ctx)->howoften = howoften; 70b51f60e1SBarry Smith (*ctx)->computeexplicitly = PETSC_FALSE; 71bbd56ea5SKarl Rupp 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_sp_eig_explicitly", &(*ctx)->computeexplicitly, NULL)); 73bbd56ea5SKarl Rupp 74b51f60e1SBarry Smith (*ctx)->comm = comm; 750d18c744SBarry Smith (*ctx)->xmin = -2.1; 76f9c1d6abSBarry Smith (*ctx)->xmax = 1.1; 770d18c744SBarry Smith (*ctx)->ymin = -1.1; 780d18c744SBarry Smith (*ctx)->ymax = 1.1; 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 804036925cSBarry Smith } 814036925cSBarry Smith 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr, PetscReal xi, PetscBool *flg) 83d71ae5a4SJacob Faibussowitsch { 84f9c1d6abSBarry Smith PetscReal yr, yi; 850d18c744SBarry Smith 86f9c1d6abSBarry Smith PetscFunctionBegin; 879566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, xr, xi, &yr, &yi)); 88f9c1d6abSBarry Smith if ((yr * yr + yi * yi) <= 1.0) *flg = PETSC_TRUE; 89f9c1d6abSBarry Smith else *flg = PETSC_FALSE; 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 910d18c744SBarry Smith } 920d18c744SBarry Smith 93d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEig(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 94d71ae5a4SJacob Faibussowitsch { 954036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx)monctx; 964036925cSBarry Smith KSP ksp = ctx->ksp; 97b51f60e1SBarry Smith PetscInt n, N, nits, neig, i, its = 200; 98b296d7d5SJed Brown PetscReal *r, *c, time_step_save; 994036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 1004036925cSBarry Smith Mat A, B; 1014036925cSBarry Smith Vec xdot; 1024036925cSBarry Smith SNES snes; 1034036925cSBarry Smith 1044036925cSBarry Smith PetscFunctionBegin; 1053ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1063ba16761SJacob Faibussowitsch if (!step) PetscFunctionReturn(PETSC_SUCCESS); 107b06615a5SLisandro Dalcin if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v, &xdot)); 1099566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1109566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &A, &B, NULL, NULL)); 1119566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B)); 112b51f60e1SBarry Smith /* 113b51f60e1SBarry Smith This doesn't work because methods keep and use internal information about the shift so it 114b51f60e1SBarry Smith seems we would need code for each method to trick the correct Jacobian in being computed. 115b296d7d5SJed Brown */ 116b296d7d5SJed Brown time_step_save = ts->time_step; 117b51f60e1SBarry Smith ts->time_step = PETSC_MAX_REAL; 118bbd56ea5SKarl Rupp 1199566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, v, A, B)); 120bbd56ea5SKarl Rupp 121b296d7d5SJed Brown ts->time_step = time_step_save; 122a174af7bSBarry Smith 1239566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, B, B)); 1249566063dSJacob Faibussowitsch PetscCall(VecGetSize(v, &n)); 125b51f60e1SBarry Smith if (n < 200) its = n; 126*09cb0f53SBarry Smith PetscCall(KSPSetTolerances(ksp, 1.e-10, PETSC_CURRENT, PETSC_CURRENT, its)); 1279566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xdot, ctx->rand)); 1289566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, xdot, xdot)); 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xdot)); 1309566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &nits)); 131b51f60e1SBarry Smith N = nits + 2; 1324036925cSBarry Smith 1334036925cSBarry Smith if (nits) { 134a174af7bSBarry Smith PetscDraw draw; 135d98bdc7eSBarry Smith PetscReal pause; 136f9c1d6abSBarry Smith PetscDrawAxis axis; 137f9c1d6abSBarry Smith PetscReal xmin, xmax, ymin, ymax; 138d98bdc7eSBarry Smith 1399566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(drawsp)); 1409566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(drawsp, ctx->xmin, ctx->xmax, ctx->ymin, ctx->ymax)); 1419566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(n, N), &r, PetscMax(n, N), &c)); 142b51f60e1SBarry Smith if (ctx->computeexplicitly) { 1439566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c)); 144b51f60e1SBarry Smith neig = n; 145b51f60e1SBarry Smith } else { 1469566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvalues(ksp, N, r, c, &neig)); 147b51f60e1SBarry Smith } 148b296d7d5SJed 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 */ 149b296d7d5SJed Brown for (i = 0; i < neig; i++) r[i] = -r[i]; 1504036925cSBarry Smith for (i = 0; i < neig; i++) { 1515f737190SBarry Smith if (ts->ops->linearstability) { 1525f737190SBarry Smith PetscReal fr, fi; 1539566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, r[i], c[i], &fr, &fi)); 15448a46eb9SPierre 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))); 1555d726895SBarry Smith } 1569566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i)); 1574036925cSBarry Smith } 1589566063dSJacob Faibussowitsch PetscCall(PetscFree2(r, c)); 1599566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(drawsp, &draw)); 1609566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 1619566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 1629566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE)); 1639566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 164f9c1d6abSBarry Smith if (ts->ops->linearstability) { 1659566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(drawsp, &axis)); 1669566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, &ymin, &ymax)); 1679566063dSJacob Faibussowitsch PetscCall(PetscDrawIndicatorFunction(draw, xmin, xmax, ymin, ymax, PETSC_DRAW_CYAN, (PetscErrorCode(*)(void *, PetscReal, PetscReal, PetscBool *))TSLinearStabilityIndicator, ts)); 1689566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_FALSE)); 1690d18c744SBarry Smith } 1709566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(drawsp)); 1714036925cSBarry Smith } 1729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1734036925cSBarry Smith } 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1754036925cSBarry Smith } 1764036925cSBarry Smith 1774036925cSBarry Smith /*@C 178bcf0153eSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with `TSMonitorSPEigCtxCreate()`. 1794036925cSBarry Smith 180c3339decSBarry Smith Collective 1814036925cSBarry Smith 1824036925cSBarry Smith Input Parameter: 1834036925cSBarry Smith . ctx - the monitor context 1844036925cSBarry Smith 1854036925cSBarry Smith Level: intermediate 1864036925cSBarry Smith 187bcf0153eSBarry Smith Note: 188bcf0153eSBarry Smith Should be passed to `TSMonitorSet()` along with `TSMonitorSPEig()` an the context created with `TSMonitorSPEigCtxCreate()` 189bcf0153eSBarry Smith 190a54bb2a9SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigCtxCreate()`, `TSMonitorSet()`, `TSMonitorSPEig()` 1914036925cSBarry Smith @*/ 192d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) 193d71ae5a4SJacob Faibussowitsch { 1944036925cSBarry Smith PetscDraw draw; 1954036925cSBarry Smith 1964036925cSBarry Smith PetscFunctionBegin; 1979566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw((*ctx)->drawsp, &draw)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->drawsp)); 2009566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&(*ctx)->ksp)); 2019566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&(*ctx)->rand)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2044036925cSBarry Smith } 205