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