Actual source code: ex10.c
1: /*
2: Simple example demonstrating creating a one sub-network DMNetwork in parallel.
4: In this example vertices 0 and 1 are not connected to any edges.
5: */
7: #include <petscdmnetwork.h>
9: int main(int argc, char **argv)
10: {
11: DM network;
12: PetscMPIInt size, rank;
13: MPI_Comm comm;
14: PetscInt e, ne, nv, v, ecompkey, vcompkey;
15: PetscInt *edgelist = NULL;
16: const PetscInt *nodes, *edges;
17: DM plex;
18: PetscSection section;
19: PetscInt Ne, Ni;
20: PetscInt nodeOffset, k = 2, nedge;
23: PetscInitialize(&argc, &argv, NULL, NULL);
24: PetscOptionsSetValue(NULL, "-petscpartitioner_use_vertex_weights", "No");
25: comm = PETSC_COMM_WORLD;
26: MPI_Comm_rank(comm, &rank);
27: MPI_Comm_size(comm, &size);
29: DMNetworkCreate(PETSC_COMM_WORLD, &network);
31: /* Register zero size components to get compkeys to be used by DMNetworkAddComponent() */
32: DMNetworkRegisterComponent(network, "ecomp", 0, &ecompkey);
33: DMNetworkRegisterComponent(network, "vcomp", 0, &vcompkey);
35: Ne = 2;
36: Ni = 1;
37: nodeOffset = (Ne + Ni) * rank; /* The global node index of the first node defined on this process */
39: /* There are three nodes on each rank and two edges. The edges only connect nodes on the given rank */
40: nedge = k * Ni;
42: if (rank == 0) {
43: nedge = 1;
44: PetscCalloc1(2 * nedge, &edgelist);
45: edgelist[0] = nodeOffset + 2;
46: edgelist[1] = nodeOffset + 3;
47: } else {
48: nedge = 2;
49: PetscCalloc1(2 * nedge, &edgelist);
50: edgelist[0] = nodeOffset + 0;
51: edgelist[1] = nodeOffset + 2;
52: edgelist[2] = nodeOffset + 1;
53: edgelist[3] = nodeOffset + 2;
54: }
56: DMNetworkSetNumSubNetworks(network, PETSC_DECIDE, 1);
57: DMNetworkAddSubnetwork(network, "Subnetwork 1", nedge, edgelist, NULL);
58: DMNetworkLayoutSetUp(network);
60: PetscPrintf(PETSC_COMM_WORLD, "Network after DMNetworkLayoutSetUp:\n");
61: DMView(network, PETSC_VIEWER_STDOUT_WORLD);
63: /* Add components and variables for the network */
64: DMNetworkGetSubnetwork(network, 0, &nv, &ne, &nodes, &edges);
65: for (e = 0; e < ne; e++) {
66: /* The edges have no degrees of freedom */
67: DMNetworkAddComponent(network, edges[e], ecompkey, NULL, 1);
68: }
69: for (v = 0; v < nv; v++) DMNetworkAddComponent(network, nodes[v], vcompkey, NULL, 2);
71: DMSetUp(network);
72: DMNetworkGetPlex(network, &plex);
73: /* DMView(plex,PETSC_VIEWER_STDOUT_WORLD); */
74: DMGetLocalSection(plex, §ion);
75: PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);
77: PetscFree(edgelist);
79: DMNetworkDistribute(&network, 0);
80: PetscPrintf(PETSC_COMM_WORLD, "\nNetwork after DMNetworkDistribute:\n");
81: DMView(network, PETSC_VIEWER_STDOUT_WORLD);
82: DMNetworkGetPlex(network, &plex);
83: /* DMView(plex,PETSC_VIEWER_STDOUT_WORLD); */
84: DMGetLocalSection(plex, §ion);
85: PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);
87: DMDestroy(&network);
88: PetscFinalize();
89: return 0;
90: }
92: /*TEST
94: build:
95: requires: !complex double
97: test:
98: nsize: 2
99: args: -petscpartitioner_type simple
101: TEST*/