Actual source code: networkmonitor.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <petscdmnetwork.h> /*I  "petscdmnetwork.h"  I*/
  2: #include <petscdraw.h>

  6: /*@
  7:   DMNetworkMonitorCreate - Creates a network monitor context

  9:   Collective on MPI_Comm

 11:   Input Parameters:
 12: . network - network to monitor

 14:   Output Parameters:
 15: . Monitorptr - Location to put network monitor context

 17:   Level: intermediate

 19: .seealso: DMNetworkMonitorDestroy(), DMNetworkMonitorAdd()
 20: @*/
 21: PetscErrorCode DMNetworkMonitorCreate(DM network,DMNetworkMonitor *monitorptr)
 22: {
 23:   PetscErrorCode   ierr;
 24:   DMNetworkMonitor monitor;
 25:   MPI_Comm         comm;
 26:   PetscMPIInt      size;

 29:   PetscObjectGetComm((PetscObject)network,&comm);
 30:   MPI_Comm_size(comm, &size);
 31:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Parallel DMNetworkMonitor is not supported yet");

 33:   PetscMalloc1(1,&monitor);
 34:   monitor->comm      = comm;
 35:   monitor->network   = network;
 36:   monitor->firstnode = PETSC_NULL;

 38:   *monitorptr = monitor;
 39:   return(0);
 40: }

 44: /*@
 45:   DMNetworkMonitorDestroy - Destroys a network monitor and all associated viewers

 47:   Collective on DMNetworkMonitor

 49:   Input Parameters:
 50: . monitor - monitor to destroy

 52:   Level: intermediate

 54: .seealso: DMNetworkMonitorCreate, DMNetworkMonitorAdd
 55: @*/
 56: PetscErrorCode DMNetworkMonitorDestroy(DMNetworkMonitor *monitor)
 57: {

 61:   while ((*monitor)->firstnode) {
 62:     DMNetworkMonitorPop(*monitor);
 63:   }

 65:   PetscFree(*monitor);
 66:   return(0);
 67: }

 71: /*@
 72:   DMNetworkMonitorPop - Removes the most recently added viewer

 74:   Collective on DMNetworkMonitor

 76:   Input Parameters:
 77: . monitor - the monitor

 79:   Level: intermediate

 81: .seealso: DMNetworkMonitorCreate(), DMNetworkMonitorDestroy()
 82: @*/
 83: PetscErrorCode DMNetworkMonitorPop(DMNetworkMonitor monitor)
 84: {
 85:   PetscErrorCode       ierr;
 86:   DMNetworkMonitorList node;

 89:   if (monitor->firstnode) {
 90:     /* Update links */
 91:     node = monitor->firstnode;
 92:     monitor->firstnode = node->next;

 94:     /* Free list node */
 95:     PetscViewerDestroy(&(node->viewer));
 96:     VecDestroy(&(node->v));
 97:     PetscFree(node);
 98:   }
 99:   return(0);
100: }

104: /*@
105:   DMNetworkMonitorAdd - Adds a new viewer to monitor

107:   Collective on DMNetworkMonitor

109:   Input Parameters:
110: + monitor - the monitor
111: . name - name of viewer
112: . element - vertex / edge number
113: . nodes - number of nodes
114: . start - variable starting offset
115: . blocksize - variable blocksize
116: . ymin - ymin for viewer
117: . ymax - ymax for viewer
118: - hold - determines if plot limits should be held

120:   Level: intermediate

122:   Notes:
123:   This is written to be independent of the semantics associated to the variables
124:   at a given network vertex / edge.

126:   Precisely, the parameters nodes, start and blocksize allow you to select a general
127:   strided subarray of the variables to monitor.

129: .seealso: DMNetworkMonitorCreate(), DMNetworkMonitorDestroy() 
130: @*/
131: PetscErrorCode DMNetworkMonitorAdd(DMNetworkMonitor monitor,const char *name,PetscInt element,PetscInt nodes,PetscInt start,PetscInt blocksize,PetscReal ymin,PetscReal ymax,PetscBool hold)
132: {
133:   PetscErrorCode       ierr;
134:   PetscDrawLG          drawlg;
135:   PetscDrawAxis        axis;
136:   PetscMPIInt          rank, size;
137:   DMNetworkMonitorList node;
138:   char                 titleBuffer[64];
139:   PetscInt             vStart,vEnd,eStart,eEnd;

142:   MPI_Comm_rank(monitor->comm, &rank);
143:   MPI_Comm_size(monitor->comm, &size);

145:   DMNetworkGetVertexRange(monitor->network, &vStart, &vEnd);
146:   DMNetworkGetEdgeRange(monitor->network, &eStart, &eEnd);

148:   /* Make window title */
149:   if (vStart <= element && element < vEnd) {
150:     PetscSNPrintf(titleBuffer, 64, "%s @ vertex %d [%d / %d]", name, element - vStart, rank, size-1);
151:   } else if (eStart <= element && element < eEnd) {
152:     PetscSNPrintf(titleBuffer, 64, "%s @ edge %d [%d / %d]", name, element - eStart, rank, size-1);
153:   } else {
154:     /* vertex / edge is not on local machine, so skip! */
155:     return(0);
156:   }

158:   PetscMalloc1(1, &node);

160:   /* Setup viewer. */
161:   PetscViewerDrawOpen(monitor->comm, PETSC_NULL, titleBuffer, PETSC_DECIDE, PETSC_DECIDE, PETSC_DRAW_QUARTER_SIZE, PETSC_DRAW_QUARTER_SIZE, &(node->viewer));
162:   PetscViewerPushFormat(node->viewer, PETSC_VIEWER_DRAW_LG);
163:   PetscViewerDrawGetDrawLG(node->viewer, 0, &drawlg);
164:   PetscDrawLGGetAxis(drawlg, &axis);
165:   PetscDrawAxisSetLimits(axis, 0, nodes-1, ymin, ymax);
166:   PetscDrawAxisSetHoldLimits(axis, hold);

168:   /* Setup vector storage for drawing. */
169:   VecCreateSeq(PETSC_COMM_SELF, nodes, &(node->v));

171:   node->element   = element;
172:   node->nodes     = nodes;
173:   node->start     = start;
174:   node->blocksize = blocksize;

176:   node->next         = monitor->firstnode;
177:   monitor->firstnode = node;
178:   return(0);
179: }

183: /*@
184:   DMNetworkMonitorView - Monitor function for TSMonitorSet.

186:   Collectiveon DMNetworkMonitor

188:   Input Parameters:
189: + monitor - DMNetworkMonitor object
190: - x - TS solution vector

192:   Level: intermediate

194: .seealso: DMNetworkMonitorCreate(), DMNetworkMonitorDestroy(), DMNetworkMonitorAdd()
195: @*/
196: PetscErrorCode DMNetworkMonitorView(DMNetworkMonitor monitor,Vec x)
197: {
198:   PetscErrorCode      ierr;
199:   PetscInt            varoffset,i,start;
200:   const PetscScalar   *xx;
201:   PetscScalar         *vv;
202:   DMNetworkMonitorList node;

205:   VecGetArrayRead(x, &xx);
206:   for (node = monitor->firstnode; node; node = node->next) {
207:     DMNetworkGetVariableGlobalOffset(monitor->network, node->element, &varoffset);
208:     VecGetArray(node->v, &vv);
209:     start = varoffset + node->start;
210:     for (i = 0; i < node->nodes; i++) {
211:       vv[i] = xx[start+i*node->blocksize];
212:     }
213:     VecRestoreArray(node->v, &vv);
214:     VecView(node->v, node->viewer);
215:   }
216:   VecRestoreArrayRead(x, &xx);
217:   return(0);
218: }