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 9*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 10*d71ae5a4SJacob Faibussowitsch { 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 22327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-petscpartitioner_use_vertex_weights", "No")); 25f11a936eSBarry Smith comm = PETSC_COMM_WORLD; 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 28f11a936eSBarry Smith 299566063dSJacob Faibussowitsch PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &network)); 30f11a936eSBarry Smith 316aad120cSJose E. Roman /* Register zero size components to get compkeys to be used by DMNetworkAddComponent() */ 329566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(network, "ecomp", 0, &ecompkey)); 339566063dSJacob Faibussowitsch PetscCall(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; 449566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * nedge, &edgelist)); 455c6496baSHong Zhang edgelist[0] = nodeOffset + 2; 465c6496baSHong Zhang edgelist[1] = nodeOffset + 3; 475c6496baSHong Zhang } else { 485c6496baSHong Zhang nedge = 2; 499566063dSJacob Faibussowitsch PetscCall(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 569566063dSJacob Faibussowitsch PetscCall(DMNetworkSetNumSubNetworks(network, PETSC_DECIDE, 1)); 579566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(network, "Subnetwork 1", nedge, edgelist, NULL)); 589566063dSJacob Faibussowitsch PetscCall(DMNetworkLayoutSetUp(network)); 59f11a936eSBarry Smith 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Network after DMNetworkLayoutSetUp:\n")); 619566063dSJacob Faibussowitsch PetscCall(DMView(network, PETSC_VIEWER_STDOUT_WORLD)); 625c6496baSHong Zhang 63f11a936eSBarry Smith /* Add components and variables for the network */ 649566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(network, 0, &nv, &ne, &nodes, &edges)); 65f11a936eSBarry Smith for (e = 0; e < ne; e++) { 66f11a936eSBarry Smith /* The edges have no degrees of freedom */ 679566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(network, edges[e], ecompkey, NULL, 1)); 68f11a936eSBarry Smith } 6948a46eb9SPierre Jolivet for (v = 0; v < nv; v++) PetscCall(DMNetworkAddComponent(network, nodes[v], vcompkey, NULL, 2)); 70f11a936eSBarry Smith 719566063dSJacob Faibussowitsch PetscCall(DMSetUp(network)); 729566063dSJacob Faibussowitsch PetscCall(DMNetworkGetPlex(network, &plex)); 739566063dSJacob Faibussowitsch /* PetscCall(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */ 749566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(plex, §ion)); 759566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD)); 76f11a936eSBarry Smith 779566063dSJacob Faibussowitsch PetscCall(PetscFree(edgelist)); 78f11a936eSBarry Smith 799566063dSJacob Faibussowitsch PetscCall(DMNetworkDistribute(&network, 0)); 809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nNetwork after DMNetworkDistribute:\n")); 819566063dSJacob Faibussowitsch PetscCall(DMView(network, PETSC_VIEWER_STDOUT_WORLD)); 829566063dSJacob Faibussowitsch PetscCall(DMNetworkGetPlex(network, &plex)); 839566063dSJacob Faibussowitsch /* PetscCall(DMView(plex,PETSC_VIEWER_STDOUT_WORLD)); */ 849566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(plex, §ion)); 859566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD)); 86f11a936eSBarry Smith 879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&network)); 889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 89b122ec5aSJacob Faibussowitsch return 0; 90f11a936eSBarry Smith } 91f11a936eSBarry Smith 92f11a936eSBarry Smith /*TEST 93f11a936eSBarry Smith 94f11a936eSBarry Smith build: 95f11a936eSBarry Smith requires: !complex double 96f11a936eSBarry Smith 97f11a936eSBarry Smith test: 98f11a936eSBarry Smith nsize: 2 995c6496baSHong Zhang args: -petscpartitioner_type simple 100f11a936eSBarry Smith 101f11a936eSBarry Smith TEST*/ 102