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