xref: /petsc/src/dm/impls/network/networkview.c (revision df1a93fe4eb31800d2f94dca3e1569c95400adb9)
1d2fd5932SHong Zhang #include <petsc/private/dmnetworkimpl.h> /*I  "petscdmnetwork.h"  I*/
2d2fd5932SHong Zhang 
3*df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_CSV(DM dm, PetscViewer viewer)
4*df1a93feSDuncan Campbell {
5*df1a93feSDuncan Campbell   DM              dmcoords;
6*df1a93feSDuncan Campbell   PetscInt        nsubnets, i, subnet, nvertices, nedges, vertex, edge;
7*df1a93feSDuncan Campbell   PetscInt        vertexOffsets[2], globalEdgeVertices[2];
8*df1a93feSDuncan Campbell   PetscScalar     vertexCoords[2];
9*df1a93feSDuncan Campbell   const PetscInt *vertices, *edges, *edgeVertices;
10*df1a93feSDuncan Campbell   Vec             allVertexCoords;
11*df1a93feSDuncan Campbell   PetscMPIInt     rank;
12*df1a93feSDuncan Campbell   MPI_Comm        comm;
13*df1a93feSDuncan Campbell 
14*df1a93feSDuncan Campbell   PetscFunctionBegin;
15*df1a93feSDuncan Campbell   // Get the network containing coordinate information
16*df1a93feSDuncan Campbell   PetscCall(DMGetCoordinateDM(dm, &dmcoords));
17*df1a93feSDuncan Campbell   // Get the coordinate vector for the network
18*df1a93feSDuncan Campbell   PetscCall(DMGetCoordinatesLocal(dm, &allVertexCoords));
19*df1a93feSDuncan Campbell   // Get the MPI communicator and this process' rank
20*df1a93feSDuncan Campbell   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
21*df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_rank(comm, &rank));
22*df1a93feSDuncan Campbell   // Start synchronized printing
23*df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
24*df1a93feSDuncan Campbell 
25*df1a93feSDuncan Campbell   // Write the header
26*df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Type,Rank,ID,X,Y,Z,Name,Color\n"));
27*df1a93feSDuncan Campbell 
28*df1a93feSDuncan Campbell   // Iterate each subnetwork (Note: We need to get the global number of subnets apparently)
29*df1a93feSDuncan Campbell   PetscCall(DMNetworkGetNumSubNetworks(dm, NULL, &nsubnets));
30*df1a93feSDuncan Campbell   for (subnet = 0; subnet < nsubnets; subnet++) {
31*df1a93feSDuncan Campbell     // Get the subnetwork's vertices and edges
32*df1a93feSDuncan Campbell     PetscCall(DMNetworkGetSubnetwork(dm, subnet, &nvertices, &nedges, &vertices, &edges));
33*df1a93feSDuncan Campbell 
34*df1a93feSDuncan Campbell     // Write out each vertex
35*df1a93feSDuncan Campbell     for (i = 0; i < nvertices; i++) {
36*df1a93feSDuncan Campbell       vertex = vertices[i];
37*df1a93feSDuncan Campbell       // Get the offset into the coordinate vector for the vertex
38*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetLocalVecOffset(dmcoords, vertex, ALL_COMPONENTS, vertexOffsets));
39*df1a93feSDuncan Campbell       vertexOffsets[1] = vertexOffsets[0] + 1;
40*df1a93feSDuncan Campbell       // Remap vertex to the global value
41*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetGlobalVertexIndex(dm, vertex, &vertex));
42*df1a93feSDuncan Campbell       // Get the vertex position from the coordinate vector
43*df1a93feSDuncan Campbell       PetscCall(VecGetValues(allVertexCoords, 2, vertexOffsets, vertexCoords));
44*df1a93feSDuncan Campbell 
45*df1a93feSDuncan Campbell       // TODO: Determine vertex color/name
46*df1a93feSDuncan Campbell       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Node,%" PetscInt_FMT ",%" PetscInt_FMT ",%lf,%lf,0,%" PetscInt_FMT "\n", (PetscInt)rank, vertex, (double)PetscRealPart(vertexCoords[0]), (double)PetscRealPart(vertexCoords[1]), vertex));
47*df1a93feSDuncan Campbell     }
48*df1a93feSDuncan Campbell 
49*df1a93feSDuncan Campbell     // Write out each edge
50*df1a93feSDuncan Campbell     for (i = 0; i < nedges; i++) {
51*df1a93feSDuncan Campbell       edge = edges[i];
52*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetConnectedVertices(dm, edge, &edgeVertices));
53*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetGlobalVertexIndex(dm, edgeVertices[0], &globalEdgeVertices[0]));
54*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetGlobalVertexIndex(dm, edgeVertices[1], &globalEdgeVertices[1]));
55*df1a93feSDuncan Campbell       PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edge, &edge));
56*df1a93feSDuncan Campbell 
57*df1a93feSDuncan Campbell       // TODO: Determine edge color/name
58*df1a93feSDuncan Campbell       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Edge,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",0,%" PetscInt_FMT "\n", (PetscInt)rank, edge, globalEdgeVertices[0], globalEdgeVertices[1], edge));
59*df1a93feSDuncan Campbell     }
60*df1a93feSDuncan Campbell   }
61*df1a93feSDuncan Campbell   // End synchronized printing
62*df1a93feSDuncan Campbell   PetscCall(PetscViewerFlush(viewer));
63*df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
64*df1a93feSDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
65*df1a93feSDuncan Campbell }
66*df1a93feSDuncan Campbell 
67*df1a93feSDuncan Campbell #include <petscdraw.h>
68*df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_Matplotlib(DM dm, PetscViewer viewer)
69*df1a93feSDuncan Campbell {
70*df1a93feSDuncan Campbell   PetscMPIInt rank, size, rank2;
71*df1a93feSDuncan Campbell   MPI_Comm    comm;
72*df1a93feSDuncan Campbell   char        filename[PETSC_MAX_PATH_LEN + 1], proccall[PETSC_MAX_PATH_LEN + 500], scriptFile[PETSC_MAX_PATH_LEN + 1], streamBuffer[256];
73*df1a93feSDuncan Campbell   PetscViewer csvViewer;
74*df1a93feSDuncan Campbell   FILE       *processFile = NULL;
75*df1a93feSDuncan Campbell   PetscBool   isnull;
76*df1a93feSDuncan Campbell   PetscDraw   draw;
77*df1a93feSDuncan Campbell 
78*df1a93feSDuncan Campbell   PetscFunctionBegin;
79*df1a93feSDuncan Campbell   // Deal with the PetscDraw we are given
80*df1a93feSDuncan Campbell   PetscCall(PetscViewerDrawGetDraw(viewer, 1, &draw));
81*df1a93feSDuncan Campbell   PetscCall(PetscDrawIsNull(draw, &isnull));
82*df1a93feSDuncan Campbell   PetscCall(PetscDrawSetVisible(draw, PETSC_FALSE));
83*df1a93feSDuncan Campbell 
84*df1a93feSDuncan Campbell   // Clear the file name buffer so all communicated bytes are well-defined
85*df1a93feSDuncan Campbell   PetscCall(PetscMemzero(filename, sizeof(filename)));
86*df1a93feSDuncan Campbell 
87*df1a93feSDuncan Campbell   // Get the MPI communicator and this process' rank
88*df1a93feSDuncan Campbell   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
89*df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_rank(comm, &rank));
90*df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_size(comm, &size));
91*df1a93feSDuncan Campbell 
92*df1a93feSDuncan Campbell   // Generate and broadcast the temporary file name from rank 0
93*df1a93feSDuncan Campbell   if (rank == 0) {
94*df1a93feSDuncan Campbell #if defined(PETSC_HAVE_TMPNAM_S)
95*df1a93feSDuncan Campbell     // Acquire a temporary file to write to and open an ASCII/CSV viewer
96*df1a93feSDuncan Campbell     PetscCheck(tmpnam_s(filename, sizeof(filename)) == 0, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
97*df1a93feSDuncan Campbell #elif defined(PETSC_HAVE_MKSTEMP) && __STDC_VERSION__ > 199901L
98*df1a93feSDuncan Campbell     size_t numChars;
99*df1a93feSDuncan Campbell     // Same thing, but for POSIX systems on which tmpnam is deprecated
100*df1a93feSDuncan Campbell     // Note: Configure may detect mkstemp but it will not be defined if compiling for C99, so check additional defines to see if we can use it
101*df1a93feSDuncan Campbell     PetscCall(PetscStrcpy(filename, "/tmp/"));
102*df1a93feSDuncan Campbell     // Mkstemp requires us to explicitly specify part of the path, but some systems may not like putting files in /tmp/ so have an option for it
103*df1a93feSDuncan Campbell     PetscCall(PetscOptionsGetString(NULL, NULL, "-dmnetwork_view_tmpdir", filename, sizeof(filename), NULL));
104*df1a93feSDuncan Campbell     // Make sure the filename ends with a '/'
105*df1a93feSDuncan Campbell     PetscCall(PetscStrlen(filename, &numChars));
106*df1a93feSDuncan Campbell     if (filename[numChars - 1] != '/') {
107*df1a93feSDuncan Campbell       filename[numChars]     = '/';
108*df1a93feSDuncan Campbell       filename[numChars + 1] = 0;
109*df1a93feSDuncan Campbell     }
110*df1a93feSDuncan Campbell     // Perform the actual temporary file creation
111*df1a93feSDuncan Campbell     PetscCall(PetscStrcat(filename, "XXXXXX"));
112*df1a93feSDuncan Campbell     PetscCheck(mkstemp(filename) != -1, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
113*df1a93feSDuncan Campbell #else
114*df1a93feSDuncan Campbell     // Same thing, but for older C versions which don't have the safe form
115*df1a93feSDuncan Campbell     PetscCheck(tmpnam(filename) != NULL, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
116*df1a93feSDuncan Campbell #endif
117*df1a93feSDuncan Campbell     // Broadcast the filename to all other MPI ranks
118*df1a93feSDuncan Campbell     for (rank2 = 1; rank2 < size; rank2++) PetscCallMPI(MPI_Send(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, rank2, 0, comm));
119*df1a93feSDuncan Campbell   } else {
120*df1a93feSDuncan Campbell     // Receive the file name
121*df1a93feSDuncan Campbell     PetscCallMPI(MPI_Recv(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, 0, 0, comm, MPI_STATUS_IGNORE));
122*df1a93feSDuncan Campbell   }
123*df1a93feSDuncan Campbell 
124*df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, filename, &csvViewer));
125*df1a93feSDuncan Campbell   PetscCall(PetscViewerPushFormat(csvViewer, PETSC_VIEWER_ASCII_CSV));
126*df1a93feSDuncan Campbell 
127*df1a93feSDuncan Campbell   // Use the CSV viewer to write out the local network
128*df1a93feSDuncan Campbell   PetscCall(DMView_Network_CSV(dm, csvViewer));
129*df1a93feSDuncan Campbell 
130*df1a93feSDuncan Campbell   // Close the viewer
131*df1a93feSDuncan Campbell   PetscCall(PetscViewerDestroy(&csvViewer));
132*df1a93feSDuncan Campbell 
133*df1a93feSDuncan Campbell   // Get the value of $PETSC_DIR
134*df1a93feSDuncan Campbell   PetscCall(PetscStrreplace(PETSC_COMM_WORLD, "${PETSC_DIR}/share/petsc/bin/dmnetwork_view.py", scriptFile, sizeof(scriptFile)));
135*df1a93feSDuncan Campbell   PetscCall(PetscFixFilename(scriptFile, scriptFile));
136*df1a93feSDuncan Campbell   // Generate the system call for 'python3 $PETSC_DIR/share/petsc/dmnetwork_view.py <file>'
137*df1a93feSDuncan Campbell   PetscCall(PetscArrayzero(proccall, sizeof(proccall)));
138*df1a93feSDuncan Campbell   PetscCall(PetscSNPrintf(proccall, sizeof(proccall), "%s %s %s %s", PETSC_PYTHON_EXE, scriptFile, (isnull ? "-tx" : ""), filename));
139*df1a93feSDuncan Campbell 
140*df1a93feSDuncan Campbell #if defined(PETSC_HAVE_POPEN)
141*df1a93feSDuncan Campbell   // Perform the call to run the python script (Note: while this is called on all ranks POpen will only run on rank 0)
142*df1a93feSDuncan Campbell   PetscCall(PetscPOpen(PETSC_COMM_WORLD, NULL, proccall, "r", &processFile));
143*df1a93feSDuncan Campbell   if (processFile != NULL) {
144*df1a93feSDuncan Campbell     while (fgets(streamBuffer, sizeof(streamBuffer), processFile) != NULL) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s", streamBuffer));
145*df1a93feSDuncan Campbell   }
146*df1a93feSDuncan Campbell   PetscCall(PetscPClose(PETSC_COMM_WORLD, processFile));
147*df1a93feSDuncan Campbell #else
148*df1a93feSDuncan Campbell   // Same thing, but using the standard library for systems that don't have POpen/PClose (only run on rank 0)
149*df1a93feSDuncan Campbell   if (rank == 0) {
150*df1a93feSDuncan Campbell     PetscCheck(system(proccall) == 0, comm, PETSC_ERR_SYS, "Failed to call viewer script");
151*df1a93feSDuncan Campbell     // Barrier so that all ranks wait until the call completes
152*df1a93feSDuncan Campbell     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
153*df1a93feSDuncan Campbell   }
154*df1a93feSDuncan Campbell #endif
155*df1a93feSDuncan Campbell   // Clean up the temporary file we used using rank 0
156*df1a93feSDuncan Campbell   if (rank == 0) PetscCheck(remove(filename) == 0, comm, PETSC_ERR_SYS, "Failed to delete temporary file");
157*df1a93feSDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
158*df1a93feSDuncan Campbell }
159*df1a93feSDuncan Campbell 
160d2fd5932SHong Zhang PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
161d2fd5932SHong Zhang {
162*df1a93feSDuncan Campbell   PetscBool         iascii, isdraw;
163*df1a93feSDuncan Campbell   PetscViewerFormat format;
164d2fd5932SHong Zhang 
165d2fd5932SHong Zhang   PetscFunctionBegin;
166d2fd5932SHong Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
167d2fd5932SHong Zhang   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
168*df1a93feSDuncan Campbell   PetscCall(PetscViewerGetFormat(viewer, &format));
169*df1a93feSDuncan Campbell 
170*df1a93feSDuncan Campbell   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
171*df1a93feSDuncan Campbell   if (isdraw) {
172*df1a93feSDuncan Campbell     PetscCall(DMView_Network_Matplotlib(dm, viewer));
173*df1a93feSDuncan Campbell     PetscFunctionReturn(PETSC_SUCCESS);
174*df1a93feSDuncan Campbell   }
175*df1a93feSDuncan Campbell 
176d2fd5932SHong Zhang   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
177d2fd5932SHong Zhang   if (iascii) {
178d2fd5932SHong Zhang     const PetscInt *cone, *vtx, *edges;
179d2fd5932SHong Zhang     PetscInt        vfrom, vto, i, j, nv, ne, nsv, p, nsubnet;
180d2fd5932SHong Zhang     DM_Network     *network = (DM_Network *)dm->data;
181*df1a93feSDuncan Campbell     PetscMPIInt     rank;
182*df1a93feSDuncan Campbell 
183*df1a93feSDuncan Campbell     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
184*df1a93feSDuncan Campbell     if (format == PETSC_VIEWER_ASCII_CSV) {
185*df1a93feSDuncan Campbell       PetscCall(DMView_Network_CSV(dm, viewer));
186*df1a93feSDuncan Campbell       PetscFunctionReturn(PETSC_SUCCESS);
187*df1a93feSDuncan Campbell     }
188d2fd5932SHong Zhang 
189d2fd5932SHong Zhang     nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */
190*df1a93feSDuncan Campbell     if (!rank) {
191d2fd5932SHong Zhang       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  NSubnets: %" PetscInt_FMT "; NEdges: %" PetscInt_FMT "; NVertices: %" PetscInt_FMT "; NSharedVertices: %" PetscInt_FMT ".\n", nsubnet, network->cloneshared->NEdges, network->cloneshared->NVertices,
192d2fd5932SHong Zhang                             network->cloneshared->Nsvtx));
193d2fd5932SHong Zhang     }
194d2fd5932SHong Zhang 
195d2fd5932SHong Zhang     PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL));
196d2fd5932SHong Zhang     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
197d2fd5932SHong Zhang     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv));
198d2fd5932SHong Zhang 
199d2fd5932SHong Zhang     for (i = 0; i < nsubnet; i++) {
200d2fd5932SHong Zhang       PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges));
201d2fd5932SHong Zhang       if (ne) {
202d2fd5932SHong Zhang         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv));
203d2fd5932SHong Zhang         for (j = 0; j < ne; j++) {
204d2fd5932SHong Zhang           p = edges[j];
205d2fd5932SHong Zhang           PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone));
206d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom));
207d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto));
208d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p));
209d2fd5932SHong Zhang           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto));
210d2fd5932SHong Zhang         }
211d2fd5932SHong Zhang       }
212d2fd5932SHong Zhang     }
213d2fd5932SHong Zhang 
214d2fd5932SHong Zhang     /* Shared vertices */
215d2fd5932SHong Zhang     PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx));
216d2fd5932SHong Zhang     if (nsv) {
217d2fd5932SHong Zhang       PetscInt        gidx;
218d2fd5932SHong Zhang       PetscBool       ghost;
219d2fd5932SHong Zhang       const PetscInt *sv = NULL;
220d2fd5932SHong Zhang 
221d2fd5932SHong Zhang       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n"));
222d2fd5932SHong Zhang       for (i = 0; i < nsv; i++) {
223d2fd5932SHong Zhang         PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost));
224d2fd5932SHong Zhang         if (ghost) continue;
225d2fd5932SHong Zhang 
226d2fd5932SHong Zhang         PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv));
227d2fd5932SHong Zhang         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1]));
228d2fd5932SHong Zhang         for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1]));
229d2fd5932SHong Zhang       }
230d2fd5932SHong Zhang     }
231d2fd5932SHong Zhang     PetscCall(PetscViewerFlush(viewer));
232d2fd5932SHong Zhang     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
233d2fd5932SHong Zhang   } else PetscCheck(iascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
234d2fd5932SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
235d2fd5932SHong Zhang }
236