xref: /petsc/src/dm/impls/network/networkview.c (revision baca6076736e385493fd849a2765fd1676cb8d6d)
112e0ef65SDuncan Campbell #include <petscconf.h>
212e0ef65SDuncan Campbell // We need to define this ahead of any other includes to make sure mkstemp is actually defined
312e0ef65SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP)
412e0ef65SDuncan Campbell   #define _XOPEN_SOURCE 600
512e0ef65SDuncan Campbell #endif
65f25b224SDuncan Campbell #include "petsc/private/petscimpl.h"
75f25b224SDuncan Campbell #include "petscerror.h"
85f25b224SDuncan Campbell #include "petscis.h"
95f25b224SDuncan Campbell #include "petscstring.h"
105f25b224SDuncan Campbell #include "petscsys.h"
115f25b224SDuncan Campbell #include "petscsystypes.h"
12d2fd5932SHong Zhang #include <petsc/private/dmnetworkimpl.h> /*I  "petscdmnetwork.h"  I*/
1312e0ef65SDuncan Campbell #include <petscdraw.h>
14d2fd5932SHong Zhang 
15df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_CSV(DM dm, PetscViewer viewer)
16df1a93feSDuncan Campbell {
17df1a93feSDuncan Campbell   DM              dmcoords;
187cd471e7SHong Zhang   PetscInt        nsubnets, i, subnet, nvertices, nedges, vertex, edge, gidx, ncomp;
19df1a93feSDuncan Campbell   PetscInt        vertexOffsets[2], globalEdgeVertices[2];
207cd471e7SHong Zhang   PetscScalar     vertexCoords[2], *color_ptr, color;
21df1a93feSDuncan Campbell   const PetscInt *vertices, *edges, *edgeVertices;
22df1a93feSDuncan Campbell   Vec             allVertexCoords;
23df1a93feSDuncan Campbell   PetscMPIInt     rank;
24df1a93feSDuncan Campbell   MPI_Comm        comm;
25df1a93feSDuncan Campbell 
26df1a93feSDuncan Campbell   PetscFunctionBegin;
277cd471e7SHong Zhang   // Get the coordinate information from dmcoords
287cd471e7SHong Zhang   PetscCheck(dm->coordinates[0].dm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "CoordinateDM not created");
29df1a93feSDuncan Campbell   PetscCall(DMGetCoordinateDM(dm, &dmcoords));
307cd471e7SHong Zhang 
317cd471e7SHong Zhang   PetscCall(DMGetCoordinateDim(dmcoords, &i));
32*baca6076SPierre Jolivet   PetscCheck(i == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim %" PetscInt_FMT " != 2 is not supported yet", i);
337cd471e7SHong Zhang 
347cd471e7SHong Zhang   // Get the coordinate vector from dm
35df1a93feSDuncan Campbell   PetscCall(DMGetCoordinatesLocal(dm, &allVertexCoords));
367cd471e7SHong Zhang 
37df1a93feSDuncan Campbell   // Get the MPI communicator and this process' rank
387cd471e7SHong Zhang   PetscCall(PetscObjectGetComm((PetscObject)dmcoords, &comm));
39df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_rank(comm, &rank));
407cd471e7SHong Zhang 
41df1a93feSDuncan Campbell   // Start synchronized printing
42df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
43df1a93feSDuncan Campbell 
44df1a93feSDuncan Campbell   // Write the header
45df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Type,Rank,ID,X,Y,Z,Name,Color\n"));
46df1a93feSDuncan Campbell 
47df1a93feSDuncan Campbell   // Iterate each subnetwork (Note: We need to get the global number of subnets apparently)
487cd471e7SHong Zhang   PetscCall(DMNetworkGetNumSubNetworks(dmcoords, NULL, &nsubnets));
49df1a93feSDuncan Campbell   for (subnet = 0; subnet < nsubnets; subnet++) {
50df1a93feSDuncan Campbell     // Get the subnetwork's vertices and edges
517cd471e7SHong Zhang     PetscCall(DMNetworkGetSubnetwork(dmcoords, subnet, &nvertices, &nedges, &vertices, &edges));
52df1a93feSDuncan Campbell 
53df1a93feSDuncan Campbell     // Write out each vertex
54df1a93feSDuncan Campbell     for (i = 0; i < nvertices; i++) {
55df1a93feSDuncan Campbell       vertex = vertices[i];
567cd471e7SHong Zhang 
57df1a93feSDuncan Campbell       // Get the offset into the coordinate vector for the vertex
58df1a93feSDuncan Campbell       PetscCall(DMNetworkGetLocalVecOffset(dmcoords, vertex, ALL_COMPONENTS, vertexOffsets));
59df1a93feSDuncan Campbell       vertexOffsets[1] = vertexOffsets[0] + 1;
60df1a93feSDuncan Campbell       // Remap vertex to the global value
617cd471e7SHong Zhang       PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, vertex, &gidx));
62df1a93feSDuncan Campbell       // Get the vertex position from the coordinate vector
63df1a93feSDuncan Campbell       PetscCall(VecGetValues(allVertexCoords, 2, vertexOffsets, vertexCoords));
64df1a93feSDuncan Campbell 
657cd471e7SHong Zhang       // Get vertex color; TODO: name
667cd471e7SHong Zhang       PetscCall(DMNetworkGetNumComponents(dmcoords, vertex, &ncomp));
677cd471e7SHong Zhang       PetscCheck(ncomp <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "num of components %" PetscInt_FMT " must be <= 1", ncomp);
687cd471e7SHong Zhang       color = 0.0;
697cd471e7SHong Zhang       if (ncomp == 1) {
707cd471e7SHong Zhang         PetscCall(DMNetworkGetComponent(dmcoords, vertex, 0, NULL, (void **)&color_ptr, NULL));
717cd471e7SHong Zhang         color = *color_ptr;
727cd471e7SHong Zhang       }
737cd471e7SHong Zhang       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Node,%" PetscInt_FMT ",%" PetscInt_FMT ",%lf,%lf,0,%" PetscInt_FMT ",%lf\n", (PetscInt)rank, gidx, (double)PetscRealPart(vertexCoords[0]), (double)PetscRealPart(vertexCoords[1]), gidx, (double)PetscRealPart(color)));
74df1a93feSDuncan Campbell     }
75df1a93feSDuncan Campbell 
76df1a93feSDuncan Campbell     // Write out each edge
77df1a93feSDuncan Campbell     for (i = 0; i < nedges; i++) {
78df1a93feSDuncan Campbell       edge = edges[i];
797cd471e7SHong Zhang       PetscCall(DMNetworkGetConnectedVertices(dmcoords, edge, &edgeVertices));
807cd471e7SHong Zhang       PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[0], &globalEdgeVertices[0]));
817cd471e7SHong Zhang       PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[1], &globalEdgeVertices[1]));
827cd471e7SHong Zhang       PetscCall(DMNetworkGetGlobalEdgeIndex(dmcoords, edge, &edge));
83df1a93feSDuncan Campbell 
84df1a93feSDuncan Campbell       // TODO: Determine edge color/name
85df1a93feSDuncan 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));
86df1a93feSDuncan Campbell     }
87df1a93feSDuncan Campbell   }
88df1a93feSDuncan Campbell   // End synchronized printing
89df1a93feSDuncan Campbell   PetscCall(PetscViewerFlush(viewer));
90df1a93feSDuncan Campbell   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
91df1a93feSDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
92df1a93feSDuncan Campbell }
93df1a93feSDuncan Campbell 
94df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_Matplotlib(DM dm, PetscViewer viewer)
95df1a93feSDuncan Campbell {
96cd2bb8e3SDuncan Campbell   PetscMPIInt rank, size;
97df1a93feSDuncan Campbell   MPI_Comm    comm;
985039956bSHong Zhang   char        filename[PETSC_MAX_PATH_LEN + 1], options[512], proccall[PETSC_MAX_PATH_LEN + 512], scriptFile[PETSC_MAX_PATH_LEN + 1], buffer[256], buffer2[256];
99df1a93feSDuncan Campbell   PetscViewer csvViewer;
100df1a93feSDuncan Campbell   FILE       *processFile = NULL;
1015039956bSHong Zhang   PetscBool   isnull, optionShowRanks = PETSC_FALSE, optionRankIsSet = PETSC_FALSE, showNoNodes = PETSC_FALSE, showNoNumbering = PETSC_FALSE, optionShowVertices = PETSC_FALSE, optionViewPadding = PETSC_FALSE;
102df1a93feSDuncan Campbell   PetscDraw   draw;
1035f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
1045039956bSHong Zhang   PetscReal   drawPause, viewPadding = 1.0;
1055f25b224SDuncan Campbell   PetscInt    i;
106cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP)
107cd2bb8e3SDuncan Campbell   PetscBool isSharedTmp;
108cd2bb8e3SDuncan Campbell #endif
109df1a93feSDuncan Campbell 
110df1a93feSDuncan Campbell   PetscFunctionBegin;
111df1a93feSDuncan Campbell   // Deal with the PetscDraw we are given
112df1a93feSDuncan Campbell   PetscCall(PetscViewerDrawGetDraw(viewer, 1, &draw));
113df1a93feSDuncan Campbell   PetscCall(PetscDrawIsNull(draw, &isnull));
114df1a93feSDuncan Campbell   PetscCall(PetscDrawSetVisible(draw, PETSC_FALSE));
115df1a93feSDuncan Campbell 
116df1a93feSDuncan Campbell   // Clear the file name buffer so all communicated bytes are well-defined
117df1a93feSDuncan Campbell   PetscCall(PetscMemzero(filename, sizeof(filename)));
118df1a93feSDuncan Campbell 
119df1a93feSDuncan Campbell   // Get the MPI communicator and this process' rank
120df1a93feSDuncan Campbell   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
121df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_rank(comm, &rank));
122df1a93feSDuncan Campbell   PetscCallMPI(MPI_Comm_size(comm, &size));
123df1a93feSDuncan Campbell 
124cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP)
125cd2bb8e3SDuncan Campbell   // Get if the temporary directory is shared
126cd2bb8e3SDuncan Campbell   // Note: This must be done collectively on every rank, it cannot be done on a single rank
127cd2bb8e3SDuncan Campbell   PetscCall(PetscSharedTmp(comm, &isSharedTmp));
128cd2bb8e3SDuncan Campbell #endif
129cd2bb8e3SDuncan Campbell 
1305f25b224SDuncan Campbell   /* Process Options */
1315f25b224SDuncan Campbell   optionShowRanks = network->vieweroptions.showallranks;
1325f25b224SDuncan Campbell   showNoNodes     = network->vieweroptions.shownovertices;
1335f25b224SDuncan Campbell   showNoNumbering = network->vieweroptions.shownonumbering;
1345f25b224SDuncan Campbell 
1355f25b224SDuncan Campbell   /*
1365f25b224SDuncan Campbell     TODO:  if the option -dmnetwork_view_tmpdir can be moved up here that would be good as well.
1375f25b224SDuncan Campbell   */
1385f25b224SDuncan Campbell   PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "MatPlotLib PetscViewer DMNetwork Options", "PetscViewer");
1395f25b224SDuncan Campbell   PetscCall(PetscOptionsBool("-dmnetwork_view_all_ranks", "View all ranks in the DMNetwork", NULL, optionShowRanks, &optionShowRanks, NULL));
1405f25b224SDuncan Campbell   PetscCall(PetscOptionsString("-dmnetwork_view_rank_range", "Set of ranks to view the DMNetwork on", NULL, buffer, buffer, sizeof(buffer), &optionRankIsSet));
1415f25b224SDuncan Campbell   PetscCall(PetscOptionsBool("-dmnetwork_view_no_vertices", "Do not view vertices", NULL, showNoNodes, &showNoNodes, NULL));
1425f25b224SDuncan Campbell   PetscCall(PetscOptionsBool("-dmnetwork_view_no_numbering", "Do not view edge and vertex numbering", NULL, showNoNumbering, &showNoNumbering, NULL));
1435039956bSHong Zhang   PetscCall(PetscOptionsString("-dmnetwork_view_zoomin_vertices", "Focus the view on the given set of vertices", NULL, buffer2, buffer2, sizeof(buffer2), &optionShowVertices));
1445039956bSHong Zhang   PetscCall(PetscOptionsReal("-dmnetwork_view_zoomin_vertices_padding", "Set the padding when viewing specific vertices", NULL, viewPadding, &viewPadding, &optionViewPadding));
1455f25b224SDuncan Campbell   PetscOptionsEnd();
1465f25b224SDuncan Campbell 
147df1a93feSDuncan Campbell   // Generate and broadcast the temporary file name from rank 0
148df1a93feSDuncan Campbell   if (rank == 0) {
149df1a93feSDuncan Campbell #if defined(PETSC_HAVE_TMPNAM_S)
150df1a93feSDuncan Campbell     // Acquire a temporary file to write to and open an ASCII/CSV viewer
151df1a93feSDuncan Campbell     PetscCheck(tmpnam_s(filename, sizeof(filename)) == 0, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
15212e0ef65SDuncan Campbell #elif defined(PETSC_HAVE_MKSTEMP)
153cd2bb8e3SDuncan Campbell     PetscBool isTmpOverridden;
154df1a93feSDuncan Campbell     size_t    numChars;
155df1a93feSDuncan Campbell     // Same thing, but for POSIX systems on which tmpnam is deprecated
156df1a93feSDuncan 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
157df1a93feSDuncan 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
1586e4289a0SDuncan Campbell     PetscCall(PetscOptionsGetString(NULL, NULL, "-dmnetwork_view_tmpdir", filename, sizeof(filename), &isTmpOverridden));
1596e4289a0SDuncan Campbell     // If not specified by option try using a shared tmp on the system
1606e4289a0SDuncan Campbell     if (!isTmpOverridden) {
1616e4289a0SDuncan Campbell       // Validate that if tmp is not overridden it is at least shared
162cd2bb8e3SDuncan 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");
163cd2bb8e3SDuncan Campbell       PetscCall(PetscGetTmp(PETSC_COMM_SELF, filename, sizeof(filename)));
164cd2bb8e3SDuncan Campbell     }
165df1a93feSDuncan Campbell     // Make sure the filename ends with a '/'
166df1a93feSDuncan Campbell     PetscCall(PetscStrlen(filename, &numChars));
167df1a93feSDuncan Campbell     if (filename[numChars - 1] != '/') {
168df1a93feSDuncan Campbell       filename[numChars]     = '/';
169df1a93feSDuncan Campbell       filename[numChars + 1] = 0;
170df1a93feSDuncan Campbell     }
171df1a93feSDuncan Campbell     // Perform the actual temporary file creation
172c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(filename, "XXXXXX", sizeof(filename)));
173df1a93feSDuncan Campbell     PetscCheck(mkstemp(filename) != -1, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
174df1a93feSDuncan Campbell #else
175df1a93feSDuncan Campbell     // Same thing, but for older C versions which don't have the safe form
176df1a93feSDuncan Campbell     PetscCheck(tmpnam(filename) != NULL, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
177df1a93feSDuncan Campbell #endif
178df1a93feSDuncan Campbell   }
179df1a93feSDuncan Campbell 
180cd2bb8e3SDuncan Campbell   // Broadcast the filename to all other MPI ranks
181cd2bb8e3SDuncan Campbell   PetscCallMPI(MPI_Bcast(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, 0, comm));
182cd2bb8e3SDuncan Campbell 
183e66d8692SHong Zhang   PetscCall(PetscViewerASCIIOpen(comm, filename, &csvViewer));
184df1a93feSDuncan Campbell   PetscCall(PetscViewerPushFormat(csvViewer, PETSC_VIEWER_ASCII_CSV));
185df1a93feSDuncan Campbell 
186df1a93feSDuncan Campbell   // Use the CSV viewer to write out the local network
187df1a93feSDuncan Campbell   PetscCall(DMView_Network_CSV(dm, csvViewer));
188df1a93feSDuncan Campbell 
189df1a93feSDuncan Campbell   // Close the viewer
190df1a93feSDuncan Campbell   PetscCall(PetscViewerDestroy(&csvViewer));
191df1a93feSDuncan Campbell 
1925f25b224SDuncan Campbell   // Generate options string
1935f25b224SDuncan Campbell   PetscCall(PetscMemzero(options, sizeof(options)));
1945f25b224SDuncan Campbell   // If the draw is null run as a "test execute" ie. do nothing just test that the script was called correctly
1955f25b224SDuncan Campbell   PetscCall(PetscStrlcat(options, isnull ? " -tx " : " ", sizeof(options)));
1965f25b224SDuncan Campbell   PetscCall(PetscDrawGetPause(draw, &drawPause));
1975f25b224SDuncan Campbell   if (drawPause > 0) {
1985f25b224SDuncan Campbell     char pausebuffer[64];
1995f25b224SDuncan Campbell     PetscCall(PetscSNPrintf(pausebuffer, sizeof(pausebuffer), "%f", (double)drawPause));
2005f25b224SDuncan Campbell     PetscCall(PetscStrlcat(options, " -dt ", sizeof(options)));
2015f25b224SDuncan Campbell     PetscCall(PetscStrlcat(options, pausebuffer, sizeof(options)));
2025f25b224SDuncan Campbell   }
2035f25b224SDuncan Campbell   if (optionShowRanks || optionRankIsSet) {
2045f25b224SDuncan Campbell     // Show all ranks only if the option is set in code or by the user AND not showing specific ranks AND there is more than one process
2055f25b224SDuncan Campbell     if (optionShowRanks && !optionRankIsSet && size != 1) PetscCall(PetscStrlcat(options, " -dar ", sizeof(options)));
2065f25b224SDuncan Campbell     // Do not show the global plot if the user requests it OR if one specific rank is requested
2075f25b224SDuncan Campbell     if (network->vieweroptions.dontshowglobal || optionRankIsSet) PetscCall(PetscStrlcat(options, " -ncp ", sizeof(options)));
2085f25b224SDuncan Campbell 
2095f25b224SDuncan Campbell     if (optionRankIsSet) {
2105f25b224SDuncan Campbell       // If a range of ranks to draw is specified append it
2115f25b224SDuncan Campbell       PetscCall(PetscStrlcat(options, " -drr ", sizeof(options)));
2125f25b224SDuncan Campbell       PetscCall(PetscStrlcat(options, buffer, sizeof(options)));
2135f25b224SDuncan Campbell     } else {
2145f25b224SDuncan Campbell       // Otherwise, use the options provided in code
2155f25b224SDuncan Campbell       if (network->vieweroptions.viewranks) {
2165f25b224SDuncan Campbell         const PetscInt *viewranks;
2175f25b224SDuncan Campbell         PetscInt        viewrankssize;
2185f25b224SDuncan Campbell         char            rankbuffer[64];
2195f25b224SDuncan Campbell         PetscCall(ISGetTotalIndices(network->vieweroptions.viewranks, &viewranks));
2205f25b224SDuncan Campbell         PetscCall(ISGetSize(network->vieweroptions.viewranks, &viewrankssize));
2215f25b224SDuncan Campbell         PetscCall(PetscStrlcat(options, " -drr ", sizeof(options)));
2225f25b224SDuncan Campbell         for (i = 0; i < viewrankssize; i++) {
2235f25b224SDuncan Campbell           PetscCall(PetscSNPrintf(rankbuffer, sizeof(rankbuffer), "%" PetscInt_FMT, viewranks[i]));
2245f25b224SDuncan Campbell           PetscCall(PetscStrlcat(options, rankbuffer, sizeof(options)));
2255f25b224SDuncan Campbell         }
2265f25b224SDuncan Campbell         PetscCall(ISRestoreTotalIndices(network->vieweroptions.viewranks, &viewranks));
2275f25b224SDuncan Campbell       } // if not provided an IS of viewing ranks, skip viewing
2285f25b224SDuncan Campbell     }
2295f25b224SDuncan Campbell   }
2305039956bSHong Zhang   if (optionShowVertices) {
2315039956bSHong Zhang     // Pass vertices to focus on if defined
2325039956bSHong Zhang     PetscCall(PetscStrlcat(options, " -vsv ", sizeof(options)));
2335039956bSHong Zhang     PetscCall(PetscStrlcat(options, buffer2, sizeof(options)));
2345039956bSHong Zhang     optionViewPadding = PETSC_TRUE;
2355039956bSHong Zhang     // Pass padding if set
2365039956bSHong Zhang     if (optionViewPadding) {
2375039956bSHong Zhang       PetscCall(PetscSNPrintf(buffer2, sizeof(buffer2), "%f", (double)viewPadding));
2385039956bSHong Zhang       PetscCall(PetscStrlcat(options, " -vp ", sizeof(options)));
2395039956bSHong Zhang       PetscCall(PetscStrlcat(options, buffer2, sizeof(options)));
2405039956bSHong Zhang     }
2415039956bSHong Zhang   }
2425f25b224SDuncan Campbell 
2435f25b224SDuncan Campbell   // Check for options for visibility...
2445f25b224SDuncan Campbell   if (showNoNodes) PetscCall(PetscStrlcat(options, " -nn ", sizeof(options)));
2455f25b224SDuncan Campbell   if (showNoNumbering) PetscCall(PetscStrlcat(options, " -nnl -nel ", sizeof(options)));
2465f25b224SDuncan Campbell 
247df1a93feSDuncan Campbell   // Get the value of $PETSC_DIR
248e66d8692SHong Zhang   PetscCall(PetscStrreplace(comm, "${PETSC_DIR}/share/petsc/bin/dmnetwork_view.py", scriptFile, sizeof(scriptFile)));
249df1a93feSDuncan Campbell   PetscCall(PetscFixFilename(scriptFile, scriptFile));
2505f25b224SDuncan Campbell   // Generate the system call for 'python3 $PETSC_DIR/share/petsc/dmnetwork_view.py <options> <file>'
251df1a93feSDuncan Campbell   PetscCall(PetscArrayzero(proccall, sizeof(proccall)));
2525f25b224SDuncan Campbell   PetscCall(PetscSNPrintf(proccall, sizeof(proccall), "%s %s %s %s", PETSC_PYTHON_EXE, scriptFile, options, filename));
253df1a93feSDuncan Campbell 
254df1a93feSDuncan Campbell #if defined(PETSC_HAVE_POPEN)
255df1a93feSDuncan Campbell   // Perform the call to run the python script (Note: while this is called on all ranks POpen will only run on rank 0)
256e66d8692SHong Zhang   PetscCall(PetscPOpen(comm, NULL, proccall, "r", &processFile));
257df1a93feSDuncan Campbell   if (processFile != NULL) {
2585f25b224SDuncan Campbell     while (fgets(buffer, sizeof(buffer), processFile) != NULL) PetscCall(PetscPrintf(comm, "%s", buffer));
259df1a93feSDuncan Campbell   }
260e66d8692SHong Zhang   PetscCall(PetscPClose(comm, processFile));
261df1a93feSDuncan Campbell #else
262df1a93feSDuncan Campbell   // Same thing, but using the standard library for systems that don't have POpen/PClose (only run on rank 0)
263e66d8692SHong Zhang   if (rank == 0) PetscCheck(system(proccall) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to call viewer script");
264df1a93feSDuncan Campbell   // Barrier so that all ranks wait until the call completes
265e66d8692SHong Zhang   PetscCallMPI(MPI_Barrier(comm));
266df1a93feSDuncan Campbell #endif
267df1a93feSDuncan Campbell   // Clean up the temporary file we used using rank 0
268e66d8692SHong Zhang   if (rank == 0) PetscCheck(remove(filename) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to delete temporary file");
269df1a93feSDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
270df1a93feSDuncan Campbell }
271df1a93feSDuncan Campbell 
272d2fd5932SHong Zhang PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
273d2fd5932SHong Zhang {
274df1a93feSDuncan Campbell   PetscBool         iascii, isdraw;
275df1a93feSDuncan Campbell   PetscViewerFormat format;
276d2fd5932SHong Zhang 
277d2fd5932SHong Zhang   PetscFunctionBegin;
278d2fd5932SHong Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
279d2fd5932SHong Zhang   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
280df1a93feSDuncan Campbell   PetscCall(PetscViewerGetFormat(viewer, &format));
281df1a93feSDuncan Campbell 
282df1a93feSDuncan Campbell   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
283df1a93feSDuncan Campbell   if (isdraw) {
284df1a93feSDuncan Campbell     PetscCall(DMView_Network_Matplotlib(dm, viewer));
285df1a93feSDuncan Campbell     PetscFunctionReturn(PETSC_SUCCESS);
286df1a93feSDuncan Campbell   }
287df1a93feSDuncan Campbell 
288d2fd5932SHong Zhang   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
289d2fd5932SHong Zhang   if (iascii) {
290d2fd5932SHong Zhang     const PetscInt *cone, *vtx, *edges;
291d2fd5932SHong Zhang     PetscInt        vfrom, vto, i, j, nv, ne, nsv, p, nsubnet;
292d2fd5932SHong Zhang     DM_Network     *network = (DM_Network *)dm->data;
293df1a93feSDuncan Campbell     PetscMPIInt     rank;
294df1a93feSDuncan Campbell 
295df1a93feSDuncan Campbell     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
296df1a93feSDuncan Campbell     if (format == PETSC_VIEWER_ASCII_CSV) {
297df1a93feSDuncan Campbell       PetscCall(DMView_Network_CSV(dm, viewer));
298df1a93feSDuncan Campbell       PetscFunctionReturn(PETSC_SUCCESS);
299df1a93feSDuncan Campbell     }
300d2fd5932SHong Zhang 
301d2fd5932SHong Zhang     nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */
302df1a93feSDuncan Campbell     if (!rank) {
303d2fd5932SHong 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,
304d2fd5932SHong Zhang                             network->cloneshared->Nsvtx));
305d2fd5932SHong Zhang     }
306d2fd5932SHong Zhang 
307d2fd5932SHong Zhang     PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL));
308d2fd5932SHong Zhang     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
309d2fd5932SHong Zhang     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv));
310d2fd5932SHong Zhang 
311d2fd5932SHong Zhang     for (i = 0; i < nsubnet; i++) {
312d2fd5932SHong Zhang       PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges));
313d2fd5932SHong Zhang       if (ne) {
314d2fd5932SHong Zhang         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv));
315d2fd5932SHong Zhang         for (j = 0; j < ne; j++) {
316d2fd5932SHong Zhang           p = edges[j];
317d2fd5932SHong Zhang           PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone));
318d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom));
319d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto));
320d2fd5932SHong Zhang           PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p));
321d2fd5932SHong Zhang           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto));
322d2fd5932SHong Zhang         }
323d2fd5932SHong Zhang       }
324d2fd5932SHong Zhang     }
325d2fd5932SHong Zhang 
326d2fd5932SHong Zhang     /* Shared vertices */
327d2fd5932SHong Zhang     PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx));
328d2fd5932SHong Zhang     if (nsv) {
329d2fd5932SHong Zhang       PetscInt        gidx;
330d2fd5932SHong Zhang       PetscBool       ghost;
331d2fd5932SHong Zhang       const PetscInt *sv = NULL;
332d2fd5932SHong Zhang 
333d2fd5932SHong Zhang       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n"));
334d2fd5932SHong Zhang       for (i = 0; i < nsv; i++) {
335d2fd5932SHong Zhang         PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost));
336d2fd5932SHong Zhang         if (ghost) continue;
337d2fd5932SHong Zhang 
338d2fd5932SHong Zhang         PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv));
339d2fd5932SHong Zhang         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1]));
340d2fd5932SHong Zhang         for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1]));
341d2fd5932SHong Zhang       }
342d2fd5932SHong Zhang     }
343d2fd5932SHong Zhang     PetscCall(PetscViewerFlush(viewer));
344d2fd5932SHong Zhang     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
345d2fd5932SHong Zhang   } else PetscCheck(iascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
346d2fd5932SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
347d2fd5932SHong Zhang }
3485f25b224SDuncan Campbell 
3495f25b224SDuncan Campbell /*@
3505f25b224SDuncan Campbell   DMNetworkViewSetShowRanks - Sets viewing the `DMETNWORK` on each rank individually.
3515f25b224SDuncan Campbell 
3525f25b224SDuncan Campbell   Logically Collective
3535f25b224SDuncan Campbell 
3545f25b224SDuncan Campbell   Input Parameter:
3555f25b224SDuncan Campbell . dm - the `DMNETWORK` object
3565f25b224SDuncan Campbell 
3575f25b224SDuncan Campbell   Output Parameter:
3585f25b224SDuncan Campbell . showranks - `PETSC_TRUE` if viewing each rank's sub network individually
3595f25b224SDuncan Campbell 
3605f25b224SDuncan Campbell   Level: beginner
3615f25b224SDuncan Campbell 
3625f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
3635f25b224SDuncan Campbell @*/
3645f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowRanks(DM dm, PetscBool showranks)
3655f25b224SDuncan Campbell {
3665f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
3675f25b224SDuncan Campbell 
3685f25b224SDuncan Campbell   PetscFunctionBegin;
3695f25b224SDuncan Campbell   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
3705f25b224SDuncan Campbell   network->vieweroptions.showallranks = showranks;
3715f25b224SDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
3725f25b224SDuncan Campbell }
3735f25b224SDuncan Campbell 
3745f25b224SDuncan Campbell /*@
3755f25b224SDuncan Campbell   DMNetworkViewSetShowGlobal - Set viewing the global network.
3765f25b224SDuncan Campbell 
3775f25b224SDuncan Campbell   Logically Collective
3785f25b224SDuncan Campbell 
3795f25b224SDuncan Campbell   Input Parameter:
3805f25b224SDuncan Campbell . dm - the `DMNETWORK` object
3815f25b224SDuncan Campbell 
3825f25b224SDuncan Campbell   Output Parameter:
3835f25b224SDuncan Campbell . showglobal - `PETSC_TRUE` if viewing the global network
3845f25b224SDuncan Campbell 
3855f25b224SDuncan Campbell   Level: beginner
3865f25b224SDuncan Campbell 
3875f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
3885f25b224SDuncan Campbell @*/
3895f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowGlobal(DM dm, PetscBool showglobal)
3905f25b224SDuncan Campbell {
3915f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
3925f25b224SDuncan Campbell 
3935f25b224SDuncan Campbell   PetscFunctionBegin;
3945f25b224SDuncan Campbell   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
3955f25b224SDuncan Campbell   network->vieweroptions.dontshowglobal = (PetscBool)(!showglobal);
3965f25b224SDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
3975f25b224SDuncan Campbell }
3985f25b224SDuncan Campbell 
3995f25b224SDuncan Campbell /*@
4005f25b224SDuncan Campbell   DMNetworkViewSetShowVertices - Sets whether to display the vertices in viewing routines.
4015f25b224SDuncan Campbell 
4025f25b224SDuncan Campbell   Logically Collective
4035f25b224SDuncan Campbell 
4045f25b224SDuncan Campbell   Input Parameter:
4055f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4065f25b224SDuncan Campbell 
4075f25b224SDuncan Campbell   Output Parameter:
4085f25b224SDuncan Campbell . showvertices - `PETSC_TRUE` if visualizing the vertices
4095f25b224SDuncan Campbell 
4105f25b224SDuncan Campbell   Level: beginner
4115f25b224SDuncan Campbell 
4125f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
4135f25b224SDuncan Campbell @*/
4145f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowVertices(DM dm, PetscBool showvertices)
4155f25b224SDuncan Campbell {
4165f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
4175f25b224SDuncan Campbell 
4185f25b224SDuncan Campbell   PetscFunctionBegin;
4195f25b224SDuncan Campbell   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4205f25b224SDuncan Campbell   network->vieweroptions.shownovertices = (PetscBool)(!showvertices);
4215f25b224SDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
4225f25b224SDuncan Campbell }
4235f25b224SDuncan Campbell 
4245f25b224SDuncan Campbell /*@
4255f25b224SDuncan Campbell   DMNetworkViewSetShowNumbering - Set displaying the numbering of edges and vertices in viewing routines.
4265f25b224SDuncan Campbell 
4275f25b224SDuncan Campbell   Logically Collective
4285f25b224SDuncan Campbell 
4295f25b224SDuncan Campbell   Input Parameter:
4305f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4315f25b224SDuncan Campbell 
4325f25b224SDuncan Campbell   Output Parameter:
4335f25b224SDuncan Campbell . shownumbering - `PETSC_TRUE` if displaying the numbering of edges and vertices
4345f25b224SDuncan Campbell 
4355f25b224SDuncan Campbell   Level: beginner
4365f25b224SDuncan Campbell 
4375f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetViewRanks()`
4385f25b224SDuncan Campbell @*/
4395f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowNumbering(DM dm, PetscBool shownumbering)
4405f25b224SDuncan Campbell {
4415f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
4425f25b224SDuncan Campbell 
4435f25b224SDuncan Campbell   PetscFunctionBegin;
4445f25b224SDuncan Campbell   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4455f25b224SDuncan Campbell   network->vieweroptions.shownonumbering = (PetscBool)(!shownumbering);
4465f25b224SDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
4475f25b224SDuncan Campbell }
4485f25b224SDuncan Campbell 
4495f25b224SDuncan Campbell /*@
4505f25b224SDuncan Campbell   DMNetworkViewSetViewRanks - View the `DMNETWORK` on each of the specified ranks individually.
4515f25b224SDuncan Campbell 
4525f25b224SDuncan Campbell   Collective
4535f25b224SDuncan Campbell 
4545f25b224SDuncan Campbell   Input Parameter:
4555f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4565f25b224SDuncan Campbell 
4575f25b224SDuncan Campbell   Output Parameter:
4585f25b224SDuncan Campbell . viewranks - set of ranks to view the `DMNETWORK` on individually
4595f25b224SDuncan Campbell 
4605f25b224SDuncan Campbell   Level: beginner
4615f25b224SDuncan Campbell 
4625f25b224SDuncan Campbell   Note:
4635f25b224SDuncan Campbell   `DMNetwork` takes ownership of the input viewranks `IS`, it should be destroyed by the caller.
4645f25b224SDuncan Campbell 
4655f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`
4665f25b224SDuncan Campbell @*/
4675f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetViewRanks(DM dm, IS viewranks)
4685f25b224SDuncan Campbell {
4695f25b224SDuncan Campbell   DM_Network *network = (DM_Network *)dm->data;
4705f25b224SDuncan Campbell 
4715f25b224SDuncan Campbell   PetscFunctionBegin;
4725f25b224SDuncan Campbell   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4735f25b224SDuncan Campbell   PetscValidHeaderSpecific(viewranks, IS_CLASSID, 2);
4745f25b224SDuncan Campbell   PetscCheckSameComm(dm, 1, viewranks, 2);
4755f25b224SDuncan Campbell   PetscCall(ISDestroy(&network->vieweroptions.viewranks));
4765f25b224SDuncan Campbell   PetscCall(PetscObjectReference((PetscObject)viewranks));
4775f25b224SDuncan Campbell   network->vieweroptions.viewranks = viewranks;
4785f25b224SDuncan Campbell   PetscFunctionReturn(PETSC_SUCCESS);
4795f25b224SDuncan Campbell }
480