xref: /petsc/src/ts/utils/dmnetworkts.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1e669de00SBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmplex.h" I*/
2e669de00SBarry Smith #include <petscts.h>
3e669de00SBarry Smith #include <petscdraw.h>
4e669de00SBarry Smith 
5e669de00SBarry Smith /*
64500feddSGetnet Betrie    TSMonitorLGCtxDestroy - Destroys  line graph contexts that where created with TSMonitorLGCtxNetworkCreate().
7e669de00SBarry Smith 
8e669de00SBarry Smith    Collective on TSMonitorLGCtx_Network
9e669de00SBarry Smith 
10e669de00SBarry Smith    Input Parameter:
11e669de00SBarry Smith .  ctx - the monitor context
12e669de00SBarry Smith 
13e669de00SBarry Smith */
149371c9d4SSatish Balay PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *ctx) {
15e669de00SBarry Smith   PetscInt i;
16e669de00SBarry Smith 
17e669de00SBarry Smith   PetscFunctionBegin;
18*48a46eb9SPierre Jolivet   for (i = 0; i < (*ctx)->nlg; i++) PetscCall(PetscDrawLGDestroy(&(*ctx)->lg[i]));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->lg));
209566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
21e669de00SBarry Smith   PetscFunctionReturn(0);
22e669de00SBarry Smith }
23e669de00SBarry Smith 
249371c9d4SSatish Balay PetscErrorCode TSMonitorLGCtxNetworkCreate(TS ts, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtxNetwork *ctx) {
25e669de00SBarry Smith   PetscDraw draw;
26e669de00SBarry Smith   MPI_Comm  comm;
27e669de00SBarry Smith   DM        dm;
28e669de00SBarry Smith   PetscInt  i, Start, End, e, nvar;
29e669de00SBarry Smith 
30e669de00SBarry Smith   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
339566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
34e669de00SBarry Smith   i = 0;
35e669de00SBarry Smith   /* loop over edges counting number of line graphs needed */
369566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
37e669de00SBarry Smith   for (e = Start; e < End; e++) {
389566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
39e669de00SBarry Smith     if (!nvar) continue;
40e669de00SBarry Smith     i++;
41e669de00SBarry Smith   }
42e669de00SBarry Smith   /* loop over vertices */
439566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
44e669de00SBarry Smith   for (e = Start; e < End; e++) {
459566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
46e669de00SBarry Smith     if (!nvar) continue;
47e669de00SBarry Smith     i++;
48e669de00SBarry Smith   }
49e669de00SBarry Smith   (*ctx)->nlg = i;
509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(i, &(*ctx)->lg));
51e669de00SBarry Smith 
52e669de00SBarry Smith   i = 0;
53e669de00SBarry Smith   /* loop over edges creating all needed line graphs*/
549566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
55e669de00SBarry Smith   for (e = Start; e < End; e++) {
569566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
57e669de00SBarry Smith     if (!nvar) continue;
589566063dSJacob Faibussowitsch     PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
599566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetFromOptions(draw));
609566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
619566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
629566063dSJacob Faibussowitsch     PetscCall(PetscDrawDestroy(&draw));
63e669de00SBarry Smith     i++;
64e669de00SBarry Smith   }
65e669de00SBarry Smith   /* loop over vertices */
669566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
67e669de00SBarry Smith   for (e = Start; e < End; e++) {
689566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
69e669de00SBarry Smith     if (!nvar) continue;
709566063dSJacob Faibussowitsch     PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
719566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetFromOptions(draw));
729566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
739566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
749566063dSJacob Faibussowitsch     PetscCall(PetscDrawDestroy(&draw));
75e669de00SBarry Smith     i++;
76e669de00SBarry Smith   }
779566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&draw));
78e669de00SBarry Smith   (*ctx)->howoften = howoften;
79e669de00SBarry Smith   PetscFunctionReturn(0);
80e669de00SBarry Smith }
81e669de00SBarry Smith 
82e669de00SBarry Smith /*
83e669de00SBarry Smith    TSMonitorLGCtxNetworkSolution - Monitors progress of the TS solvers for a DMNetwork solution with one window for each vertex and each edge
84e669de00SBarry Smith 
85e669de00SBarry Smith    Collective on TS
86e669de00SBarry Smith 
87e669de00SBarry Smith    Input Parameters:
88e669de00SBarry Smith +  ts - the TS context
89e669de00SBarry Smith .  step - current time-step
90e669de00SBarry Smith .  ptime - current time
91e669de00SBarry Smith .  u - current solution
924500feddSGetnet Betrie -  dctx - the TSMonitorLGCtxNetwork object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreateNetwork()
93e669de00SBarry Smith 
94e669de00SBarry Smith    Options Database:
95e669de00SBarry Smith .   -ts_monitor_lg_solution_variables
96e669de00SBarry Smith 
97e669de00SBarry Smith    Level: intermediate
98e669de00SBarry Smith 
99e669de00SBarry Smith    Notes:
100e669de00SBarry Smith     Each process in a parallel run displays its component solutions in a separate window
101e669de00SBarry Smith 
102e669de00SBarry Smith */
1039371c9d4SSatish Balay PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx) {
104e669de00SBarry Smith   TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx;
105e669de00SBarry Smith   const PetscScalar    *xv;
106e669de00SBarry Smith   PetscScalar          *yv;
107e669de00SBarry Smith   PetscInt              i, v, Start, End, offset, nvar, e;
108e669de00SBarry Smith   TSConvergedReason     reason;
109e669de00SBarry Smith   DM                    dm;
110e669de00SBarry Smith   Vec                   uv;
111e669de00SBarry Smith 
112e669de00SBarry Smith   PetscFunctionBegin;
113e669de00SBarry Smith   if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */
114e669de00SBarry Smith   if (!step) {
115e669de00SBarry Smith     PetscDrawAxis axis;
116e669de00SBarry Smith 
117e669de00SBarry Smith     for (i = 0; i < ctx->nlg; i++) {
1189566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGGetAxis(ctx->lg[i], &axis));
1199566063dSJacob Faibussowitsch       PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
1209566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGReset(ctx->lg[i]));
121e669de00SBarry Smith     }
122e669de00SBarry Smith   }
123e669de00SBarry Smith 
124e669de00SBarry Smith   if (ctx->semilogy) {
125e669de00SBarry Smith     PetscInt n, j;
126e669de00SBarry Smith 
1279566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uv));
1289566063dSJacob Faibussowitsch     PetscCall(VecCopy(u, uv));
1299566063dSJacob Faibussowitsch     PetscCall(VecGetArray(uv, &yv));
1309566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(uv, &n));
131e669de00SBarry Smith     for (j = 0; j < n; j++) {
1324500feddSGetnet Betrie       if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12;
1334e936a74SSatish Balay       else yv[j] = PetscLog10Real(PetscRealPart(yv[j]));
134e669de00SBarry Smith     }
135e669de00SBarry Smith     xv = yv;
136e669de00SBarry Smith   } else {
1379566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(u, &xv));
138e669de00SBarry Smith   }
139e669de00SBarry Smith   /* iterate over edges */
1409566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
141e669de00SBarry Smith   i = 0;
1429566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
143e669de00SBarry Smith   for (e = Start; e < End; e++) {
1449566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
145e669de00SBarry Smith     if (!nvar) continue;
146e669de00SBarry Smith 
1479566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(dm, e, ALL_COMPONENTS, &offset));
1489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
149e669de00SBarry Smith     i++;
150e669de00SBarry Smith   }
151e669de00SBarry Smith 
152e669de00SBarry Smith   /* iterate over vertices */
1539566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
154e669de00SBarry Smith   for (v = Start; v < End; v++) {
1559566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(dm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
156e669de00SBarry Smith     if (!nvar) continue;
157e669de00SBarry Smith 
1589566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(dm, v, ALL_COMPONENTS, &offset));
1599566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
160e669de00SBarry Smith     i++;
161e669de00SBarry Smith   }
162e669de00SBarry Smith   if (ctx->semilogy) {
1639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(uv, &yv));
1649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uv));
165e669de00SBarry Smith   } else {
1669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(u, &xv));
167e669de00SBarry Smith   }
168e669de00SBarry Smith 
1699566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
170e669de00SBarry Smith   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) {
171e669de00SBarry Smith     for (i = 0; i < ctx->nlg; i++) {
1729566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGDraw(ctx->lg[i]));
1739566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGSave(ctx->lg[i]));
174e669de00SBarry Smith     }
175e669de00SBarry Smith   }
176e669de00SBarry Smith   PetscFunctionReturn(0);
177e669de00SBarry Smith }
178