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 17bcf0153eSBarry Smith TSMonitorSPEigCtxCreate - Creates a context for use with `TS` to monitor the eigenvalues of the linearized operator 184036925cSBarry Smith 19bcf0153eSBarry Smith Collective 204036925cSBarry Smith 214036925cSBarry Smith Input Parameters: 222fe279fdSBarry Smith + comm - the communicator to share the monitor 232fe279fdSBarry Smith . host - the X display to open, or `NULL` for the local machine 244036925cSBarry Smith . label - the title to put in the title bar 252fe279fdSBarry Smith . x - the horizontal screen coordinates of the upper left coordinate of the window 262fe279fdSBarry Smith . y - the vertical coordinates of the upper left coordinate of the window 272fe279fdSBarry Smith . m - the screen width in pixels 282fe279fdSBarry Smith . n - the screen height in pixels 294036925cSBarry Smith - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time 304036925cSBarry Smith 314036925cSBarry Smith Output Parameter: 324036925cSBarry Smith . ctx - the context 334036925cSBarry Smith 344036925cSBarry Smith Options Database Key: 35b51f60e1SBarry Smith . -ts_monitor_sp_eig - plot egienvalues of linearized right hand side 364036925cSBarry Smith 37bcf0153eSBarry Smith Level: intermediate 38bcf0153eSBarry Smith 394036925cSBarry Smith Notes: 40bcf0153eSBarry Smith Use `TSMonitorSPEigCtxDestroy()` to destroy the context 414036925cSBarry Smith 42b51f60e1SBarry Smith Currently only works if the Jacobian is provided explicitly. 43b51f60e1SBarry Smith 44b51f60e1SBarry Smith Currently only works for ODEs u_t - F(t,u) = 0; that is with no mass matrix. 45b51f60e1SBarry Smith 46*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigTimeStep()`, `TSMonitorSet()`, `TSMonitorLGSolution()`, `TSMonitorLGError()` 474036925cSBarry Smith @*/ 48d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorSPEigCtx *ctx) 49d71ae5a4SJacob Faibussowitsch { 504036925cSBarry Smith PetscDraw win; 514036925cSBarry Smith PC pc; 524036925cSBarry Smith 534036925cSBarry Smith PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 559566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &(*ctx)->rand)); 569566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions((*ctx)->rand)); 579566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &win)); 589566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(win)); 599566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(win, 1, &(*ctx)->drawsp)); 609566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &(*ctx)->ksp)); 619566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix((*ctx)->ksp, "ts_monitor_sp_eig_")); /* this is wrong, used use also prefix from the TS */ 629566063dSJacob Faibussowitsch PetscCall(KSPSetType((*ctx)->ksp, KSPGMRES)); 639566063dSJacob Faibussowitsch PetscCall(KSPGMRESSetRestart((*ctx)->ksp, 200)); 649566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances((*ctx)->ksp, 1.e-10, PETSC_DEFAULT, PETSC_DEFAULT, 200)); 659566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues((*ctx)->ksp, PETSC_TRUE)); 669566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions((*ctx)->ksp)); 679566063dSJacob Faibussowitsch PetscCall(KSPGetPC((*ctx)->ksp, &pc)); 689566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 69bbd56ea5SKarl Rupp 704036925cSBarry Smith (*ctx)->howoften = howoften; 71b51f60e1SBarry Smith (*ctx)->computeexplicitly = PETSC_FALSE; 72bbd56ea5SKarl Rupp 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_sp_eig_explicitly", &(*ctx)->computeexplicitly, NULL)); 74bbd56ea5SKarl Rupp 75b51f60e1SBarry Smith (*ctx)->comm = comm; 760d18c744SBarry Smith (*ctx)->xmin = -2.1; 77f9c1d6abSBarry Smith (*ctx)->xmax = 1.1; 780d18c744SBarry Smith (*ctx)->ymin = -1.1; 790d18c744SBarry Smith (*ctx)->ymax = 1.1; 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 814036925cSBarry Smith } 824036925cSBarry Smith 83d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLinearStabilityIndicator(TS ts, PetscReal xr, PetscReal xi, PetscBool *flg) 84d71ae5a4SJacob Faibussowitsch { 85f9c1d6abSBarry Smith PetscReal yr, yi; 860d18c744SBarry Smith 87f9c1d6abSBarry Smith PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, xr, xi, &yr, &yi)); 89f9c1d6abSBarry Smith if ((yr * yr + yi * yi) <= 1.0) *flg = PETSC_TRUE; 90f9c1d6abSBarry Smith else *flg = PETSC_FALSE; 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 920d18c744SBarry Smith } 930d18c744SBarry Smith 94d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEig(TS ts, PetscInt step, PetscReal ptime, Vec v, void *monctx) 95d71ae5a4SJacob Faibussowitsch { 964036925cSBarry Smith TSMonitorSPEigCtx ctx = (TSMonitorSPEigCtx)monctx; 974036925cSBarry Smith KSP ksp = ctx->ksp; 98b51f60e1SBarry Smith PetscInt n, N, nits, neig, i, its = 200; 99b296d7d5SJed Brown PetscReal *r, *c, time_step_save; 1004036925cSBarry Smith PetscDrawSP drawsp = ctx->drawsp; 1014036925cSBarry Smith Mat A, B; 1024036925cSBarry Smith Vec xdot; 1034036925cSBarry Smith SNES snes; 1044036925cSBarry Smith 1054036925cSBarry Smith PetscFunctionBegin; 1063ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 1073ba16761SJacob Faibussowitsch if (!step) PetscFunctionReturn(PETSC_SUCCESS); 108b06615a5SLisandro Dalcin if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) { 1099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v, &xdot)); 1109566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1119566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &A, &B, NULL, NULL)); 1129566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B)); 113b51f60e1SBarry Smith /* 114b51f60e1SBarry Smith This doesn't work because methods keep and use internal information about the shift so it 115b51f60e1SBarry Smith seems we would need code for each method to trick the correct Jacobian in being computed. 116b296d7d5SJed Brown */ 117b296d7d5SJed Brown time_step_save = ts->time_step; 118b51f60e1SBarry Smith ts->time_step = PETSC_MAX_REAL; 119bbd56ea5SKarl Rupp 1209566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, v, A, B)); 121bbd56ea5SKarl Rupp 122b296d7d5SJed Brown ts->time_step = time_step_save; 123a174af7bSBarry Smith 1249566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, B, B)); 1259566063dSJacob Faibussowitsch PetscCall(VecGetSize(v, &n)); 126b51f60e1SBarry Smith if (n < 200) its = n; 1279566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp, 1.e-10, PETSC_DEFAULT, PETSC_DEFAULT, its)); 1289566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xdot, ctx->rand)); 1299566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, xdot, xdot)); 1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xdot)); 1319566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &nits)); 132b51f60e1SBarry Smith N = nits + 2; 1334036925cSBarry Smith 1344036925cSBarry Smith if (nits) { 135a174af7bSBarry Smith PetscDraw draw; 136d98bdc7eSBarry Smith PetscReal pause; 137f9c1d6abSBarry Smith PetscDrawAxis axis; 138f9c1d6abSBarry Smith PetscReal xmin, xmax, ymin, ymax; 139d98bdc7eSBarry Smith 1409566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(drawsp)); 1419566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(drawsp, ctx->xmin, ctx->xmax, ctx->ymin, ctx->ymax)); 1429566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(n, N), &r, PetscMax(n, N), &c)); 143b51f60e1SBarry Smith if (ctx->computeexplicitly) { 1449566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c)); 145b51f60e1SBarry Smith neig = n; 146b51f60e1SBarry Smith } else { 1479566063dSJacob Faibussowitsch PetscCall(KSPComputeEigenvalues(ksp, N, r, c, &neig)); 148b51f60e1SBarry Smith } 149b296d7d5SJed 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 */ 150b296d7d5SJed Brown for (i = 0; i < neig; i++) r[i] = -r[i]; 1514036925cSBarry Smith for (i = 0; i < neig; i++) { 1525f737190SBarry Smith if (ts->ops->linearstability) { 1535f737190SBarry Smith PetscReal fr, fi; 1549566063dSJacob Faibussowitsch PetscCall(TSComputeLinearStability(ts, r[i], c[i], &fr, &fi)); 15548a46eb9SPierre 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))); 1565d726895SBarry Smith } 1579566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i)); 1584036925cSBarry Smith } 1599566063dSJacob Faibussowitsch PetscCall(PetscFree2(r, c)); 1609566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw(drawsp, &draw)); 1619566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 1629566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 1639566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE)); 1649566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 165f9c1d6abSBarry Smith if (ts->ops->linearstability) { 1669566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(drawsp, &axis)); 1679566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, &ymin, &ymax)); 1689566063dSJacob Faibussowitsch PetscCall(PetscDrawIndicatorFunction(draw, xmin, xmax, ymin, ymax, PETSC_DRAW_CYAN, (PetscErrorCode(*)(void *, PetscReal, PetscReal, PetscBool *))TSLinearStabilityIndicator, ts)); 1699566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(drawsp, PETSC_FALSE)); 1700d18c744SBarry Smith } 1719566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSave(drawsp)); 1724036925cSBarry Smith } 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1744036925cSBarry Smith } 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1764036925cSBarry Smith } 1774036925cSBarry Smith 1784036925cSBarry Smith /*@C 179bcf0153eSBarry Smith TSMonitorSPEigCtxDestroy - Destroys a scatter plot context that was created with `TSMonitorSPEigCtxCreate()`. 1804036925cSBarry Smith 181c3339decSBarry Smith Collective 1824036925cSBarry Smith 1834036925cSBarry Smith Input Parameter: 1844036925cSBarry Smith . ctx - the monitor context 1854036925cSBarry Smith 1864036925cSBarry Smith Level: intermediate 1874036925cSBarry Smith 188bcf0153eSBarry Smith Note: 189bcf0153eSBarry Smith Should be passed to `TSMonitorSet()` along with `TSMonitorSPEig()` an the context created with `TSMonitorSPEigCtxCreate()` 190bcf0153eSBarry Smith 191*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSMonitorSPEigCtxCreate()`, `TSMonitorSet()`, `TSMonitorSPEig();` 1924036925cSBarry Smith @*/ 193d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *ctx) 194d71ae5a4SJacob Faibussowitsch { 1954036925cSBarry Smith PetscDraw draw; 1964036925cSBarry Smith 1974036925cSBarry Smith PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetDraw((*ctx)->drawsp, &draw)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&(*ctx)->drawsp)); 2019566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&(*ctx)->ksp)); 2029566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&(*ctx)->rand)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2054036925cSBarry Smith } 206