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