1e669de00SBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmplex.h" I*/
2e669de00SBarry Smith #include <petscts.h>
3e669de00SBarry Smith #include <petscdraw.h>
4e669de00SBarry Smith
514d0ab18SJacob Faibussowitsch /*@C
614d0ab18SJacob Faibussowitsch TSMonitorLGCtxNetworkDestroy - Destroys line graph contexts that where created with `TSMonitorLGCtxNetworkCreate()`.
7e669de00SBarry Smith
8c3339decSBarry Smith Collective
9e669de00SBarry Smith
10e669de00SBarry Smith Input Parameter:
11e669de00SBarry Smith . ctx - the monitor context
12e669de00SBarry Smith
1314d0ab18SJacob Faibussowitsch Level: intermediate
1414d0ab18SJacob Faibussowitsch
1514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkSolution()`
1614d0ab18SJacob Faibussowitsch @*/
TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork * ctx)17d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *ctx)
18d71ae5a4SJacob Faibussowitsch {
19e669de00SBarry Smith PetscInt i;
20e669de00SBarry Smith
21e669de00SBarry Smith PetscFunctionBegin;
2248a46eb9SPierre Jolivet for (i = 0; i < (*ctx)->nlg; i++) PetscCall(PetscDrawLGDestroy(&(*ctx)->lg[i]));
239566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->lg));
249566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx));
253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
26e669de00SBarry Smith }
27e669de00SBarry Smith
TSMonitorLGCtxNetworkCreate(TS ts,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtxNetwork * ctx)28d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxNetworkCreate(TS ts, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtxNetwork *ctx)
29d71ae5a4SJacob Faibussowitsch {
30e669de00SBarry Smith PetscDraw draw;
31e669de00SBarry Smith MPI_Comm comm;
32e669de00SBarry Smith DM dm;
33e669de00SBarry Smith PetscInt i, Start, End, e, nvar;
34e669de00SBarry Smith
35e669de00SBarry Smith PetscFunctionBegin;
369566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
379566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
389566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx));
39e669de00SBarry Smith i = 0;
40e669de00SBarry Smith /* loop over edges counting number of line graphs needed */
419566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
42e669de00SBarry Smith for (e = Start; e < End; e++) {
439566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
44e669de00SBarry Smith if (!nvar) continue;
45e669de00SBarry Smith i++;
46e669de00SBarry Smith }
47e669de00SBarry Smith /* loop over vertices */
489566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
49e669de00SBarry Smith for (e = Start; e < End; e++) {
509566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
51e669de00SBarry Smith if (!nvar) continue;
52e669de00SBarry Smith i++;
53e669de00SBarry Smith }
54e669de00SBarry Smith (*ctx)->nlg = i;
559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i, &(*ctx)->lg));
56e669de00SBarry Smith
57e669de00SBarry Smith i = 0;
58e669de00SBarry Smith /* loop over edges creating all needed line graphs*/
599566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
60e669de00SBarry Smith for (e = Start; e < End; e++) {
619566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
62e669de00SBarry Smith if (!nvar) continue;
639566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
649566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw));
659566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
669566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
679566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw));
68e669de00SBarry Smith i++;
69e669de00SBarry Smith }
70e669de00SBarry Smith /* loop over vertices */
719566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
72e669de00SBarry Smith for (e = Start; e < End; e++) {
739566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
74e669de00SBarry Smith if (!nvar) continue;
759566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
769566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw));
779566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
789566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
799566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw));
80e669de00SBarry Smith i++;
81e669de00SBarry Smith }
829566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw));
83e669de00SBarry Smith (*ctx)->howoften = howoften;
843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
85e669de00SBarry Smith }
86e669de00SBarry Smith
8714d0ab18SJacob Faibussowitsch /*@C
88bcf0153eSBarry Smith TSMonitorLGCtxNetworkSolution - Monitors progress of the `TS` solvers for a `DMNETWORK` solution with one window for each vertex and each edge
89e669de00SBarry Smith
90c3339decSBarry Smith Collective
91e669de00SBarry Smith
92e669de00SBarry Smith Input Parameters:
93bcf0153eSBarry Smith + ts - the `TS` context
94e669de00SBarry Smith . step - current time-step
95e669de00SBarry Smith . ptime - current time
96e669de00SBarry Smith . u - current solution
97bcf0153eSBarry Smith - dctx - the `TSMonitorLGCtxNetwork` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreateNetwork()`
98e669de00SBarry Smith
99bcf0153eSBarry Smith Options Database Key:
100b43aa488SJacob Faibussowitsch . -ts_monitor_lg_solution_variables - monitor solution variables
101e669de00SBarry Smith
102e669de00SBarry Smith Level: intermediate
103e669de00SBarry Smith
104bcf0153eSBarry Smith Note:
10514d0ab18SJacob Faibussowitsch Each process in a parallel run displays its component solutions in a separate graphics window
106e669de00SBarry Smith
10714d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkDestroy()`
10814d0ab18SJacob Faibussowitsch @*/
TSMonitorLGCtxNetworkSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void * dctx)109d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
110d71ae5a4SJacob Faibussowitsch {
111*835f2295SStefano Zampini #if defined(PETSC_USE_COMPLEX)
112*835f2295SStefano Zampini PetscFunctionBegin;
113*835f2295SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
114*835f2295SStefano Zampini #else
115e669de00SBarry Smith TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx;
116e669de00SBarry Smith const PetscScalar *xv;
117e669de00SBarry Smith PetscScalar *yv;
118e669de00SBarry Smith PetscInt i, v, Start, End, offset, nvar, e;
119e669de00SBarry Smith TSConvergedReason reason;
120e669de00SBarry Smith DM dm;
121e669de00SBarry Smith Vec uv;
122e669de00SBarry Smith
123e669de00SBarry Smith PetscFunctionBegin;
1243ba16761SJacob Faibussowitsch if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
125e669de00SBarry Smith if (!step) {
126e669de00SBarry Smith PetscDrawAxis axis;
127e669de00SBarry Smith
128e669de00SBarry Smith for (i = 0; i < ctx->nlg; i++) {
1299566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg[i], &axis));
1309566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
1319566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg[i]));
132e669de00SBarry Smith }
133e669de00SBarry Smith }
134e669de00SBarry Smith
135e669de00SBarry Smith if (ctx->semilogy) {
136e669de00SBarry Smith PetscInt n, j;
137e669de00SBarry Smith
1389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &uv));
1399566063dSJacob Faibussowitsch PetscCall(VecCopy(u, uv));
1409566063dSJacob Faibussowitsch PetscCall(VecGetArray(uv, &yv));
1419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(uv, &n));
142e669de00SBarry Smith for (j = 0; j < n; j++) {
1434500feddSGetnet Betrie if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12;
1444e936a74SSatish Balay else yv[j] = PetscLog10Real(PetscRealPart(yv[j]));
145e669de00SBarry Smith }
146e669de00SBarry Smith xv = yv;
147e669de00SBarry Smith } else {
1489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &xv));
149e669de00SBarry Smith }
150e669de00SBarry Smith /* iterate over edges */
1519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
152e669de00SBarry Smith i = 0;
1539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
154e669de00SBarry Smith for (e = Start; e < End; e++) {
1559566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
156e669de00SBarry Smith if (!nvar) continue;
157e669de00SBarry Smith
1589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(dm, e, ALL_COMPONENTS, &offset));
159*835f2295SStefano Zampini PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, xv + offset));
160e669de00SBarry Smith i++;
161e669de00SBarry Smith }
162e669de00SBarry Smith
163e669de00SBarry Smith /* iterate over vertices */
1649566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
165e669de00SBarry Smith for (v = Start; v < End; v++) {
1669566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
167e669de00SBarry Smith if (!nvar) continue;
168e669de00SBarry Smith
1699566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(dm, v, ALL_COMPONENTS, &offset));
170*835f2295SStefano Zampini PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, xv + offset));
171e669de00SBarry Smith i++;
172e669de00SBarry Smith }
173e669de00SBarry Smith if (ctx->semilogy) {
1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(uv, &yv));
1759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uv));
176e669de00SBarry Smith } else {
1779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &xv));
178e669de00SBarry Smith }
179e669de00SBarry Smith
1809566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason));
181e669de00SBarry Smith if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) {
182e669de00SBarry Smith for (i = 0; i < ctx->nlg; i++) {
1839566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg[i]));
1849566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg[i]));
185e669de00SBarry Smith }
186e669de00SBarry Smith }
1873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
188*835f2295SStefano Zampini #endif
189e669de00SBarry Smith }
190