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