#include /*I "petscdmnetwork.h" I*/ PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) { PetscBool iascii; PetscMPIInt rank; PetscFunctionBegin; PetscValidHeaderSpecific(dm, DM_CLASSID, 1); PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); if (iascii) { const PetscInt *cone, *vtx, *edges; PetscInt vfrom, vto, i, j, nv, ne, nsv, p, nsubnet; DM_Network *network = (DM_Network *)dm->data; nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */ if (rank == 0) { 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, network->cloneshared->Nsvtx)); } PetscCall(DMNetworkGetSharedVertices(dm, &nsv, NULL)); PetscCall(PetscViewerASCIIPushSynchronized(viewer)); PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv)); for (i = 0; i < nsubnet; i++) { PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges)); if (ne) { PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv)); for (j = 0; j < ne; j++) { p = edges[j]; PetscCall(DMNetworkGetConnectedVertices(dm, p, &cone)); PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom)); PetscCall(DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto)); PetscCall(DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p)); PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto)); } } } /* Shared vertices */ PetscCall(DMNetworkGetSharedVertices(dm, NULL, &vtx)); if (nsv) { PetscInt gidx; PetscBool ghost; const PetscInt *sv = NULL; PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n")); for (i = 0; i < nsv; i++) { PetscCall(DMNetworkIsGhostVertex(dm, vtx[i], &ghost)); if (ghost) continue; PetscCall(DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv)); PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1])); for (j = 1; j < nv; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1])); } } PetscCall(PetscViewerFlush(viewer)); PetscCall(PetscViewerASCIIPopSynchronized(viewer)); } else PetscCheck(iascii, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name); PetscFunctionReturn(PETSC_SUCCESS); }