1c4762a1bSJed Brown static char help[] = "Example of using graph partitioning to partition a graph\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*T 4c4762a1bSJed Brown Concepts: Mat^mat partitioning 5c4762a1bSJed Brown Concepts: Mat^image segmentation 6c4762a1bSJed Brown Processors: n 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc, char **args) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown Mat A; 14c4762a1bSJed Brown MatPartitioning part; 15c4762a1bSJed Brown IS is; 16c4762a1bSJed Brown PetscInt r,N = 10, start, end, *vweights; 17c4762a1bSJed Brown PetscBool set_vweights=PETSC_FALSE,use_edge_weights=PETSC_FALSE; 18c4762a1bSJed Brown PetscMPIInt rank; 19c4762a1bSJed Brown MPI_Comm comm; 20c4762a1bSJed Brown PetscErrorCode ierr; 21c4762a1bSJed Brown 22c4762a1bSJed Brown ierr = PetscInitialize(&argc, &args, (char*) 0, help);if (ierr) return ierr; 23c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL)); 255f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm, &A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A, 3, NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_vertex_weights",&set_vweights,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_use_edge_weights",&use_edge_weights,NULL)); 33c4762a1bSJed Brown /* Create a linear mesh */ 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A, &start, &end)); 35c4762a1bSJed Brown if (set_vweights) { 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(end-start,&vweights)); 37c4762a1bSJed Brown for (r = start; r < end; ++r) 38c4762a1bSJed Brown vweights[r-start] = rank+1; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown for (r = start; r < end; ++r) { 41c4762a1bSJed Brown if (r == 0) { 42c4762a1bSJed Brown PetscInt cols[2]; 43c4762a1bSJed Brown PetscScalar vals[2]; 44c4762a1bSJed Brown 45c4762a1bSJed Brown cols[0] = r; cols[1] = r+1; 46c4762a1bSJed Brown vals[0] = 1.0; vals[1] = use_edge_weights? 2.0: 1.0; 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES)); 49c4762a1bSJed Brown } else if (r == N-1) { 50c4762a1bSJed Brown PetscInt cols[2]; 51c4762a1bSJed Brown PetscScalar vals[2]; 52c4762a1bSJed Brown 53c4762a1bSJed Brown cols[0] = r-1; cols[1] = r; 54c4762a1bSJed Brown vals[0] = use_edge_weights? 3.0:1.0; vals[1] = 1.0; 55c4762a1bSJed Brown 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES)); 57c4762a1bSJed Brown } else { 58c4762a1bSJed Brown PetscInt cols[3]; 59c4762a1bSJed Brown PetscScalar vals[3]; 60c4762a1bSJed Brown 61c4762a1bSJed Brown cols[0] = r-1; cols[1] = r; cols[2] = r+1; 62c4762a1bSJed Brown /* ADJ matrix needs to be symmetric */ 63c4762a1bSJed Brown vals[0] = use_edge_weights? (cols[0]==0? 2.0:5.0):1.0; 64c4762a1bSJed Brown vals[1] = 1.0; 65c4762a1bSJed Brown vals[2] = use_edge_weights? (cols[2]==N-1? 3.0:5.0):1.0; 66c4762a1bSJed Brown 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES)); 68c4762a1bSJed Brown } 69c4762a1bSJed Brown } 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningCreate(comm, &part)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetAdjacency(part, A)); 75c4762a1bSJed Brown if (set_vweights) { 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetVertexWeights(part,vweights)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown if (use_edge_weights) { 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetUseEdgeWeights(part,use_edge_weights)); 80c4762a1bSJed Brown 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningGetUseEdgeWeights(part,&use_edge_weights)); 82*28b400f6SJacob Faibussowitsch PetscCheck(use_edge_weights,comm,PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly "); 83c4762a1bSJed Brown } 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetFromOptions(part)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningApply(part, &is)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(is, PETSC_VIEWER_STDOUT_WORLD)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningDestroy(&part)); 89c4762a1bSJed Brown 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 91c4762a1bSJed Brown ierr = PetscFinalize(); 92c4762a1bSJed Brown return ierr; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /*TEST 96c4762a1bSJed Brown 97c4762a1bSJed Brown test: 98c4762a1bSJed Brown nsize: 3 99c4762a1bSJed Brown requires: parmetis 100c4762a1bSJed Brown args: -mat_partitioning_type parmetis 101c4762a1bSJed Brown 102c4762a1bSJed Brown test: 103c4762a1bSJed Brown suffix: 2 104c4762a1bSJed Brown nsize: 3 105c4762a1bSJed Brown requires: ptscotch 106c4762a1bSJed Brown args: -mat_partitioning_type ptscotch 107c4762a1bSJed Brown 108c4762a1bSJed Brown test: 109c4762a1bSJed Brown suffix: 3 110c4762a1bSJed Brown nsize: 4 111c4762a1bSJed Brown requires: party 112c4762a1bSJed Brown args: -mat_partitioning_type party 113c4762a1bSJed Brown 114c4762a1bSJed Brown test: 115c4762a1bSJed Brown suffix: 4 116c4762a1bSJed Brown nsize: 3 117c4762a1bSJed Brown requires: chaco 118c4762a1bSJed Brown args: -mat_partitioning_type chaco 119c4762a1bSJed Brown 120c4762a1bSJed Brown test: 121c4762a1bSJed Brown suffix: 5 122c4762a1bSJed Brown nsize: 3 123c4762a1bSJed Brown requires: parmetis 124c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown suffix: 6 128c4762a1bSJed Brown nsize: 3 129c4762a1bSJed Brown requires: parmetis 130c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100 -test_vertex_weights 1 -mat_partitioning_use_edge_weights 1 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown suffix: 7 134c4762a1bSJed Brown nsize: 2 135c4762a1bSJed Brown requires: parmetis 136c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 2 -mat_partitioning_nparts 10 -mat_partitioning_hierarchical_fineparttype hierarch -malloc_dump -N 100 -mat_partitioning_improve 1 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: 8 140c4762a1bSJed Brown nsize: 2 141c4762a1bSJed Brown requires: parmetis 142c4762a1bSJed Brown args: -mat_partitioning_type parmetis -mat_partitioning_nparts 3 -test_use_edge_weights 1 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown suffix: 9 146c4762a1bSJed Brown nsize: 2 147c4762a1bSJed Brown requires: ptscotch 148c4762a1bSJed Brown args: -mat_partitioning_type ptscotch -mat_partitioning_nparts 3 -test_use_edge_weights 1 -mat_partitioning_ptscotch_proc_weight 0 149c4762a1bSJed Brown 150c4762a1bSJed Brown TEST*/ 151