xref: /petsc/src/ts/utils/dmnetworkts.c (revision 14d0ab18b70ee075d422f67a0c1395817de67fab)
1e669de00SBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmplex.h" I*/
2e669de00SBarry Smith #include <petscts.h>
3e669de00SBarry Smith #include <petscdraw.h>
4e669de00SBarry Smith 
5*14d0ab18SJacob Faibussowitsch /*@C
6*14d0ab18SJacob 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 
13*14d0ab18SJacob Faibussowitsch   Level: intermediate
14*14d0ab18SJacob Faibussowitsch 
15*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkSolution()`
16*14d0ab18SJacob Faibussowitsch @*/
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 
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 
87*14d0ab18SJacob 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:
105*14d0ab18SJacob Faibussowitsch   Each process in a parallel run displays its component solutions in a separate graphics window
106e669de00SBarry Smith 
107*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkDestroy()`
108*14d0ab18SJacob Faibussowitsch @*/
109d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
110d71ae5a4SJacob Faibussowitsch {
111e669de00SBarry Smith   TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx;
112e669de00SBarry Smith   const PetscScalar    *xv;
113e669de00SBarry Smith   PetscScalar          *yv;
114e669de00SBarry Smith   PetscInt              i, v, Start, End, offset, nvar, e;
115e669de00SBarry Smith   TSConvergedReason     reason;
116e669de00SBarry Smith   DM                    dm;
117e669de00SBarry Smith   Vec                   uv;
118e669de00SBarry Smith 
119e669de00SBarry Smith   PetscFunctionBegin;
1203ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
121e669de00SBarry Smith   if (!step) {
122e669de00SBarry Smith     PetscDrawAxis axis;
123e669de00SBarry Smith 
124e669de00SBarry Smith     for (i = 0; i < ctx->nlg; i++) {
1259566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGGetAxis(ctx->lg[i], &axis));
1269566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
1279566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGReset(ctx->lg[i]));
128e669de00SBarry Smith     }
129e669de00SBarry Smith   }
130e669de00SBarry Smith 
131e669de00SBarry Smith   if (ctx->semilogy) {
132e669de00SBarry Smith     PetscInt n, j;
133e669de00SBarry Smith 
1349566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uv));
1359566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, uv));
1369566063dSJacob Faibussowitsch     PetscCall(VecGetArray(uv, &yv));
1379566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(uv, &n));
138e669de00SBarry Smith     for (j = 0; j < n; j++) {
1394500feddSGetnet Betrie       if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12;
1404e936a74SSatish Balay       else yv[j] = PetscLog10Real(PetscRealPart(yv[j]));
141e669de00SBarry Smith     }
142e669de00SBarry Smith     xv = yv;
143e669de00SBarry Smith   } else {
1449566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(u, &xv));
145e669de00SBarry Smith   }
146e669de00SBarry Smith   /* iterate over edges */
1479566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
148e669de00SBarry Smith   i = 0;
1499566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
150e669de00SBarry Smith   for (e = Start; e < End; e++) {
1519566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
152e669de00SBarry Smith     if (!nvar) continue;
153e669de00SBarry Smith 
1549566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(dm, e, ALL_COMPONENTS, &offset));
1559566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
156e669de00SBarry Smith     i++;
157e669de00SBarry Smith   }
158e669de00SBarry Smith 
159e669de00SBarry Smith   /* iterate over vertices */
1609566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
161e669de00SBarry Smith   for (v = Start; v < End; v++) {
1629566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
163e669de00SBarry Smith     if (!nvar) continue;
164e669de00SBarry Smith 
1659566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(dm, v, ALL_COMPONENTS, &offset));
1669566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
167e669de00SBarry Smith     i++;
168e669de00SBarry Smith   }
169e669de00SBarry Smith   if (ctx->semilogy) {
1709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(uv, &yv));
1719566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uv));
172e669de00SBarry Smith   } else {
1739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(u, &xv));
174e669de00SBarry Smith   }
175e669de00SBarry Smith 
1769566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
177e669de00SBarry Smith   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) {
178e669de00SBarry Smith     for (i = 0; i < ctx->nlg; i++) {
1799566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGDraw(ctx->lg[i]));
1809566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSave(ctx->lg[i]));
181e669de00SBarry Smith     }
182e669de00SBarry Smith   }
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184e669de00SBarry Smith }
185