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; 18*7cd471e7SHong Zhang PetscInt nsubnets, i, subnet, nvertices, nedges, vertex, edge, gidx, ncomp; 19df1a93feSDuncan Campbell PetscInt vertexOffsets[2], globalEdgeVertices[2]; 20*7cd471e7SHong 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; 27*7cd471e7SHong Zhang // Get the coordinate information from dmcoords 28*7cd471e7SHong Zhang PetscCheck(dm->coordinates[0].dm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "CoordinateDM not created"); 29df1a93feSDuncan Campbell PetscCall(DMGetCoordinateDM(dm, &dmcoords)); 30*7cd471e7SHong Zhang 31*7cd471e7SHong Zhang PetscCall(DMGetCoordinateDim(dmcoords, &i)); 32*7cd471e7SHong Zhang PetscCheck(i == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim %" PetscInt_FMT " != 2 is not supporte yet", i); 33*7cd471e7SHong Zhang 34*7cd471e7SHong Zhang // Get the coordinate vector from dm 35df1a93feSDuncan Campbell PetscCall(DMGetCoordinatesLocal(dm, &allVertexCoords)); 36*7cd471e7SHong Zhang 37df1a93feSDuncan Campbell // Get the MPI communicator and this process' rank 38*7cd471e7SHong Zhang PetscCall(PetscObjectGetComm((PetscObject)dmcoords, &comm)); 39df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(comm, &rank)); 40*7cd471e7SHong 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) 48*7cd471e7SHong Zhang PetscCall(DMNetworkGetNumSubNetworks(dmcoords, NULL, &nsubnets)); 49df1a93feSDuncan Campbell for (subnet = 0; subnet < nsubnets; subnet++) { 50df1a93feSDuncan Campbell // Get the subnetwork's vertices and edges 51*7cd471e7SHong 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]; 56*7cd471e7SHong 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 61*7cd471e7SHong 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 65*7cd471e7SHong Zhang // Get vertex color; TODO: name 66*7cd471e7SHong Zhang PetscCall(DMNetworkGetNumComponents(dmcoords, vertex, &ncomp)); 67*7cd471e7SHong Zhang PetscCheck(ncomp <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "num of components %" PetscInt_FMT " must be <= 1", ncomp); 68*7cd471e7SHong Zhang color = 0.0; 69*7cd471e7SHong Zhang if (ncomp == 1) { 70*7cd471e7SHong Zhang PetscCall(DMNetworkGetComponent(dmcoords, vertex, 0, NULL, (void **)&color_ptr, NULL)); 71*7cd471e7SHong Zhang color = *color_ptr; 72*7cd471e7SHong Zhang } 73*7cd471e7SHong 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]; 79*7cd471e7SHong Zhang PetscCall(DMNetworkGetConnectedVertices(dmcoords, edge, &edgeVertices)); 80*7cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[0], &globalEdgeVertices[0])); 81*7cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[1], &globalEdgeVertices[1])); 82*7cd471e7SHong 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; 985f25b224SDuncan Campbell char filename[PETSC_MAX_PATH_LEN + 1], options[512], proccall[PETSC_MAX_PATH_LEN + 512], scriptFile[PETSC_MAX_PATH_LEN + 1], buffer[256]; 99df1a93feSDuncan Campbell PetscViewer csvViewer; 100df1a93feSDuncan Campbell FILE *processFile = NULL; 1015f25b224SDuncan Campbell PetscBool isnull, optionShowRanks = PETSC_FALSE, optionRankIsSet = PETSC_FALSE, showNoNodes = PETSC_FALSE, showNoNumbering = PETSC_FALSE; 102df1a93feSDuncan Campbell PetscDraw draw; 1035f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 1045f25b224SDuncan Campbell PetscReal drawPause; 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)); 1435f25b224SDuncan Campbell PetscOptionsEnd(); 1445f25b224SDuncan Campbell 145df1a93feSDuncan Campbell // Generate and broadcast the temporary file name from rank 0 146df1a93feSDuncan Campbell if (rank == 0) { 147df1a93feSDuncan Campbell #if defined(PETSC_HAVE_TMPNAM_S) 148df1a93feSDuncan Campbell // Acquire a temporary file to write to and open an ASCII/CSV viewer 149df1a93feSDuncan Campbell PetscCheck(tmpnam_s(filename, sizeof(filename)) == 0, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 15012e0ef65SDuncan Campbell #elif defined(PETSC_HAVE_MKSTEMP) 151cd2bb8e3SDuncan Campbell PetscBool isTmpOverridden; 152df1a93feSDuncan Campbell size_t numChars; 153df1a93feSDuncan Campbell // Same thing, but for POSIX systems on which tmpnam is deprecated 154df1a93feSDuncan 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 155df1a93feSDuncan 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 1566e4289a0SDuncan Campbell PetscCall(PetscOptionsGetString(NULL, NULL, "-dmnetwork_view_tmpdir", filename, sizeof(filename), &isTmpOverridden)); 1576e4289a0SDuncan Campbell // If not specified by option try using a shared tmp on the system 1586e4289a0SDuncan Campbell if (!isTmpOverridden) { 1596e4289a0SDuncan Campbell // Validate that if tmp is not overridden it is at least shared 160cd2bb8e3SDuncan 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"); 161cd2bb8e3SDuncan Campbell PetscCall(PetscGetTmp(PETSC_COMM_SELF, filename, sizeof(filename))); 162cd2bb8e3SDuncan Campbell } 163df1a93feSDuncan Campbell // Make sure the filename ends with a '/' 164df1a93feSDuncan Campbell PetscCall(PetscStrlen(filename, &numChars)); 165df1a93feSDuncan Campbell if (filename[numChars - 1] != '/') { 166df1a93feSDuncan Campbell filename[numChars] = '/'; 167df1a93feSDuncan Campbell filename[numChars + 1] = 0; 168df1a93feSDuncan Campbell } 169df1a93feSDuncan Campbell // Perform the actual temporary file creation 170c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(filename, "XXXXXX", sizeof(filename))); 171df1a93feSDuncan Campbell PetscCheck(mkstemp(filename) != -1, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 172df1a93feSDuncan Campbell #else 173df1a93feSDuncan Campbell // Same thing, but for older C versions which don't have the safe form 174df1a93feSDuncan Campbell PetscCheck(tmpnam(filename) != NULL, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 175df1a93feSDuncan Campbell #endif 176df1a93feSDuncan Campbell } 177df1a93feSDuncan Campbell 178cd2bb8e3SDuncan Campbell // Broadcast the filename to all other MPI ranks 179cd2bb8e3SDuncan Campbell PetscCallMPI(MPI_Bcast(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, 0, comm)); 180cd2bb8e3SDuncan Campbell 181e66d8692SHong Zhang PetscCall(PetscViewerASCIIOpen(comm, filename, &csvViewer)); 182df1a93feSDuncan Campbell PetscCall(PetscViewerPushFormat(csvViewer, PETSC_VIEWER_ASCII_CSV)); 183df1a93feSDuncan Campbell 184df1a93feSDuncan Campbell // Use the CSV viewer to write out the local network 185df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, csvViewer)); 186df1a93feSDuncan Campbell 187df1a93feSDuncan Campbell // Close the viewer 188df1a93feSDuncan Campbell PetscCall(PetscViewerDestroy(&csvViewer)); 189df1a93feSDuncan Campbell 1905f25b224SDuncan Campbell // Generate options string 1915f25b224SDuncan Campbell PetscCall(PetscMemzero(options, sizeof(options))); 1925f25b224SDuncan Campbell // If the draw is null run as a "test execute" ie. do nothing just test that the script was called correctly 1935f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, isnull ? " -tx " : " ", sizeof(options))); 1945f25b224SDuncan Campbell PetscCall(PetscDrawGetPause(draw, &drawPause)); 1955f25b224SDuncan Campbell if (drawPause > 0) { 1965f25b224SDuncan Campbell char pausebuffer[64]; 1975f25b224SDuncan Campbell PetscCall(PetscSNPrintf(pausebuffer, sizeof(pausebuffer), "%f", (double)drawPause)); 1985f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -dt ", sizeof(options))); 1995f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, pausebuffer, sizeof(options))); 2005f25b224SDuncan Campbell } 2015f25b224SDuncan Campbell if (optionShowRanks || optionRankIsSet) { 2025f25b224SDuncan 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 2035f25b224SDuncan Campbell if (optionShowRanks && !optionRankIsSet && size != 1) PetscCall(PetscStrlcat(options, " -dar ", sizeof(options))); 2045f25b224SDuncan Campbell // Do not show the global plot if the user requests it OR if one specific rank is requested 2055f25b224SDuncan Campbell if (network->vieweroptions.dontshowglobal || optionRankIsSet) PetscCall(PetscStrlcat(options, " -ncp ", sizeof(options))); 2065f25b224SDuncan Campbell 2075f25b224SDuncan Campbell if (optionRankIsSet) { 2085f25b224SDuncan Campbell // If a range of ranks to draw is specified append it 2095f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options))); 2105f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, buffer, sizeof(options))); 2115f25b224SDuncan Campbell } else { 2125f25b224SDuncan Campbell // Otherwise, use the options provided in code 2135f25b224SDuncan Campbell if (network->vieweroptions.viewranks) { 2145f25b224SDuncan Campbell const PetscInt *viewranks; 2155f25b224SDuncan Campbell PetscInt viewrankssize; 2165f25b224SDuncan Campbell char rankbuffer[64]; 2175f25b224SDuncan Campbell PetscCall(ISGetTotalIndices(network->vieweroptions.viewranks, &viewranks)); 2185f25b224SDuncan Campbell PetscCall(ISGetSize(network->vieweroptions.viewranks, &viewrankssize)); 2195f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options))); 2205f25b224SDuncan Campbell for (i = 0; i < viewrankssize; i++) { 2215f25b224SDuncan Campbell PetscCall(PetscSNPrintf(rankbuffer, sizeof(rankbuffer), "%" PetscInt_FMT, viewranks[i])); 2225f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, rankbuffer, sizeof(options))); 2235f25b224SDuncan Campbell } 2245f25b224SDuncan Campbell PetscCall(ISRestoreTotalIndices(network->vieweroptions.viewranks, &viewranks)); 2255f25b224SDuncan Campbell } // if not provided an IS of viewing ranks, skip viewing 2265f25b224SDuncan Campbell } 2275f25b224SDuncan Campbell } 2285f25b224SDuncan Campbell 2295f25b224SDuncan Campbell // Check for options for visibility... 2305f25b224SDuncan Campbell if (showNoNodes) PetscCall(PetscStrlcat(options, " -nn ", sizeof(options))); 2315f25b224SDuncan Campbell if (showNoNumbering) PetscCall(PetscStrlcat(options, " -nnl -nel ", sizeof(options))); 2325f25b224SDuncan Campbell 233df1a93feSDuncan Campbell // Get the value of $PETSC_DIR 234e66d8692SHong Zhang PetscCall(PetscStrreplace(comm, "${PETSC_DIR}/share/petsc/bin/dmnetwork_view.py", scriptFile, sizeof(scriptFile))); 235df1a93feSDuncan Campbell PetscCall(PetscFixFilename(scriptFile, scriptFile)); 2365f25b224SDuncan Campbell // Generate the system call for 'python3 $PETSC_DIR/share/petsc/dmnetwork_view.py <options> <file>' 237df1a93feSDuncan Campbell PetscCall(PetscArrayzero(proccall, sizeof(proccall))); 2385f25b224SDuncan Campbell PetscCall(PetscSNPrintf(proccall, sizeof(proccall), "%s %s %s %s", PETSC_PYTHON_EXE, scriptFile, options, filename)); 239df1a93feSDuncan Campbell 240df1a93feSDuncan Campbell #if defined(PETSC_HAVE_POPEN) 241df1a93feSDuncan Campbell // Perform the call to run the python script (Note: while this is called on all ranks POpen will only run on rank 0) 242e66d8692SHong Zhang PetscCall(PetscPOpen(comm, NULL, proccall, "r", &processFile)); 243df1a93feSDuncan Campbell if (processFile != NULL) { 2445f25b224SDuncan Campbell while (fgets(buffer, sizeof(buffer), processFile) != NULL) PetscCall(PetscPrintf(comm, "%s", buffer)); 245df1a93feSDuncan Campbell } 246e66d8692SHong Zhang PetscCall(PetscPClose(comm, processFile)); 247df1a93feSDuncan Campbell #else 248df1a93feSDuncan Campbell // Same thing, but using the standard library for systems that don't have POpen/PClose (only run on rank 0) 249e66d8692SHong Zhang if (rank == 0) PetscCheck(system(proccall) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to call viewer script"); 250df1a93feSDuncan Campbell // Barrier so that all ranks wait until the call completes 251e66d8692SHong Zhang PetscCallMPI(MPI_Barrier(comm)); 252df1a93feSDuncan Campbell #endif 253df1a93feSDuncan Campbell // Clean up the temporary file we used using rank 0 254e66d8692SHong Zhang if (rank == 0) PetscCheck(remove(filename) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to delete temporary file"); 255df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 256df1a93feSDuncan Campbell } 257df1a93feSDuncan Campbell 258d2fd5932SHong Zhang PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 259d2fd5932SHong Zhang { 260df1a93feSDuncan Campbell PetscBool iascii, isdraw; 261df1a93feSDuncan Campbell PetscViewerFormat format; 262d2fd5932SHong Zhang 263d2fd5932SHong Zhang PetscFunctionBegin; 264d2fd5932SHong Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 265d2fd5932SHong Zhang PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 266df1a93feSDuncan Campbell PetscCall(PetscViewerGetFormat(viewer, &format)); 267df1a93feSDuncan Campbell 268df1a93feSDuncan Campbell PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 269df1a93feSDuncan Campbell if (isdraw) { 270df1a93feSDuncan Campbell PetscCall(DMView_Network_Matplotlib(dm, viewer)); 271df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 272df1a93feSDuncan Campbell } 273df1a93feSDuncan Campbell 274d2fd5932SHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 275d2fd5932SHong Zhang if (iascii) { 276d2fd5932SHong Zhang const PetscInt *cone, *vtx, *edges; 277d2fd5932SHong Zhang PetscInt vfrom, vto, i, j, nv, ne, nsv, p, nsubnet; 278d2fd5932SHong Zhang DM_Network *network = (DM_Network *)dm->data; 279df1a93feSDuncan Campbell PetscMPIInt rank; 280df1a93feSDuncan Campbell 281df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 282df1a93feSDuncan Campbell if (format == PETSC_VIEWER_ASCII_CSV) { 283df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, viewer)); 284df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 285df1a93feSDuncan Campbell } 286d2fd5932SHong Zhang 287d2fd5932SHong Zhang nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */ 288df1a93feSDuncan Campbell if (!rank) { 289d2fd5932SHong 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, 290d2fd5932SHong Zhang network->cloneshared->Nsvtx)); 291d2fd5932SHong Zhang } 292d2fd5932SHong Zhang 293d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL)); 294d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 295d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv)); 296d2fd5932SHong Zhang 297d2fd5932SHong Zhang for (i = 0; i < nsubnet; i++) { 298d2fd5932SHong Zhang PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges)); 299d2fd5932SHong Zhang if (ne) { 300d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv)); 301d2fd5932SHong Zhang for (j = 0; j < ne; j++) { 302d2fd5932SHong Zhang p = edges[j]; 303d2fd5932SHong Zhang PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone)); 304d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom)); 305d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto)); 306d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p)); 307d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto)); 308d2fd5932SHong Zhang } 309d2fd5932SHong Zhang } 310d2fd5932SHong Zhang } 311d2fd5932SHong Zhang 312d2fd5932SHong Zhang /* Shared vertices */ 313d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx)); 314d2fd5932SHong Zhang if (nsv) { 315d2fd5932SHong Zhang PetscInt gidx; 316d2fd5932SHong Zhang PetscBool ghost; 317d2fd5932SHong Zhang const PetscInt *sv = NULL; 318d2fd5932SHong Zhang 319d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n")); 320d2fd5932SHong Zhang for (i = 0; i < nsv; i++) { 321d2fd5932SHong Zhang PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost)); 322d2fd5932SHong Zhang if (ghost) continue; 323d2fd5932SHong Zhang 324d2fd5932SHong Zhang PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv)); 325d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1])); 326d2fd5932SHong Zhang for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1])); 327d2fd5932SHong Zhang } 328d2fd5932SHong Zhang } 329d2fd5932SHong Zhang PetscCall(PetscViewerFlush(viewer)); 330d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 331d2fd5932SHong Zhang } else PetscCheck(iascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name); 332d2fd5932SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 333d2fd5932SHong Zhang } 3345f25b224SDuncan Campbell 3355f25b224SDuncan Campbell /*@ 3365f25b224SDuncan Campbell DMNetworkViewSetShowRanks - Sets viewing the `DMETNWORK` on each rank individually. 3375f25b224SDuncan Campbell 3385f25b224SDuncan Campbell Logically Collective 3395f25b224SDuncan Campbell 3405f25b224SDuncan Campbell Input Parameter: 3415f25b224SDuncan Campbell . dm - the `DMNETWORK` object 3425f25b224SDuncan Campbell 3435f25b224SDuncan Campbell Output Parameter: 3445f25b224SDuncan Campbell . showranks - `PETSC_TRUE` if viewing each rank's sub network individually 3455f25b224SDuncan Campbell 3465f25b224SDuncan Campbell Level: beginner 3475f25b224SDuncan Campbell 3485f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 3495f25b224SDuncan Campbell @*/ 3505f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowRanks(DM dm, PetscBool showranks) 3515f25b224SDuncan Campbell { 3525f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 3535f25b224SDuncan Campbell 3545f25b224SDuncan Campbell PetscFunctionBegin; 3555f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 3565f25b224SDuncan Campbell network->vieweroptions.showallranks = showranks; 3575f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 3585f25b224SDuncan Campbell } 3595f25b224SDuncan Campbell 3605f25b224SDuncan Campbell /*@ 3615f25b224SDuncan Campbell DMNetworkViewSetShowGlobal - Set viewing the global network. 3625f25b224SDuncan Campbell 3635f25b224SDuncan Campbell Logically Collective 3645f25b224SDuncan Campbell 3655f25b224SDuncan Campbell Input Parameter: 3665f25b224SDuncan Campbell . dm - the `DMNETWORK` object 3675f25b224SDuncan Campbell 3685f25b224SDuncan Campbell Output Parameter: 3695f25b224SDuncan Campbell . showglobal - `PETSC_TRUE` if viewing the global network 3705f25b224SDuncan Campbell 3715f25b224SDuncan Campbell Level: beginner 3725f25b224SDuncan Campbell 3735f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 3745f25b224SDuncan Campbell @*/ 3755f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowGlobal(DM dm, PetscBool showglobal) 3765f25b224SDuncan Campbell { 3775f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 3785f25b224SDuncan Campbell 3795f25b224SDuncan Campbell PetscFunctionBegin; 3805f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 3815f25b224SDuncan Campbell network->vieweroptions.dontshowglobal = (PetscBool)(!showglobal); 3825f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 3835f25b224SDuncan Campbell } 3845f25b224SDuncan Campbell 3855f25b224SDuncan Campbell /*@ 3865f25b224SDuncan Campbell DMNetworkViewSetShowVertices - Sets whether to display the vertices in viewing routines. 3875f25b224SDuncan Campbell 3885f25b224SDuncan Campbell Logically Collective 3895f25b224SDuncan Campbell 3905f25b224SDuncan Campbell Input Parameter: 3915f25b224SDuncan Campbell . dm - the `DMNETWORK` object 3925f25b224SDuncan Campbell 3935f25b224SDuncan Campbell Output Parameter: 3945f25b224SDuncan Campbell . showvertices - `PETSC_TRUE` if visualizing the vertices 3955f25b224SDuncan Campbell 3965f25b224SDuncan Campbell Level: beginner 3975f25b224SDuncan Campbell 3985f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 3995f25b224SDuncan Campbell @*/ 4005f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowVertices(DM dm, PetscBool showvertices) 4015f25b224SDuncan Campbell { 4025f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 4035f25b224SDuncan Campbell 4045f25b224SDuncan Campbell PetscFunctionBegin; 4055f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 4065f25b224SDuncan Campbell network->vieweroptions.shownovertices = (PetscBool)(!showvertices); 4075f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 4085f25b224SDuncan Campbell } 4095f25b224SDuncan Campbell 4105f25b224SDuncan Campbell /*@ 4115f25b224SDuncan Campbell DMNetworkViewSetShowNumbering - Set displaying the numbering of edges and vertices in viewing routines. 4125f25b224SDuncan Campbell 4135f25b224SDuncan Campbell Logically Collective 4145f25b224SDuncan Campbell 4155f25b224SDuncan Campbell Input Parameter: 4165f25b224SDuncan Campbell . dm - the `DMNETWORK` object 4175f25b224SDuncan Campbell 4185f25b224SDuncan Campbell Output Parameter: 4195f25b224SDuncan Campbell . shownumbering - `PETSC_TRUE` if displaying the numbering of edges and vertices 4205f25b224SDuncan Campbell 4215f25b224SDuncan Campbell Level: beginner 4225f25b224SDuncan Campbell 4235f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetViewRanks()` 4245f25b224SDuncan Campbell @*/ 4255f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowNumbering(DM dm, PetscBool shownumbering) 4265f25b224SDuncan Campbell { 4275f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 4285f25b224SDuncan Campbell 4295f25b224SDuncan Campbell PetscFunctionBegin; 4305f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 4315f25b224SDuncan Campbell network->vieweroptions.shownonumbering = (PetscBool)(!shownumbering); 4325f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 4335f25b224SDuncan Campbell } 4345f25b224SDuncan Campbell 4355f25b224SDuncan Campbell /*@ 4365f25b224SDuncan Campbell DMNetworkViewSetViewRanks - View the `DMNETWORK` on each of the specified ranks individually. 4375f25b224SDuncan Campbell 4385f25b224SDuncan Campbell Collective 4395f25b224SDuncan Campbell 4405f25b224SDuncan Campbell Input Parameter: 4415f25b224SDuncan Campbell . dm - the `DMNETWORK` object 4425f25b224SDuncan Campbell 4435f25b224SDuncan Campbell Output Parameter: 4445f25b224SDuncan Campbell . viewranks - set of ranks to view the `DMNETWORK` on individually 4455f25b224SDuncan Campbell 4465f25b224SDuncan Campbell Level: beginner 4475f25b224SDuncan Campbell 4485f25b224SDuncan Campbell Note: 4495f25b224SDuncan Campbell `DMNetwork` takes ownership of the input viewranks `IS`, it should be destroyed by the caller. 4505f25b224SDuncan Campbell 4515f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()` 4525f25b224SDuncan Campbell @*/ 4535f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetViewRanks(DM dm, IS viewranks) 4545f25b224SDuncan Campbell { 4555f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 4565f25b224SDuncan Campbell 4575f25b224SDuncan Campbell PetscFunctionBegin; 4585f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 4595f25b224SDuncan Campbell PetscValidHeaderSpecific(viewranks, IS_CLASSID, 2); 4605f25b224SDuncan Campbell PetscCheckSameComm(dm, 1, viewranks, 2); 4615f25b224SDuncan Campbell PetscCall(ISDestroy(&network->vieweroptions.viewranks)); 4625f25b224SDuncan Campbell PetscCall(PetscObjectReference((PetscObject)viewranks)); 4635f25b224SDuncan Campbell network->vieweroptions.viewranks = viewranks; 4645f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 4655f25b224SDuncan Campbell } 466