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