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 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; 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsSetValue(NULL,"-petscpartitioner_use_vertex_weights","No")); 25f11a936eSBarry Smith comm = PETSC_COMM_WORLD; 26*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 27*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 28f11a936eSBarry Smith 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkCreate(PETSC_COMM_WORLD,&network)); 30f11a936eSBarry Smith 31eac198afSGetnet /* Register zero size componets to get compkeys to be used by DMNetworkAddComponent() */ 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(network,"ecomp",0,&ecompkey)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(network,"vcomp",0,&vcompkey)); 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 425c6496baSHong Zhang if (rank == 0) { 435c6496baSHong Zhang nedge = 1; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(2*nedge,&edgelist)); 455c6496baSHong Zhang edgelist[0] = nodeOffset + 2; 465c6496baSHong Zhang edgelist[1] = nodeOffset + 3; 475c6496baSHong Zhang } else { 485c6496baSHong Zhang nedge = 2; 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(2*nedge,&edgelist)); 50f11a936eSBarry Smith edgelist[0] = nodeOffset + 0; 51f11a936eSBarry Smith edgelist[1] = nodeOffset + 2; 52f11a936eSBarry Smith edgelist[2] = nodeOffset + 1; 53f11a936eSBarry Smith edgelist[3] = nodeOffset + 2; 545c6496baSHong Zhang } 55f11a936eSBarry Smith 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkSetNumSubNetworks(network,PETSC_DECIDE,1)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddSubnetwork(network,"Subnetwork 1",nedge,edgelist,NULL)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkLayoutSetUp(network)); 59f11a936eSBarry Smith 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Network after DMNetworkLayoutSetUp:\n")); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(network,PETSC_VIEWER_STDOUT_WORLD)); 625c6496baSHong Zhang 63f11a936eSBarry Smith /* Add components and variables for the network */ 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(network,0,&nv,&ne,&nodes,&edges)); 65f11a936eSBarry Smith for (e = 0; e < ne; e++) { 66f11a936eSBarry Smith /* The edges have no degrees of freedom */ 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(network,edges[e],ecompkey,NULL,1)); 68f11a936eSBarry Smith } 69f11a936eSBarry Smith for (v = 0; v < nv; v++) { 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(network,nodes[v],vcompkey,NULL,2)); 71f11a936eSBarry Smith } 72f11a936eSBarry Smith 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(network)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetPlex(network,&plex)); 75*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */ 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(plex,§ion)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD)); 78f11a936eSBarry Smith 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(edgelist)); 80f11a936eSBarry Smith 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkDistribute(&network,0)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nNetwork after DMNetworkDistribute:\n")); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(network,PETSC_VIEWER_STDOUT_WORLD)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetPlex(network,&plex)); 85*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */ 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(plex,§ion)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD)); 88f11a936eSBarry Smith 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&network)); 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 1015c6496baSHong Zhang args: -petscpartitioner_type simple 102f11a936eSBarry Smith 103f11a936eSBarry Smith TEST*/ 104