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 6*5f25b224SDuncan Campbell #include "petsc/private/petscimpl.h" 7*5f25b224SDuncan Campbell #include "petscerror.h" 8*5f25b224SDuncan Campbell #include "petscis.h" 9*5f25b224SDuncan Campbell #include "petscstring.h" 10*5f25b224SDuncan Campbell #include "petscsys.h" 11*5f25b224SDuncan 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; 18df1a93feSDuncan Campbell PetscInt nsubnets, i, subnet, nvertices, nedges, vertex, edge; 19df1a93feSDuncan Campbell PetscInt vertexOffsets[2], globalEdgeVertices[2]; 20df1a93feSDuncan Campbell PetscScalar vertexCoords[2]; 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; 27df1a93feSDuncan Campbell // Get the network containing coordinate information 28df1a93feSDuncan Campbell PetscCall(DMGetCoordinateDM(dm, &dmcoords)); 29df1a93feSDuncan Campbell // Get the coordinate vector for the network 30df1a93feSDuncan Campbell PetscCall(DMGetCoordinatesLocal(dm, &allVertexCoords)); 31df1a93feSDuncan Campbell // Get the MPI communicator and this process' rank 32df1a93feSDuncan Campbell PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 33df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(comm, &rank)); 34df1a93feSDuncan Campbell // Start synchronized printing 35df1a93feSDuncan Campbell PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 36df1a93feSDuncan Campbell 37df1a93feSDuncan Campbell // Write the header 38df1a93feSDuncan Campbell PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Type,Rank,ID,X,Y,Z,Name,Color\n")); 39df1a93feSDuncan Campbell 40df1a93feSDuncan Campbell // Iterate each subnetwork (Note: We need to get the global number of subnets apparently) 41df1a93feSDuncan Campbell PetscCall(DMNetworkGetNumSubNetworks(dm, NULL, &nsubnets)); 42df1a93feSDuncan Campbell for (subnet = 0; subnet < nsubnets; subnet++) { 43df1a93feSDuncan Campbell // Get the subnetwork's vertices and edges 44df1a93feSDuncan Campbell PetscCall(DMNetworkGetSubnetwork(dm, subnet, &nvertices, &nedges, &vertices, &edges)); 45df1a93feSDuncan Campbell 46df1a93feSDuncan Campbell // Write out each vertex 47df1a93feSDuncan Campbell for (i = 0; i < nvertices; i++) { 48df1a93feSDuncan Campbell vertex = vertices[i]; 49df1a93feSDuncan Campbell // Get the offset into the coordinate vector for the vertex 50df1a93feSDuncan Campbell PetscCall(DMNetworkGetLocalVecOffset(dmcoords, vertex, ALL_COMPONENTS, vertexOffsets)); 51df1a93feSDuncan Campbell vertexOffsets[1] = vertexOffsets[0] + 1; 52df1a93feSDuncan Campbell // Remap vertex to the global value 53df1a93feSDuncan Campbell PetscCall(DMNetworkGetGlobalVertexIndex(dm, vertex, &vertex)); 54df1a93feSDuncan Campbell // Get the vertex position from the coordinate vector 55df1a93feSDuncan Campbell PetscCall(VecGetValues(allVertexCoords, 2, vertexOffsets, vertexCoords)); 56df1a93feSDuncan Campbell 57df1a93feSDuncan Campbell // TODO: Determine vertex color/name 58df1a93feSDuncan 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)); 59df1a93feSDuncan Campbell } 60df1a93feSDuncan Campbell 61df1a93feSDuncan Campbell // Write out each edge 62df1a93feSDuncan Campbell for (i = 0; i < nedges; i++) { 63df1a93feSDuncan Campbell edge = edges[i]; 64df1a93feSDuncan Campbell PetscCall(DMNetworkGetConnectedVertices(dm, edge, &edgeVertices)); 65df1a93feSDuncan Campbell PetscCall(DMNetworkGetGlobalVertexIndex(dm, edgeVertices[0], &globalEdgeVertices[0])); 66df1a93feSDuncan Campbell PetscCall(DMNetworkGetGlobalVertexIndex(dm, edgeVertices[1], &globalEdgeVertices[1])); 67df1a93feSDuncan Campbell PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edge, &edge)); 68df1a93feSDuncan Campbell 69df1a93feSDuncan Campbell // TODO: Determine edge color/name 70df1a93feSDuncan 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)); 71df1a93feSDuncan Campbell } 72df1a93feSDuncan Campbell } 73df1a93feSDuncan Campbell // End synchronized printing 74df1a93feSDuncan Campbell PetscCall(PetscViewerFlush(viewer)); 75df1a93feSDuncan Campbell PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 76df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 77df1a93feSDuncan Campbell } 78df1a93feSDuncan Campbell 79df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_Matplotlib(DM dm, PetscViewer viewer) 80df1a93feSDuncan Campbell { 81cd2bb8e3SDuncan Campbell PetscMPIInt rank, size; 82df1a93feSDuncan Campbell MPI_Comm comm; 83*5f25b224SDuncan Campbell char filename[PETSC_MAX_PATH_LEN + 1], options[512], proccall[PETSC_MAX_PATH_LEN + 512], scriptFile[PETSC_MAX_PATH_LEN + 1], buffer[256]; 84df1a93feSDuncan Campbell PetscViewer csvViewer; 85df1a93feSDuncan Campbell FILE *processFile = NULL; 86*5f25b224SDuncan Campbell PetscBool isnull, optionShowRanks = PETSC_FALSE, optionRankIsSet = PETSC_FALSE, showNoNodes = PETSC_FALSE, showNoNumbering = PETSC_FALSE; 87df1a93feSDuncan Campbell PetscDraw draw; 88*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 89*5f25b224SDuncan Campbell PetscReal drawPause; 90*5f25b224SDuncan Campbell PetscInt i; 91cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP) 92cd2bb8e3SDuncan Campbell PetscBool isSharedTmp; 93cd2bb8e3SDuncan Campbell #endif 94df1a93feSDuncan Campbell 95df1a93feSDuncan Campbell PetscFunctionBegin; 96df1a93feSDuncan Campbell // Deal with the PetscDraw we are given 97df1a93feSDuncan Campbell PetscCall(PetscViewerDrawGetDraw(viewer, 1, &draw)); 98df1a93feSDuncan Campbell PetscCall(PetscDrawIsNull(draw, &isnull)); 99df1a93feSDuncan Campbell PetscCall(PetscDrawSetVisible(draw, PETSC_FALSE)); 100df1a93feSDuncan Campbell 101df1a93feSDuncan Campbell // Clear the file name buffer so all communicated bytes are well-defined 102df1a93feSDuncan Campbell PetscCall(PetscMemzero(filename, sizeof(filename))); 103df1a93feSDuncan Campbell 104df1a93feSDuncan Campbell // Get the MPI communicator and this process' rank 105df1a93feSDuncan Campbell PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 106df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(comm, &rank)); 107df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_size(comm, &size)); 108df1a93feSDuncan Campbell 109cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP) 110cd2bb8e3SDuncan Campbell // Get if the temporary directory is shared 111cd2bb8e3SDuncan Campbell // Note: This must be done collectively on every rank, it cannot be done on a single rank 112cd2bb8e3SDuncan Campbell PetscCall(PetscSharedTmp(comm, &isSharedTmp)); 113cd2bb8e3SDuncan Campbell #endif 114cd2bb8e3SDuncan Campbell 115*5f25b224SDuncan Campbell /* Process Options */ 116*5f25b224SDuncan Campbell optionShowRanks = network->vieweroptions.showallranks; 117*5f25b224SDuncan Campbell showNoNodes = network->vieweroptions.shownovertices; 118*5f25b224SDuncan Campbell showNoNumbering = network->vieweroptions.shownonumbering; 119*5f25b224SDuncan Campbell 120*5f25b224SDuncan Campbell /* 121*5f25b224SDuncan Campbell TODO: if the option -dmnetwork_view_tmpdir can be moved up here that would be good as well. 122*5f25b224SDuncan Campbell */ 123*5f25b224SDuncan Campbell PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "MatPlotLib PetscViewer DMNetwork Options", "PetscViewer"); 124*5f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_all_ranks", "View all ranks in the DMNetwork", NULL, optionShowRanks, &optionShowRanks, NULL)); 125*5f25b224SDuncan Campbell PetscCall(PetscOptionsString("-dmnetwork_view_rank_range", "Set of ranks to view the DMNetwork on", NULL, buffer, buffer, sizeof(buffer), &optionRankIsSet)); 126*5f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_no_vertices", "Do not view vertices", NULL, showNoNodes, &showNoNodes, NULL)); 127*5f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_no_numbering", "Do not view edge and vertex numbering", NULL, showNoNumbering, &showNoNumbering, NULL)); 128*5f25b224SDuncan Campbell PetscOptionsEnd(); 129*5f25b224SDuncan Campbell 130df1a93feSDuncan Campbell // Generate and broadcast the temporary file name from rank 0 131df1a93feSDuncan Campbell if (rank == 0) { 132df1a93feSDuncan Campbell #if defined(PETSC_HAVE_TMPNAM_S) 133df1a93feSDuncan Campbell // Acquire a temporary file to write to and open an ASCII/CSV viewer 134df1a93feSDuncan Campbell PetscCheck(tmpnam_s(filename, sizeof(filename)) == 0, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 13512e0ef65SDuncan Campbell #elif defined(PETSC_HAVE_MKSTEMP) 136cd2bb8e3SDuncan Campbell PetscBool isTmpOverridden; 137df1a93feSDuncan Campbell size_t numChars; 138df1a93feSDuncan Campbell // Same thing, but for POSIX systems on which tmpnam is deprecated 139df1a93feSDuncan 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 140df1a93feSDuncan 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 1416e4289a0SDuncan Campbell PetscCall(PetscOptionsGetString(NULL, NULL, "-dmnetwork_view_tmpdir", filename, sizeof(filename), &isTmpOverridden)); 1426e4289a0SDuncan Campbell // If not specified by option try using a shared tmp on the system 1436e4289a0SDuncan Campbell if (!isTmpOverridden) { 1446e4289a0SDuncan Campbell // Validate that if tmp is not overridden it is at least shared 145cd2bb8e3SDuncan 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"); 146cd2bb8e3SDuncan Campbell PetscCall(PetscGetTmp(PETSC_COMM_SELF, filename, sizeof(filename))); 147cd2bb8e3SDuncan Campbell } 148df1a93feSDuncan Campbell // Make sure the filename ends with a '/' 149df1a93feSDuncan Campbell PetscCall(PetscStrlen(filename, &numChars)); 150df1a93feSDuncan Campbell if (filename[numChars - 1] != '/') { 151df1a93feSDuncan Campbell filename[numChars] = '/'; 152df1a93feSDuncan Campbell filename[numChars + 1] = 0; 153df1a93feSDuncan Campbell } 154df1a93feSDuncan Campbell // Perform the actual temporary file creation 155c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(filename, "XXXXXX", sizeof(filename))); 156df1a93feSDuncan Campbell PetscCheck(mkstemp(filename) != -1, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 157df1a93feSDuncan Campbell #else 158df1a93feSDuncan Campbell // Same thing, but for older C versions which don't have the safe form 159df1a93feSDuncan Campbell PetscCheck(tmpnam(filename) != NULL, comm, PETSC_ERR_SYS, "Could not acquire temporary file"); 160df1a93feSDuncan Campbell #endif 161df1a93feSDuncan Campbell } 162df1a93feSDuncan Campbell 163cd2bb8e3SDuncan Campbell // Broadcast the filename to all other MPI ranks 164cd2bb8e3SDuncan Campbell PetscCallMPI(MPI_Bcast(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, 0, comm)); 165cd2bb8e3SDuncan Campbell 166e66d8692SHong Zhang PetscCall(PetscViewerASCIIOpen(comm, filename, &csvViewer)); 167df1a93feSDuncan Campbell PetscCall(PetscViewerPushFormat(csvViewer, PETSC_VIEWER_ASCII_CSV)); 168df1a93feSDuncan Campbell 169df1a93feSDuncan Campbell // Use the CSV viewer to write out the local network 170df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, csvViewer)); 171df1a93feSDuncan Campbell 172df1a93feSDuncan Campbell // Close the viewer 173df1a93feSDuncan Campbell PetscCall(PetscViewerDestroy(&csvViewer)); 174df1a93feSDuncan Campbell 175*5f25b224SDuncan Campbell // Generate options string 176*5f25b224SDuncan Campbell PetscCall(PetscMemzero(options, sizeof(options))); 177*5f25b224SDuncan Campbell // If the draw is null run as a "test execute" ie. do nothing just test that the script was called correctly 178*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, isnull ? " -tx " : " ", sizeof(options))); 179*5f25b224SDuncan Campbell PetscCall(PetscDrawGetPause(draw, &drawPause)); 180*5f25b224SDuncan Campbell if (drawPause > 0) { 181*5f25b224SDuncan Campbell char pausebuffer[64]; 182*5f25b224SDuncan Campbell PetscCall(PetscSNPrintf(pausebuffer, sizeof(pausebuffer), "%f", (double)drawPause)); 183*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -dt ", sizeof(options))); 184*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, pausebuffer, sizeof(options))); 185*5f25b224SDuncan Campbell } 186*5f25b224SDuncan Campbell if (optionShowRanks || optionRankIsSet) { 187*5f25b224SDuncan 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 188*5f25b224SDuncan Campbell if (optionShowRanks && !optionRankIsSet && size != 1) PetscCall(PetscStrlcat(options, " -dar ", sizeof(options))); 189*5f25b224SDuncan Campbell // Do not show the global plot if the user requests it OR if one specific rank is requested 190*5f25b224SDuncan Campbell if (network->vieweroptions.dontshowglobal || optionRankIsSet) PetscCall(PetscStrlcat(options, " -ncp ", sizeof(options))); 191*5f25b224SDuncan Campbell 192*5f25b224SDuncan Campbell if (optionRankIsSet) { 193*5f25b224SDuncan Campbell // If a range of ranks to draw is specified append it 194*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options))); 195*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, buffer, sizeof(options))); 196*5f25b224SDuncan Campbell } else { 197*5f25b224SDuncan Campbell // Otherwise, use the options provided in code 198*5f25b224SDuncan Campbell if (network->vieweroptions.viewranks) { 199*5f25b224SDuncan Campbell const PetscInt *viewranks; 200*5f25b224SDuncan Campbell PetscInt viewrankssize; 201*5f25b224SDuncan Campbell char rankbuffer[64]; 202*5f25b224SDuncan Campbell PetscCall(ISGetTotalIndices(network->vieweroptions.viewranks, &viewranks)); 203*5f25b224SDuncan Campbell PetscCall(ISGetSize(network->vieweroptions.viewranks, &viewrankssize)); 204*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options))); 205*5f25b224SDuncan Campbell for (i = 0; i < viewrankssize; i++) { 206*5f25b224SDuncan Campbell PetscCall(PetscSNPrintf(rankbuffer, sizeof(rankbuffer), "%" PetscInt_FMT, viewranks[i])); 207*5f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, rankbuffer, sizeof(options))); 208*5f25b224SDuncan Campbell } 209*5f25b224SDuncan Campbell PetscCall(ISRestoreTotalIndices(network->vieweroptions.viewranks, &viewranks)); 210*5f25b224SDuncan Campbell } // if not provided an IS of viewing ranks, skip viewing 211*5f25b224SDuncan Campbell } 212*5f25b224SDuncan Campbell } 213*5f25b224SDuncan Campbell 214*5f25b224SDuncan Campbell // Check for options for visibility... 215*5f25b224SDuncan Campbell if (showNoNodes) PetscCall(PetscStrlcat(options, " -nn ", sizeof(options))); 216*5f25b224SDuncan Campbell if (showNoNumbering) PetscCall(PetscStrlcat(options, " -nnl -nel ", sizeof(options))); 217*5f25b224SDuncan Campbell 218df1a93feSDuncan Campbell // Get the value of $PETSC_DIR 219e66d8692SHong Zhang PetscCall(PetscStrreplace(comm, "${PETSC_DIR}/share/petsc/bin/dmnetwork_view.py", scriptFile, sizeof(scriptFile))); 220df1a93feSDuncan Campbell PetscCall(PetscFixFilename(scriptFile, scriptFile)); 221*5f25b224SDuncan Campbell // Generate the system call for 'python3 $PETSC_DIR/share/petsc/dmnetwork_view.py <options> <file>' 222df1a93feSDuncan Campbell PetscCall(PetscArrayzero(proccall, sizeof(proccall))); 223*5f25b224SDuncan Campbell PetscCall(PetscSNPrintf(proccall, sizeof(proccall), "%s %s %s %s", PETSC_PYTHON_EXE, scriptFile, options, filename)); 224df1a93feSDuncan Campbell 225df1a93feSDuncan Campbell #if defined(PETSC_HAVE_POPEN) 226df1a93feSDuncan Campbell // Perform the call to run the python script (Note: while this is called on all ranks POpen will only run on rank 0) 227e66d8692SHong Zhang PetscCall(PetscPOpen(comm, NULL, proccall, "r", &processFile)); 228df1a93feSDuncan Campbell if (processFile != NULL) { 229*5f25b224SDuncan Campbell while (fgets(buffer, sizeof(buffer), processFile) != NULL) PetscCall(PetscPrintf(comm, "%s", buffer)); 230df1a93feSDuncan Campbell } 231e66d8692SHong Zhang PetscCall(PetscPClose(comm, processFile)); 232df1a93feSDuncan Campbell #else 233df1a93feSDuncan Campbell // Same thing, but using the standard library for systems that don't have POpen/PClose (only run on rank 0) 234e66d8692SHong Zhang if (rank == 0) PetscCheck(system(proccall) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to call viewer script"); 235df1a93feSDuncan Campbell // Barrier so that all ranks wait until the call completes 236e66d8692SHong Zhang PetscCallMPI(MPI_Barrier(comm)); 237df1a93feSDuncan Campbell #endif 238df1a93feSDuncan Campbell // Clean up the temporary file we used using rank 0 239e66d8692SHong Zhang if (rank == 0) PetscCheck(remove(filename) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to delete temporary file"); 240df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 241df1a93feSDuncan Campbell } 242df1a93feSDuncan Campbell 243d2fd5932SHong Zhang PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 244d2fd5932SHong Zhang { 245df1a93feSDuncan Campbell PetscBool iascii, isdraw; 246df1a93feSDuncan Campbell PetscViewerFormat format; 247d2fd5932SHong Zhang 248d2fd5932SHong Zhang PetscFunctionBegin; 249d2fd5932SHong Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 250d2fd5932SHong Zhang PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 251df1a93feSDuncan Campbell PetscCall(PetscViewerGetFormat(viewer, &format)); 252df1a93feSDuncan Campbell 253df1a93feSDuncan Campbell PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 254df1a93feSDuncan Campbell if (isdraw) { 255df1a93feSDuncan Campbell PetscCall(DMView_Network_Matplotlib(dm, viewer)); 256df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 257df1a93feSDuncan Campbell } 258df1a93feSDuncan Campbell 259d2fd5932SHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 260d2fd5932SHong Zhang if (iascii) { 261d2fd5932SHong Zhang const PetscInt *cone, *vtx, *edges; 262d2fd5932SHong Zhang PetscInt vfrom, vto, i, j, nv, ne, nsv, p, nsubnet; 263d2fd5932SHong Zhang DM_Network *network = (DM_Network *)dm->data; 264df1a93feSDuncan Campbell PetscMPIInt rank; 265df1a93feSDuncan Campbell 266df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 267df1a93feSDuncan Campbell if (format == PETSC_VIEWER_ASCII_CSV) { 268df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, viewer)); 269df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 270df1a93feSDuncan Campbell } 271d2fd5932SHong Zhang 272d2fd5932SHong Zhang nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */ 273df1a93feSDuncan Campbell if (!rank) { 274d2fd5932SHong 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, 275d2fd5932SHong Zhang network->cloneshared->Nsvtx)); 276d2fd5932SHong Zhang } 277d2fd5932SHong Zhang 278d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL)); 279d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 280d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv)); 281d2fd5932SHong Zhang 282d2fd5932SHong Zhang for (i = 0; i < nsubnet; i++) { 283d2fd5932SHong Zhang PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges)); 284d2fd5932SHong Zhang if (ne) { 285d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv)); 286d2fd5932SHong Zhang for (j = 0; j < ne; j++) { 287d2fd5932SHong Zhang p = edges[j]; 288d2fd5932SHong Zhang PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone)); 289d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom)); 290d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto)); 291d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p)); 292d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto)); 293d2fd5932SHong Zhang } 294d2fd5932SHong Zhang } 295d2fd5932SHong Zhang } 296d2fd5932SHong Zhang 297d2fd5932SHong Zhang /* Shared vertices */ 298d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx)); 299d2fd5932SHong Zhang if (nsv) { 300d2fd5932SHong Zhang PetscInt gidx; 301d2fd5932SHong Zhang PetscBool ghost; 302d2fd5932SHong Zhang const PetscInt *sv = NULL; 303d2fd5932SHong Zhang 304d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n")); 305d2fd5932SHong Zhang for (i = 0; i < nsv; i++) { 306d2fd5932SHong Zhang PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost)); 307d2fd5932SHong Zhang if (ghost) continue; 308d2fd5932SHong Zhang 309d2fd5932SHong Zhang PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv)); 310d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1])); 311d2fd5932SHong Zhang for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1])); 312d2fd5932SHong Zhang } 313d2fd5932SHong Zhang } 314d2fd5932SHong Zhang PetscCall(PetscViewerFlush(viewer)); 315d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 316d2fd5932SHong Zhang } else PetscCheck(iascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name); 317d2fd5932SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 318d2fd5932SHong Zhang } 319*5f25b224SDuncan Campbell 320*5f25b224SDuncan Campbell /*@ 321*5f25b224SDuncan Campbell DMNetworkViewSetShowRanks - Sets viewing the `DMETNWORK` on each rank individually. 322*5f25b224SDuncan Campbell 323*5f25b224SDuncan Campbell Logically Collective 324*5f25b224SDuncan Campbell 325*5f25b224SDuncan Campbell Input Parameter: 326*5f25b224SDuncan Campbell . dm - the `DMNETWORK` object 327*5f25b224SDuncan Campbell 328*5f25b224SDuncan Campbell Output Parameter: 329*5f25b224SDuncan Campbell . showranks - `PETSC_TRUE` if viewing each rank's sub network individually 330*5f25b224SDuncan Campbell 331*5f25b224SDuncan Campbell Level: beginner 332*5f25b224SDuncan Campbell 333*5f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 334*5f25b224SDuncan Campbell @*/ 335*5f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowRanks(DM dm, PetscBool showranks) 336*5f25b224SDuncan Campbell { 337*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 338*5f25b224SDuncan Campbell 339*5f25b224SDuncan Campbell PetscFunctionBegin; 340*5f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 341*5f25b224SDuncan Campbell network->vieweroptions.showallranks = showranks; 342*5f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 343*5f25b224SDuncan Campbell } 344*5f25b224SDuncan Campbell 345*5f25b224SDuncan Campbell /*@ 346*5f25b224SDuncan Campbell DMNetworkViewSetShowGlobal - Set viewing the global network. 347*5f25b224SDuncan Campbell 348*5f25b224SDuncan Campbell Logically Collective 349*5f25b224SDuncan Campbell 350*5f25b224SDuncan Campbell Input Parameter: 351*5f25b224SDuncan Campbell . dm - the `DMNETWORK` object 352*5f25b224SDuncan Campbell 353*5f25b224SDuncan Campbell Output Parameter: 354*5f25b224SDuncan Campbell . showglobal - `PETSC_TRUE` if viewing the global network 355*5f25b224SDuncan Campbell 356*5f25b224SDuncan Campbell Level: beginner 357*5f25b224SDuncan Campbell 358*5f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 359*5f25b224SDuncan Campbell @*/ 360*5f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowGlobal(DM dm, PetscBool showglobal) 361*5f25b224SDuncan Campbell { 362*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 363*5f25b224SDuncan Campbell 364*5f25b224SDuncan Campbell PetscFunctionBegin; 365*5f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 366*5f25b224SDuncan Campbell network->vieweroptions.dontshowglobal = (PetscBool)(!showglobal); 367*5f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 368*5f25b224SDuncan Campbell } 369*5f25b224SDuncan Campbell 370*5f25b224SDuncan Campbell /*@ 371*5f25b224SDuncan Campbell DMNetworkViewSetShowVertices - Sets whether to display the vertices in viewing routines. 372*5f25b224SDuncan Campbell 373*5f25b224SDuncan Campbell Logically Collective 374*5f25b224SDuncan Campbell 375*5f25b224SDuncan Campbell Input Parameter: 376*5f25b224SDuncan Campbell . dm - the `DMNETWORK` object 377*5f25b224SDuncan Campbell 378*5f25b224SDuncan Campbell Output Parameter: 379*5f25b224SDuncan Campbell . showvertices - `PETSC_TRUE` if visualizing the vertices 380*5f25b224SDuncan Campbell 381*5f25b224SDuncan Campbell Level: beginner 382*5f25b224SDuncan Campbell 383*5f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()` 384*5f25b224SDuncan Campbell @*/ 385*5f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowVertices(DM dm, PetscBool showvertices) 386*5f25b224SDuncan Campbell { 387*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 388*5f25b224SDuncan Campbell 389*5f25b224SDuncan Campbell PetscFunctionBegin; 390*5f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 391*5f25b224SDuncan Campbell network->vieweroptions.shownovertices = (PetscBool)(!showvertices); 392*5f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 393*5f25b224SDuncan Campbell } 394*5f25b224SDuncan Campbell 395*5f25b224SDuncan Campbell /*@ 396*5f25b224SDuncan Campbell DMNetworkViewSetShowNumbering - Set displaying the numbering of edges and vertices in viewing routines. 397*5f25b224SDuncan Campbell 398*5f25b224SDuncan Campbell Logically Collective 399*5f25b224SDuncan Campbell 400*5f25b224SDuncan Campbell Input Parameter: 401*5f25b224SDuncan Campbell . dm - the `DMNETWORK` object 402*5f25b224SDuncan Campbell 403*5f25b224SDuncan Campbell Output Parameter: 404*5f25b224SDuncan Campbell . shownumbering - `PETSC_TRUE` if displaying the numbering of edges and vertices 405*5f25b224SDuncan Campbell 406*5f25b224SDuncan Campbell Level: beginner 407*5f25b224SDuncan Campbell 408*5f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetViewRanks()` 409*5f25b224SDuncan Campbell @*/ 410*5f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowNumbering(DM dm, PetscBool shownumbering) 411*5f25b224SDuncan Campbell { 412*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 413*5f25b224SDuncan Campbell 414*5f25b224SDuncan Campbell PetscFunctionBegin; 415*5f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 416*5f25b224SDuncan Campbell network->vieweroptions.shownonumbering = (PetscBool)(!shownumbering); 417*5f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 418*5f25b224SDuncan Campbell } 419*5f25b224SDuncan Campbell 420*5f25b224SDuncan Campbell /*@ 421*5f25b224SDuncan Campbell DMNetworkViewSetViewRanks - View the `DMNETWORK` on each of the specified ranks individually. 422*5f25b224SDuncan Campbell 423*5f25b224SDuncan Campbell Collective 424*5f25b224SDuncan Campbell 425*5f25b224SDuncan Campbell Input Parameter: 426*5f25b224SDuncan Campbell . dm - the `DMNETWORK` object 427*5f25b224SDuncan Campbell 428*5f25b224SDuncan Campbell Output Parameter: 429*5f25b224SDuncan Campbell . viewranks - set of ranks to view the `DMNETWORK` on individually 430*5f25b224SDuncan Campbell 431*5f25b224SDuncan Campbell Level: beginner 432*5f25b224SDuncan Campbell 433*5f25b224SDuncan Campbell Note: 434*5f25b224SDuncan Campbell `DMNetwork` takes ownership of the input viewranks `IS`, it should be destroyed by the caller. 435*5f25b224SDuncan Campbell 436*5f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()` 437*5f25b224SDuncan Campbell @*/ 438*5f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetViewRanks(DM dm, IS viewranks) 439*5f25b224SDuncan Campbell { 440*5f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data; 441*5f25b224SDuncan Campbell 442*5f25b224SDuncan Campbell PetscFunctionBegin; 443*5f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK); 444*5f25b224SDuncan Campbell PetscValidHeaderSpecific(viewranks, IS_CLASSID, 2); 445*5f25b224SDuncan Campbell PetscCheckSameComm(dm, 1, viewranks, 2); 446*5f25b224SDuncan Campbell PetscCall(ISDestroy(&network->vieweroptions.viewranks)); 447*5f25b224SDuncan Campbell PetscCall(PetscObjectReference((PetscObject)viewranks)); 448*5f25b224SDuncan Campbell network->vieweroptions.viewranks = viewranks; 449*5f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS); 450*5f25b224SDuncan Campbell } 451