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