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: #if defined(PETSC_USE_COMPLEX)
112:   PetscFunctionBegin;
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: #else
115:   TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx;
116:   const PetscScalar    *xv;
117:   PetscScalar          *yv;
118:   PetscInt              i, v, Start, End, offset, nvar, e;
119:   TSConvergedReason     reason;
120:   DM                    dm;
121:   Vec                   uv;

123:   PetscFunctionBegin;
124:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
125:   if (!step) {
126:     PetscDrawAxis axis;

128:     for (i = 0; i < ctx->nlg; i++) {
129:       PetscCall(PetscDrawLGGetAxis(ctx->lg[i], &axis));
130:       PetscCall(PetscDrawAxisSetLabels(axis, "Solution as function of time", "Time", "Solution"));
131:       PetscCall(PetscDrawLGReset(ctx->lg[i]));
132:     }
133:   }

135:   if (ctx->semilogy) {
136:     PetscInt n, j;

138:     PetscCall(VecDuplicate(u, &uv));
139:     PetscCall(VecCopy(u, uv));
140:     PetscCall(VecGetArray(uv, &yv));
141:     PetscCall(VecGetLocalSize(uv, &n));
142:     for (j = 0; j < n; j++) {
143:       if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12;
144:       else yv[j] = PetscLog10Real(PetscRealPart(yv[j]));
145:     }
146:     xv = yv;
147:   } else {
148:     PetscCall(VecGetArrayRead(u, &xv));
149:   }
150:   /* iterate over edges */
151:   PetscCall(TSGetDM(ts, &dm));
152:   i = 0;
153:   PetscCall(DMNetworkGetEdgeRange(dm, &Start, &End));
154:   for (e = Start; e < End; e++) {
155:     PetscCall(DMNetworkGetComponent(dm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
156:     if (!nvar) continue;

158:     PetscCall(DMNetworkGetLocalVecOffset(dm, e, ALL_COMPONENTS, &offset));
159:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, xv + offset));
160:     i++;
161:   }

163:   /* iterate over vertices */
164:   PetscCall(DMNetworkGetVertexRange(dm, &Start, &End));
165:   for (v = Start; v < End; v++) {
166:     PetscCall(DMNetworkGetComponent(dm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
167:     if (!nvar) continue;

169:     PetscCall(DMNetworkGetLocalVecOffset(dm, v, ALL_COMPONENTS, &offset));
170:     PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i], ptime, xv + offset));
171:     i++;
172:   }
173:   if (ctx->semilogy) {
174:     PetscCall(VecRestoreArray(uv, &yv));
175:     PetscCall(VecDestroy(&uv));
176:   } else {
177:     PetscCall(VecRestoreArrayRead(u, &xv));
178:   }

180:   PetscCall(TSGetConvergedReason(ts, &reason));
181:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) {
182:     for (i = 0; i < ctx->nlg; i++) {
183:       PetscCall(PetscDrawLGDraw(ctx->lg[i]));
184:       PetscCall(PetscDrawLGSave(ctx->lg[i]));
185:     }
186:   }
187:   PetscFunctionReturn(PETSC_SUCCESS);
188: #endif
189: }