xref: /petsc/src/dm/tests/ex10.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1f11a936eSBarry Smith /*
25c6496baSHong Zhang     Simple example demonstrating creating a one sub-network DMNetwork in parallel.
35c6496baSHong Zhang 
45c6496baSHong 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   DM                network;
12f11a936eSBarry Smith   PetscMPIInt       size,rank;
13f11a936eSBarry Smith   MPI_Comm          comm;
14eac198afSGetnet   PetscInt          e,ne,nv,v,ecompkey,vcompkey;
15f11a936eSBarry Smith   PetscInt          *edgelist = NULL;
16f11a936eSBarry Smith   const PetscInt    *nodes,*edges;
17f11a936eSBarry Smith   DM                plex;
18f11a936eSBarry Smith   PetscSection      section;
19f11a936eSBarry Smith   PetscInt          Ne,Ni;
20f11a936eSBarry Smith   PetscInt          nodeOffset,k = 2,nedge;
21f11a936eSBarry Smith 
22*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,NULL));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsSetValue(NULL,"-petscpartitioner_use_vertex_weights","No"));
24f11a936eSBarry Smith   comm = PETSC_COMM_WORLD;
255f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
265f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
27f11a936eSBarry Smith 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreate(PETSC_COMM_WORLD,&network));
29f11a936eSBarry Smith 
30eac198afSGetnet   /* Register zero size componets to get compkeys to be used by DMNetworkAddComponent() */
315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkRegisterComponent(network,"ecomp",0,&ecompkey));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkRegisterComponent(network,"vcomp",0,&vcompkey));
33eac198afSGetnet 
34f11a936eSBarry Smith   Ne = 2;
35f11a936eSBarry Smith   Ni = 1;
36f11a936eSBarry Smith   nodeOffset = (Ne+Ni)*rank;   /* The global node index of the first node defined on this process */
37f11a936eSBarry Smith 
38f11a936eSBarry Smith   /* There are three nodes on each rank and two edges. The edges only connect nodes on the given rank */
39f11a936eSBarry Smith   nedge = k * Ni;
40f11a936eSBarry Smith 
415c6496baSHong Zhang   if (rank == 0) {
425c6496baSHong Zhang     nedge = 1;
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(2*nedge,&edgelist));
445c6496baSHong Zhang     edgelist[0] = nodeOffset + 2;
455c6496baSHong Zhang     edgelist[1] = nodeOffset + 3;
465c6496baSHong Zhang   } else {
475c6496baSHong Zhang     nedge = 2;
485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(2*nedge,&edgelist));
49f11a936eSBarry Smith     edgelist[0] = nodeOffset + 0;
50f11a936eSBarry Smith     edgelist[1] = nodeOffset + 2;
51f11a936eSBarry Smith     edgelist[2] = nodeOffset + 1;
52f11a936eSBarry Smith     edgelist[3] = nodeOffset + 2;
535c6496baSHong Zhang   }
54f11a936eSBarry Smith 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkSetNumSubNetworks(network,PETSC_DECIDE,1));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkAddSubnetwork(network,"Subnetwork 1",nedge,edgelist,NULL));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkLayoutSetUp(network));
58f11a936eSBarry Smith 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Network after DMNetworkLayoutSetUp:\n"));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(network,PETSC_VIEWER_STDOUT_WORLD));
615c6496baSHong Zhang 
62f11a936eSBarry Smith   /* Add components and variables for the network */
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetSubnetwork(network,0,&nv,&ne,&nodes,&edges));
64f11a936eSBarry Smith   for (e = 0; e < ne; e++) {
65f11a936eSBarry Smith     /* The edges have no degrees of freedom */
665f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkAddComponent(network,edges[e],ecompkey,NULL,1));
67f11a936eSBarry Smith   }
68f11a936eSBarry Smith   for (v = 0; v < nv; v++) {
695f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkAddComponent(network,nodes[v],vcompkey,NULL,2));
70f11a936eSBarry Smith   }
71f11a936eSBarry Smith 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(network));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetPlex(network,&plex));
745f80ce2aSJacob Faibussowitsch   /* CHKERRQ(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */
755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(plex,&section));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD));
77f11a936eSBarry Smith 
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(edgelist));
79f11a936eSBarry Smith 
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkDistribute(&network,0));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nNetwork after DMNetworkDistribute:\n"));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(network,PETSC_VIEWER_STDOUT_WORLD));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetPlex(network,&plex));
845f80ce2aSJacob Faibussowitsch   /* CHKERRQ(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(plex,&section));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD));
87f11a936eSBarry Smith 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&network));
89*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
90*b122ec5aSJacob Faibussowitsch   return 0;
91f11a936eSBarry Smith }
92f11a936eSBarry Smith 
93f11a936eSBarry Smith /*TEST
94f11a936eSBarry Smith 
95f11a936eSBarry Smith    build:
96f11a936eSBarry Smith       requires: !complex double
97f11a936eSBarry Smith 
98f11a936eSBarry Smith    test:
99f11a936eSBarry Smith       nsize: 2
1005c6496baSHong Zhang       args: -petscpartitioner_type simple
101f11a936eSBarry Smith 
102f11a936eSBarry Smith TEST*/
103