Actual source code: dmnetworkts.c
1: #include <petsc/private/dmnetworkimpl.h>
2: #include <petscts.h>
3: #include <petscdraw.h>
5: /*@C
6: TSMonitorLGCtxNetworkDestroy - Destroys line graph contexts that where created with `TSMonitorLGCtxNetworkCreate()`.
8: Collective
10: Input Parameter:
11: . ctx - the monitor context
13: Level: intermediate
15: .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkSolution()`
16: @*/
17: PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *ctx)
18: {
19: PetscInt i;
21: PetscFunctionBegin;
22: for (i = 0; i < (*ctx)->nlg; i++) PetscCall(PetscDrawLGDestroy(&(*ctx)->lg[i]));
23: PetscCall(PetscFree((*ctx)->lg));
24: PetscCall(PetscFree(*ctx));
25: PetscFunctionReturn(PETSC_SUCCESS);
26: }
28: PetscErrorCode TSMonitorLGCtxNetworkCreate(TS ts, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TSMonitorLGCtxNetwork *ctx)
29: {
30: PetscDraw draw;
31: MPI_Comm comm;
32: DM dm;
33: PetscInt i, Start, End, e, nvar;
35: PetscFunctionBegin;
36: PetscCall(TSGetDM(ts, &dm));
37: PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
38: PetscCall(PetscNew(ctx));
39: i = 0;
40: /* loop over edges counting number of line graphs needed */
41: PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
42: for (e = Start; e < End; e++) {
43: PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
44: if (!nvar) continue;
45: i++;
46: }
47: /* loop over vertices */
48: PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
49: for (e = Start; e < End; e++) {
50: PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
51: if (!nvar) continue;
52: i++;
53: }
54: (*ctx)->nlg = i;
55: PetscCall(PetscMalloc1(i, &(*ctx)->lg));
57: i = 0;
58: /* loop over edges creating all needed line graphs*/
59: PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
60: for (e = Start; e < End; e++) {
61: PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
62: if (!nvar) continue;
63: PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
64: PetscCall(PetscDrawSetFromOptions(draw));
65: PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
66: PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
67: PetscCall(PetscDrawDestroy(&draw));
68: i++;
69: }
70: /* loop over vertices */
71: PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
72: for (e = Start; e < End; e++) {
73: PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
74: if (!nvar) continue;
75: PetscCall(PetscDrawCreate(comm, host, label, x, y, m, n, &draw));
76: PetscCall(PetscDrawSetFromOptions(draw));
77: PetscCall(PetscDrawLGCreate(draw, nvar, &(*ctx)->lg[i]));
78: PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i]));
79: PetscCall(PetscDrawDestroy(&draw));
80: i++;
81: }
82: PetscCall(PetscDrawDestroy(&draw));
83: (*ctx)->howoften = howoften;
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: /*@C
88: TSMonitorLGCtxNetworkSolution - Monitors progress of the `TS` solvers for a `DMNETWORK` solution with one window for each vertex and each edge
90: Collective
92: Input Parameters:
93: + ts - the `TS` context
94: . step - current time-step
95: . ptime - current time
96: . u - current solution
97: - dctx - the `TSMonitorLGCtxNetwork` object that contains all the options for the monitoring, this is created with `TSMonitorLGCtxCreateNetwork()`
99: Options Database Key:
100: . -ts_monitor_lg_solution_variables - monitor solution variables
102: Level: intermediate
104: Note:
105: Each process in a parallel run displays its component solutions in a separate graphics window
107: .seealso: [](ch_ts), `TS`, `TSMonitorLGCtxNetworkDestroy()`
108: @*/
109: PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
110: {
111: TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx;
112: const PetscScalar *xv;
113: PetscScalar *yv;
114: PetscInt i, v, Start, End, offset, nvar, e;
115: TSConvergedReason reason;
116: DM dm;
117: Vec uv;
119: PetscFunctionBegin;
120: if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
121: if (!step) {
122: PetscDrawAxis axis;
124: for (i = 0; i < ctx->nlg; i++) {
125: PetscCall(PetscDrawLGGetAxis(ctx->lg[i], &axis));
126: PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
127: PetscCall(PetscDrawLGReset(ctx->lg[i]));
128: }
129: }
131: if (ctx->semilogy) {
132: PetscInt n, j;
134: PetscCall(VecDuplicate(u, &uv));
135: PetscCall(VecCopy(u, uv));
136: PetscCall(VecGetArray(uv, &yv));
137: PetscCall(VecGetLocalSize(uv, &n));
138: for (j = 0; j < n; j++) {
139: if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12;
140: else yv[j] = PetscLog10Real(PetscRealPart(yv[j]));
141: }
142: xv = yv;
143: } else {
144: PetscCall(VecGetArrayRead(u, &xv));
145: }
146: /* iterate over edges */
147: PetscCall(TSGetDM(ts, &dm));
148: i = 0;
149: PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
150: for (e = Start; e < End; e++) {
151: PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
152: if (!nvar) continue;
154: PetscCall(DMNetworkGetLocalVecOffset(dm, e, ALL_COMPONENTS, &offset));
155: PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
156: i++;
157: }
159: /* iterate over vertices */
160: PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
161: for (v = Start; v < End; v++) {
162: PetscCall(DMNetworkGetComponent(dm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
163: if (!nvar) continue;
165: PetscCall(DMNetworkGetLocalVecOffset(dm, v, ALL_COMPONENTS, &offset));
166: PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, (const PetscReal *)(xv + offset)));
167: i++;
168: }
169: if (ctx->semilogy) {
170: PetscCall(VecRestoreArray(uv, &yv));
171: PetscCall(VecDestroy(&uv));
172: } else {
173: PetscCall(VecRestoreArrayRead(u, &xv));
174: }
176: PetscCall(TSGetConvergedReason(ts, &reason));
177: if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) {
178: for (i = 0; i < ctx->nlg; i++) {
179: PetscCall(PetscDrawLGDraw(ctx->lg[i]));
180: PetscCall(PetscDrawLGSave(ctx->lg[i]));
181: }
182: }
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }