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