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)
41777c8a5SBarry Smith #if !defined(_XOPEN_SOURCE)
512e0ef65SDuncan Campbell #define _XOPEN_SOURCE 600
612e0ef65SDuncan Campbell #endif
71777c8a5SBarry Smith #endif
8d2fd5932SHong Zhang #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/
912e0ef65SDuncan Campbell #include <petscdraw.h>
10d2fd5932SHong Zhang
DMView_Network_CSV(DM dm,PetscViewer viewer)11df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_CSV(DM dm, PetscViewer viewer)
12df1a93feSDuncan Campbell {
13df1a93feSDuncan Campbell DM dmcoords;
147cd471e7SHong Zhang PetscInt nsubnets, i, subnet, nvertices, nedges, vertex, edge, gidx, ncomp;
15df1a93feSDuncan Campbell PetscInt vertexOffsets[2], globalEdgeVertices[2];
167cd471e7SHong Zhang PetscScalar vertexCoords[2], *color_ptr, color;
17df1a93feSDuncan Campbell const PetscInt *vertices, *edges, *edgeVertices;
18df1a93feSDuncan Campbell Vec allVertexCoords;
19df1a93feSDuncan Campbell PetscMPIInt rank;
20df1a93feSDuncan Campbell MPI_Comm comm;
21df1a93feSDuncan Campbell
22df1a93feSDuncan Campbell PetscFunctionBegin;
237cd471e7SHong Zhang // Get the coordinate information from dmcoords
247cd471e7SHong Zhang PetscCheck(dm->coordinates[0].dm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "CoordinateDM not created");
25df1a93feSDuncan Campbell PetscCall(DMGetCoordinateDM(dm, &dmcoords));
267cd471e7SHong Zhang
277cd471e7SHong Zhang PetscCall(DMGetCoordinateDim(dmcoords, &i));
28baca6076SPierre Jolivet PetscCheck(i == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim %" PetscInt_FMT " != 2 is not supported yet", i);
297cd471e7SHong Zhang
307cd471e7SHong Zhang // Get the coordinate vector from dm
31df1a93feSDuncan Campbell PetscCall(DMGetCoordinatesLocal(dm, &allVertexCoords));
327cd471e7SHong Zhang
33df1a93feSDuncan Campbell // Get the MPI communicator and this process' rank
347cd471e7SHong Zhang PetscCall(PetscObjectGetComm((PetscObject)dmcoords, &comm));
35df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(comm, &rank));
367cd471e7SHong Zhang
37df1a93feSDuncan Campbell // Start synchronized printing
38df1a93feSDuncan Campbell PetscCall(PetscViewerASCIIPushSynchronized(viewer));
39df1a93feSDuncan Campbell
40df1a93feSDuncan Campbell // Write the header
41df1a93feSDuncan Campbell PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Type,Rank,ID,X,Y,Z,Name,Color\n"));
42df1a93feSDuncan Campbell
43df1a93feSDuncan Campbell // Iterate each subnetwork (Note: We need to get the global number of subnets apparently)
447cd471e7SHong Zhang PetscCall(DMNetworkGetNumSubNetworks(dmcoords, NULL, &nsubnets));
45df1a93feSDuncan Campbell for (subnet = 0; subnet < nsubnets; subnet++) {
46df1a93feSDuncan Campbell // Get the subnetwork's vertices and edges
477cd471e7SHong Zhang PetscCall(DMNetworkGetSubnetwork(dmcoords, subnet, &nvertices, &nedges, &vertices, &edges));
48df1a93feSDuncan Campbell
49df1a93feSDuncan Campbell // Write out each vertex
50df1a93feSDuncan Campbell for (i = 0; i < nvertices; i++) {
51df1a93feSDuncan Campbell vertex = vertices[i];
527cd471e7SHong Zhang
53df1a93feSDuncan Campbell // Get the offset into the coordinate vector for the vertex
54df1a93feSDuncan Campbell PetscCall(DMNetworkGetLocalVecOffset(dmcoords, vertex, ALL_COMPONENTS, vertexOffsets));
55df1a93feSDuncan Campbell vertexOffsets[1] = vertexOffsets[0] + 1;
56df1a93feSDuncan Campbell // Remap vertex to the global value
577cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, vertex, &gidx));
58df1a93feSDuncan Campbell // Get the vertex position from the coordinate vector
59df1a93feSDuncan Campbell PetscCall(VecGetValues(allVertexCoords, 2, vertexOffsets, vertexCoords));
60df1a93feSDuncan Campbell
617cd471e7SHong Zhang // Get vertex color; TODO: name
627cd471e7SHong Zhang PetscCall(DMNetworkGetNumComponents(dmcoords, vertex, &ncomp));
637cd471e7SHong Zhang PetscCheck(ncomp <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "num of components %" PetscInt_FMT " must be <= 1", ncomp);
647cd471e7SHong Zhang color = 0.0;
657cd471e7SHong Zhang if (ncomp == 1) {
667cd471e7SHong Zhang PetscCall(DMNetworkGetComponent(dmcoords, vertex, 0, NULL, (void **)&color_ptr, NULL));
677cd471e7SHong Zhang color = *color_ptr;
687cd471e7SHong Zhang }
69835f2295SStefano Zampini PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Node,%d,%" PetscInt_FMT ",%lf,%lf,0,%" PetscInt_FMT ",%lf\n", rank, gidx, (double)PetscRealPart(vertexCoords[0]), (double)PetscRealPart(vertexCoords[1]), gidx, (double)PetscRealPart(color)));
70df1a93feSDuncan Campbell }
71df1a93feSDuncan Campbell
72df1a93feSDuncan Campbell // Write out each edge
73df1a93feSDuncan Campbell for (i = 0; i < nedges; i++) {
74df1a93feSDuncan Campbell edge = edges[i];
757cd471e7SHong Zhang PetscCall(DMNetworkGetConnectedVertices(dmcoords, edge, &edgeVertices));
767cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[0], &globalEdgeVertices[0]));
777cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, edgeVertices[1], &globalEdgeVertices[1]));
787cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalEdgeIndex(dmcoords, edge, &edge));
79df1a93feSDuncan Campbell
80df1a93feSDuncan Campbell // TODO: Determine edge color/name
81835f2295SStefano Zampini PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Edge,%d,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",0,%" PetscInt_FMT "\n", rank, edge, globalEdgeVertices[0], globalEdgeVertices[1], edge));
82df1a93feSDuncan Campbell }
83df1a93feSDuncan Campbell }
84df1a93feSDuncan Campbell // End synchronized printing
85df1a93feSDuncan Campbell PetscCall(PetscViewerFlush(viewer));
86df1a93feSDuncan Campbell PetscCall(PetscViewerASCIIPopSynchronized(viewer));
87df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
88df1a93feSDuncan Campbell }
89df1a93feSDuncan Campbell
DMView_Network_Matplotlib(DM dm,PetscViewer viewer)90df1a93feSDuncan Campbell static PetscErrorCode DMView_Network_Matplotlib(DM dm, PetscViewer viewer)
91df1a93feSDuncan Campbell {
92cd2bb8e3SDuncan Campbell PetscMPIInt rank, size;
93df1a93feSDuncan Campbell MPI_Comm comm;
945039956bSHong 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];
95df1a93feSDuncan Campbell PetscViewer csvViewer;
96df1a93feSDuncan Campbell FILE *processFile = NULL;
975039956bSHong Zhang PetscBool isnull, optionShowRanks = PETSC_FALSE, optionRankIsSet = PETSC_FALSE, showNoNodes = PETSC_FALSE, showNoNumbering = PETSC_FALSE, optionShowVertices = PETSC_FALSE, optionViewPadding = PETSC_FALSE;
98df1a93feSDuncan Campbell PetscDraw draw;
995f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
1005039956bSHong Zhang PetscReal drawPause, viewPadding = 1.0;
1015f25b224SDuncan Campbell PetscInt i;
102cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP)
103cd2bb8e3SDuncan Campbell PetscBool isSharedTmp;
104cd2bb8e3SDuncan Campbell #endif
105df1a93feSDuncan Campbell
106df1a93feSDuncan Campbell PetscFunctionBegin;
107df1a93feSDuncan Campbell // Deal with the PetscDraw we are given
108df1a93feSDuncan Campbell PetscCall(PetscViewerDrawGetDraw(viewer, 1, &draw));
109df1a93feSDuncan Campbell PetscCall(PetscDrawIsNull(draw, &isnull));
110df1a93feSDuncan Campbell PetscCall(PetscDrawSetVisible(draw, PETSC_FALSE));
111df1a93feSDuncan Campbell
112df1a93feSDuncan Campbell // Clear the file name buffer so all communicated bytes are well-defined
113df1a93feSDuncan Campbell PetscCall(PetscMemzero(filename, sizeof(filename)));
114df1a93feSDuncan Campbell
115df1a93feSDuncan Campbell // Get the MPI communicator and this process' rank
116df1a93feSDuncan Campbell PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
117df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(comm, &rank));
118df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_size(comm, &size));
119df1a93feSDuncan Campbell
120cd2bb8e3SDuncan Campbell #if defined(PETSC_HAVE_MKSTEMP)
121cd2bb8e3SDuncan Campbell // Get if the temporary directory is shared
122cd2bb8e3SDuncan Campbell // Note: This must be done collectively on every rank, it cannot be done on a single rank
123cd2bb8e3SDuncan Campbell PetscCall(PetscSharedTmp(comm, &isSharedTmp));
124cd2bb8e3SDuncan Campbell #endif
125cd2bb8e3SDuncan Campbell
1265f25b224SDuncan Campbell /* Process Options */
1275f25b224SDuncan Campbell optionShowRanks = network->vieweroptions.showallranks;
1285f25b224SDuncan Campbell showNoNodes = network->vieweroptions.shownovertices;
1295f25b224SDuncan Campbell showNoNumbering = network->vieweroptions.shownonumbering;
1305f25b224SDuncan Campbell
1315f25b224SDuncan Campbell /*
1325f25b224SDuncan Campbell TODO: if the option -dmnetwork_view_tmpdir can be moved up here that would be good as well.
1335f25b224SDuncan Campbell */
1345f25b224SDuncan Campbell PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "MatPlotLib PetscViewer DMNetwork Options", "PetscViewer");
1355f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_all_ranks", "View all ranks in the DMNetwork", NULL, optionShowRanks, &optionShowRanks, NULL));
136*0b595cf8SStefano Zampini PetscCall(PetscOptionsString("-dmnetwork_view_rank_range", "Set of ranks to view the DMNetwork on", NULL, NULL, buffer, sizeof(buffer), &optionRankIsSet));
1375f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_no_vertices", "Do not view vertices", NULL, showNoNodes, &showNoNodes, NULL));
1385f25b224SDuncan Campbell PetscCall(PetscOptionsBool("-dmnetwork_view_no_numbering", "Do not view edge and vertex numbering", NULL, showNoNumbering, &showNoNumbering, NULL));
139*0b595cf8SStefano Zampini PetscCall(PetscOptionsString("-dmnetwork_view_zoomin_vertices", "Focus the view on the given set of vertices", NULL, NULL, buffer2, sizeof(buffer2), &optionShowVertices));
1405039956bSHong Zhang PetscCall(PetscOptionsReal("-dmnetwork_view_zoomin_vertices_padding", "Set the padding when viewing specific vertices", NULL, viewPadding, &viewPadding, &optionViewPadding));
1415f25b224SDuncan Campbell PetscOptionsEnd();
1425f25b224SDuncan Campbell
143df1a93feSDuncan Campbell // Generate and broadcast the temporary file name from rank 0
144df1a93feSDuncan Campbell if (rank == 0) {
145df1a93feSDuncan Campbell #if defined(PETSC_HAVE_TMPNAM_S)
146df1a93feSDuncan Campbell // Acquire a temporary file to write to and open an ASCII/CSV viewer
147df1a93feSDuncan Campbell PetscCheck(tmpnam_s(filename, sizeof(filename)) == 0, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
14812e0ef65SDuncan Campbell #elif defined(PETSC_HAVE_MKSTEMP)
149cd2bb8e3SDuncan Campbell PetscBool isTmpOverridden;
150df1a93feSDuncan Campbell size_t numChars;
151df1a93feSDuncan Campbell // Same thing, but for POSIX systems on which tmpnam is deprecated
152df1a93feSDuncan 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
153df1a93feSDuncan 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
1546e4289a0SDuncan Campbell PetscCall(PetscOptionsGetString(NULL, NULL, "-dmnetwork_view_tmpdir", filename, sizeof(filename), &isTmpOverridden));
1556e4289a0SDuncan Campbell // If not specified by option try using a shared tmp on the system
1566e4289a0SDuncan Campbell if (!isTmpOverridden) {
1576e4289a0SDuncan Campbell // Validate that if tmp is not overridden it is at least shared
158cd2bb8e3SDuncan 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");
159cd2bb8e3SDuncan Campbell PetscCall(PetscGetTmp(PETSC_COMM_SELF, filename, sizeof(filename)));
160cd2bb8e3SDuncan Campbell }
161df1a93feSDuncan Campbell // Make sure the filename ends with a '/'
162df1a93feSDuncan Campbell PetscCall(PetscStrlen(filename, &numChars));
163df1a93feSDuncan Campbell if (filename[numChars - 1] != '/') {
164df1a93feSDuncan Campbell filename[numChars] = '/';
165df1a93feSDuncan Campbell filename[numChars + 1] = 0;
166df1a93feSDuncan Campbell }
167df1a93feSDuncan Campbell // Perform the actual temporary file creation
168c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(filename, "XXXXXX", sizeof(filename)));
169df1a93feSDuncan Campbell PetscCheck(mkstemp(filename) != -1, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
170df1a93feSDuncan Campbell #else
171df1a93feSDuncan Campbell // Same thing, but for older C versions which don't have the safe form
172df1a93feSDuncan Campbell PetscCheck(tmpnam(filename) != NULL, comm, PETSC_ERR_SYS, "Could not acquire temporary file");
173df1a93feSDuncan Campbell #endif
174df1a93feSDuncan Campbell }
175df1a93feSDuncan Campbell
176cd2bb8e3SDuncan Campbell // Broadcast the filename to all other MPI ranks
177cd2bb8e3SDuncan Campbell PetscCallMPI(MPI_Bcast(filename, PETSC_MAX_PATH_LEN, MPI_BYTE, 0, comm));
178cd2bb8e3SDuncan Campbell
179e66d8692SHong Zhang PetscCall(PetscViewerASCIIOpen(comm, filename, &csvViewer));
180df1a93feSDuncan Campbell PetscCall(PetscViewerPushFormat(csvViewer, PETSC_VIEWER_ASCII_CSV));
181df1a93feSDuncan Campbell
182df1a93feSDuncan Campbell // Use the CSV viewer to write out the local network
183df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, csvViewer));
184df1a93feSDuncan Campbell
185df1a93feSDuncan Campbell // Close the viewer
186df1a93feSDuncan Campbell PetscCall(PetscViewerDestroy(&csvViewer));
187df1a93feSDuncan Campbell
1885f25b224SDuncan Campbell // Generate options string
1895f25b224SDuncan Campbell PetscCall(PetscMemzero(options, sizeof(options)));
1905f25b224SDuncan Campbell // If the draw is null run as a "test execute" ie. do nothing just test that the script was called correctly
1915f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, isnull ? " -tx " : " ", sizeof(options)));
1925f25b224SDuncan Campbell PetscCall(PetscDrawGetPause(draw, &drawPause));
1935f25b224SDuncan Campbell if (drawPause > 0) {
1945f25b224SDuncan Campbell char pausebuffer[64];
1955f25b224SDuncan Campbell PetscCall(PetscSNPrintf(pausebuffer, sizeof(pausebuffer), "%f", (double)drawPause));
1965f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -dt ", sizeof(options)));
1975f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, pausebuffer, sizeof(options)));
1985f25b224SDuncan Campbell }
1995f25b224SDuncan Campbell if (optionShowRanks || optionRankIsSet) {
2005f25b224SDuncan 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
2015f25b224SDuncan Campbell if (optionShowRanks && !optionRankIsSet && size != 1) PetscCall(PetscStrlcat(options, " -dar ", sizeof(options)));
2025f25b224SDuncan Campbell // Do not show the global plot if the user requests it OR if one specific rank is requested
2035f25b224SDuncan Campbell if (network->vieweroptions.dontshowglobal || optionRankIsSet) PetscCall(PetscStrlcat(options, " -ncp ", sizeof(options)));
2045f25b224SDuncan Campbell
2055f25b224SDuncan Campbell if (optionRankIsSet) {
2065f25b224SDuncan Campbell // If a range of ranks to draw is specified append it
2075f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options)));
2085f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, buffer, sizeof(options)));
2095f25b224SDuncan Campbell } else {
2105f25b224SDuncan Campbell // Otherwise, use the options provided in code
2115f25b224SDuncan Campbell if (network->vieweroptions.viewranks) {
2125f25b224SDuncan Campbell const PetscInt *viewranks;
2135f25b224SDuncan Campbell PetscInt viewrankssize;
2145f25b224SDuncan Campbell char rankbuffer[64];
2155f25b224SDuncan Campbell PetscCall(ISGetTotalIndices(network->vieweroptions.viewranks, &viewranks));
2165f25b224SDuncan Campbell PetscCall(ISGetSize(network->vieweroptions.viewranks, &viewrankssize));
2175f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, " -drr ", sizeof(options)));
2185f25b224SDuncan Campbell for (i = 0; i < viewrankssize; i++) {
2195f25b224SDuncan Campbell PetscCall(PetscSNPrintf(rankbuffer, sizeof(rankbuffer), "%" PetscInt_FMT, viewranks[i]));
2205f25b224SDuncan Campbell PetscCall(PetscStrlcat(options, rankbuffer, sizeof(options)));
2215f25b224SDuncan Campbell }
2225f25b224SDuncan Campbell PetscCall(ISRestoreTotalIndices(network->vieweroptions.viewranks, &viewranks));
2235f25b224SDuncan Campbell } // if not provided an IS of viewing ranks, skip viewing
2245f25b224SDuncan Campbell }
2255f25b224SDuncan Campbell }
2265039956bSHong Zhang if (optionShowVertices) {
2275039956bSHong Zhang // Pass vertices to focus on if defined
2285039956bSHong Zhang PetscCall(PetscStrlcat(options, " -vsv ", sizeof(options)));
2295039956bSHong Zhang PetscCall(PetscStrlcat(options, buffer2, sizeof(options)));
2305039956bSHong Zhang optionViewPadding = PETSC_TRUE;
2315039956bSHong Zhang // Pass padding if set
2325039956bSHong Zhang if (optionViewPadding) {
2335039956bSHong Zhang PetscCall(PetscSNPrintf(buffer2, sizeof(buffer2), "%f", (double)viewPadding));
2345039956bSHong Zhang PetscCall(PetscStrlcat(options, " -vp ", sizeof(options)));
2355039956bSHong Zhang PetscCall(PetscStrlcat(options, buffer2, sizeof(options)));
2365039956bSHong Zhang }
2375039956bSHong Zhang }
2385f25b224SDuncan Campbell
2395f25b224SDuncan Campbell // Check for options for visibility...
2405f25b224SDuncan Campbell if (showNoNodes) PetscCall(PetscStrlcat(options, " -nn ", sizeof(options)));
2415f25b224SDuncan Campbell if (showNoNumbering) PetscCall(PetscStrlcat(options, " -nnl -nel ", sizeof(options)));
2425f25b224SDuncan Campbell
243df1a93feSDuncan Campbell // Get the value of $PETSC_DIR
244e66d8692SHong Zhang PetscCall(PetscStrreplace(comm, "${PETSC_DIR}/share/petsc/bin/dmnetwork_view.py", scriptFile, sizeof(scriptFile)));
245df1a93feSDuncan Campbell PetscCall(PetscFixFilename(scriptFile, scriptFile));
2465f25b224SDuncan Campbell // Generate the system call for 'python3 $PETSC_DIR/share/petsc/dmnetwork_view.py <options> <file>'
247df1a93feSDuncan Campbell PetscCall(PetscArrayzero(proccall, sizeof(proccall)));
2485f25b224SDuncan Campbell PetscCall(PetscSNPrintf(proccall, sizeof(proccall), "%s %s %s %s", PETSC_PYTHON_EXE, scriptFile, options, filename));
249df1a93feSDuncan Campbell
250df1a93feSDuncan Campbell #if defined(PETSC_HAVE_POPEN)
251df1a93feSDuncan Campbell // Perform the call to run the python script (Note: while this is called on all ranks POpen will only run on rank 0)
252e66d8692SHong Zhang PetscCall(PetscPOpen(comm, NULL, proccall, "r", &processFile));
253df1a93feSDuncan Campbell if (processFile != NULL) {
2545f25b224SDuncan Campbell while (fgets(buffer, sizeof(buffer), processFile) != NULL) PetscCall(PetscPrintf(comm, "%s", buffer));
255df1a93feSDuncan Campbell }
256e66d8692SHong Zhang PetscCall(PetscPClose(comm, processFile));
257df1a93feSDuncan Campbell #else
258df1a93feSDuncan Campbell // Same thing, but using the standard library for systems that don't have POpen/PClose (only run on rank 0)
259e66d8692SHong Zhang if (rank == 0) PetscCheck(system(proccall) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to call viewer script");
260df1a93feSDuncan Campbell // Barrier so that all ranks wait until the call completes
261e66d8692SHong Zhang PetscCallMPI(MPI_Barrier(comm));
262df1a93feSDuncan Campbell #endif
263df1a93feSDuncan Campbell // Clean up the temporary file we used using rank 0
264e66d8692SHong Zhang if (rank == 0) PetscCheck(remove(filename) == 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "Failed to delete temporary file");
265df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
266df1a93feSDuncan Campbell }
267df1a93feSDuncan Campbell
DMView_Network(DM dm,PetscViewer viewer)268d2fd5932SHong Zhang PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
269d2fd5932SHong Zhang {
2709f196a02SMartin Diehl PetscBool isascii, isdraw;
271df1a93feSDuncan Campbell PetscViewerFormat format;
272d2fd5932SHong Zhang
273d2fd5932SHong Zhang PetscFunctionBegin;
274d2fd5932SHong Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
275d2fd5932SHong Zhang PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
276df1a93feSDuncan Campbell PetscCall(PetscViewerGetFormat(viewer, &format));
277df1a93feSDuncan Campbell
278df1a93feSDuncan Campbell PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
279df1a93feSDuncan Campbell if (isdraw) {
280df1a93feSDuncan Campbell PetscCall(DMView_Network_Matplotlib(dm, viewer));
281df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
282df1a93feSDuncan Campbell }
283df1a93feSDuncan Campbell
2849f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2859f196a02SMartin Diehl if (isascii) {
286d2fd5932SHong Zhang const PetscInt *cone, *vtx, *edges;
287d2fd5932SHong Zhang PetscInt vfrom, vto, i, j, nv, ne, nsv, p, nsubnet;
288d2fd5932SHong Zhang DM_Network *network = (DM_Network *)dm->data;
289df1a93feSDuncan Campbell PetscMPIInt rank;
290df1a93feSDuncan Campbell
291df1a93feSDuncan Campbell PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
292df1a93feSDuncan Campbell if (format == PETSC_VIEWER_ASCII_CSV) {
293df1a93feSDuncan Campbell PetscCall(DMView_Network_CSV(dm, viewer));
294df1a93feSDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
295df1a93feSDuncan Campbell }
296d2fd5932SHong Zhang
297d2fd5932SHong Zhang nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */
298df1a93feSDuncan Campbell if (!rank) {
299d2fd5932SHong 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,
300d2fd5932SHong Zhang network->cloneshared->Nsvtx));
301d2fd5932SHong Zhang }
302d2fd5932SHong Zhang
303d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL));
304d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPushSynchronized(viewer));
305d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv));
306d2fd5932SHong Zhang
307d2fd5932SHong Zhang for (i = 0; i < nsubnet; i++) {
308d2fd5932SHong Zhang PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges));
309d2fd5932SHong Zhang if (ne) {
310d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv));
311d2fd5932SHong Zhang for (j = 0; j < ne; j++) {
312d2fd5932SHong Zhang p = edges[j];
313d2fd5932SHong Zhang PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone));
314d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom));
315d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto));
316d2fd5932SHong Zhang PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p));
317d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto));
318d2fd5932SHong Zhang }
319d2fd5932SHong Zhang }
320d2fd5932SHong Zhang }
321d2fd5932SHong Zhang
322d2fd5932SHong Zhang /* Shared vertices */
323d2fd5932SHong Zhang PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx));
324d2fd5932SHong Zhang if (nsv) {
325d2fd5932SHong Zhang PetscInt gidx;
326d2fd5932SHong Zhang PetscBool ghost;
327d2fd5932SHong Zhang const PetscInt *sv = NULL;
328d2fd5932SHong Zhang
329d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n"));
330d2fd5932SHong Zhang for (i = 0; i < nsv; i++) {
331d2fd5932SHong Zhang PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost));
332d2fd5932SHong Zhang if (ghost) continue;
333d2fd5932SHong Zhang
334d2fd5932SHong Zhang PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv));
335d2fd5932SHong Zhang PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1]));
336d2fd5932SHong Zhang for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1]));
337d2fd5932SHong Zhang }
338d2fd5932SHong Zhang }
339d2fd5932SHong Zhang PetscCall(PetscViewerFlush(viewer));
340d2fd5932SHong Zhang PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3419f196a02SMartin Diehl } else PetscCheck(isascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
342d2fd5932SHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
343d2fd5932SHong Zhang }
3445f25b224SDuncan Campbell
3455f25b224SDuncan Campbell /*@
3465f25b224SDuncan Campbell DMNetworkViewSetShowRanks - Sets viewing the `DMETNWORK` on each rank individually.
3475f25b224SDuncan Campbell
3485f25b224SDuncan Campbell Logically Collective
3495f25b224SDuncan Campbell
3505f25b224SDuncan Campbell Input Parameter:
3515f25b224SDuncan Campbell . dm - the `DMNETWORK` object
3525f25b224SDuncan Campbell
3535f25b224SDuncan Campbell Output Parameter:
3545f25b224SDuncan Campbell . showranks - `PETSC_TRUE` if viewing each rank's sub network individually
3555f25b224SDuncan Campbell
3565f25b224SDuncan Campbell Level: beginner
3575f25b224SDuncan Campbell
3585f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
3595f25b224SDuncan Campbell @*/
DMNetworkViewSetShowRanks(DM dm,PetscBool showranks)3605f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowRanks(DM dm, PetscBool showranks)
3615f25b224SDuncan Campbell {
3625f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
3635f25b224SDuncan Campbell
3645f25b224SDuncan Campbell PetscFunctionBegin;
3655f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
3665f25b224SDuncan Campbell network->vieweroptions.showallranks = showranks;
3675f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
3685f25b224SDuncan Campbell }
3695f25b224SDuncan Campbell
3705f25b224SDuncan Campbell /*@
3715f25b224SDuncan Campbell DMNetworkViewSetShowGlobal - Set viewing the global network.
3725f25b224SDuncan Campbell
3735f25b224SDuncan Campbell Logically Collective
3745f25b224SDuncan Campbell
3755f25b224SDuncan Campbell Input Parameter:
3765f25b224SDuncan Campbell . dm - the `DMNETWORK` object
3775f25b224SDuncan Campbell
3785f25b224SDuncan Campbell Output Parameter:
3795f25b224SDuncan Campbell . showglobal - `PETSC_TRUE` if viewing the global network
3805f25b224SDuncan Campbell
3815f25b224SDuncan Campbell Level: beginner
3825f25b224SDuncan Campbell
3835f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
3845f25b224SDuncan Campbell @*/
DMNetworkViewSetShowGlobal(DM dm,PetscBool showglobal)3855f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowGlobal(DM dm, PetscBool showglobal)
3865f25b224SDuncan Campbell {
3875f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
3885f25b224SDuncan Campbell
3895f25b224SDuncan Campbell PetscFunctionBegin;
3905f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
3915f25b224SDuncan Campbell network->vieweroptions.dontshowglobal = (PetscBool)(!showglobal);
3925f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
3935f25b224SDuncan Campbell }
3945f25b224SDuncan Campbell
3955f25b224SDuncan Campbell /*@
3965f25b224SDuncan Campbell DMNetworkViewSetShowVertices - Sets whether to display the vertices in viewing routines.
3975f25b224SDuncan Campbell
3985f25b224SDuncan Campbell Logically Collective
3995f25b224SDuncan Campbell
4005f25b224SDuncan Campbell Input Parameter:
4015f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4025f25b224SDuncan Campbell
4035f25b224SDuncan Campbell Output Parameter:
4045f25b224SDuncan Campbell . showvertices - `PETSC_TRUE` if visualizing the vertices
4055f25b224SDuncan Campbell
4065f25b224SDuncan Campbell Level: beginner
4075f25b224SDuncan Campbell
4085f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowNumbering()`, `DMNetworkViewSetViewRanks()`
4095f25b224SDuncan Campbell @*/
DMNetworkViewSetShowVertices(DM dm,PetscBool showvertices)4105f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowVertices(DM dm, PetscBool showvertices)
4115f25b224SDuncan Campbell {
4125f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
4135f25b224SDuncan Campbell
4145f25b224SDuncan Campbell PetscFunctionBegin;
4155f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4165f25b224SDuncan Campbell network->vieweroptions.shownovertices = (PetscBool)(!showvertices);
4175f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
4185f25b224SDuncan Campbell }
4195f25b224SDuncan Campbell
4205f25b224SDuncan Campbell /*@
4215f25b224SDuncan Campbell DMNetworkViewSetShowNumbering - Set displaying the numbering of edges and vertices in viewing routines.
4225f25b224SDuncan Campbell
4235f25b224SDuncan Campbell Logically Collective
4245f25b224SDuncan Campbell
4255f25b224SDuncan Campbell Input Parameter:
4265f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4275f25b224SDuncan Campbell
4285f25b224SDuncan Campbell Output Parameter:
4295f25b224SDuncan Campbell . shownumbering - `PETSC_TRUE` if displaying the numbering of edges and vertices
4305f25b224SDuncan Campbell
4315f25b224SDuncan Campbell Level: beginner
4325f25b224SDuncan Campbell
4335f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetViewRanks()`
4345f25b224SDuncan Campbell @*/
DMNetworkViewSetShowNumbering(DM dm,PetscBool shownumbering)4355f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetShowNumbering(DM dm, PetscBool shownumbering)
4365f25b224SDuncan Campbell {
4375f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
4385f25b224SDuncan Campbell
4395f25b224SDuncan Campbell PetscFunctionBegin;
4405f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4415f25b224SDuncan Campbell network->vieweroptions.shownonumbering = (PetscBool)(!shownumbering);
4425f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
4435f25b224SDuncan Campbell }
4445f25b224SDuncan Campbell
4455f25b224SDuncan Campbell /*@
4465f25b224SDuncan Campbell DMNetworkViewSetViewRanks - View the `DMNETWORK` on each of the specified ranks individually.
4475f25b224SDuncan Campbell
4485f25b224SDuncan Campbell Collective
4495f25b224SDuncan Campbell
4505f25b224SDuncan Campbell Input Parameter:
4515f25b224SDuncan Campbell . dm - the `DMNETWORK` object
4525f25b224SDuncan Campbell
4535f25b224SDuncan Campbell Output Parameter:
4545f25b224SDuncan Campbell . viewranks - set of ranks to view the `DMNETWORK` on individually
4555f25b224SDuncan Campbell
4565f25b224SDuncan Campbell Level: beginner
4575f25b224SDuncan Campbell
4585f25b224SDuncan Campbell Note:
4595f25b224SDuncan Campbell `DMNetwork` takes ownership of the input viewranks `IS`, it should be destroyed by the caller.
4605f25b224SDuncan Campbell
4615f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkViewSetShowRanks()`, `DMNetworkViewSetShowGlobal()`, `DMNetworkViewSetShowVertices()`, `DMNetworkViewSetShowNumbering()`
4625f25b224SDuncan Campbell @*/
DMNetworkViewSetViewRanks(DM dm,IS viewranks)4635f25b224SDuncan Campbell PetscErrorCode DMNetworkViewSetViewRanks(DM dm, IS viewranks)
4645f25b224SDuncan Campbell {
4655f25b224SDuncan Campbell DM_Network *network = (DM_Network *)dm->data;
4665f25b224SDuncan Campbell
4675f25b224SDuncan Campbell PetscFunctionBegin;
4685f25b224SDuncan Campbell PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
4695f25b224SDuncan Campbell PetscValidHeaderSpecific(viewranks, IS_CLASSID, 2);
4705f25b224SDuncan Campbell PetscCheckSameComm(dm, 1, viewranks, 2);
4715f25b224SDuncan Campbell PetscCall(ISDestroy(&network->vieweroptions.viewranks));
4725f25b224SDuncan Campbell PetscCall(PetscObjectReference((PetscObject)viewranks));
4735f25b224SDuncan Campbell network->vieweroptions.viewranks = viewranks;
4745f25b224SDuncan Campbell PetscFunctionReturn(PETSC_SUCCESS);
4755f25b224SDuncan Campbell }
476