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; 24c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL);CHKERRQ(ierr); 25ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 26c4762a1bSJed Brown ierr = MatCreate(comm, &A);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A, 3, NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_vertex_weights",&set_vweights,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_use_edge_weights",&use_edge_weights,NULL);CHKERRQ(ierr); 33c4762a1bSJed Brown /* Create a linear mesh */ 34c4762a1bSJed Brown ierr = MatGetOwnershipRange(A, &start, &end);CHKERRQ(ierr); 35c4762a1bSJed Brown if (set_vweights) { 36c4762a1bSJed Brown ierr = PetscMalloc1(end-start,&vweights);CHKERRQ(ierr); 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 48c4762a1bSJed Brown ierr = MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr); 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 56c4762a1bSJed Brown ierr = MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr); 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 67c4762a1bSJed Brown ierr = MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES);CHKERRQ(ierr); 68c4762a1bSJed Brown } 69c4762a1bSJed Brown } 70c4762a1bSJed Brown ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = MatPartitioningCreate(comm, &part);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatPartitioningSetAdjacency(part, A);CHKERRQ(ierr); 75c4762a1bSJed Brown if (set_vweights) { 76c4762a1bSJed Brown ierr = MatPartitioningSetVertexWeights(part,vweights);CHKERRQ(ierr); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown if (use_edge_weights) { 79c4762a1bSJed Brown ierr = MatPartitioningSetUseEdgeWeights(part,use_edge_weights);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown ierr = MatPartitioningGetUseEdgeWeights(part,&use_edge_weights);CHKERRQ(ierr); 82*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!use_edge_weights,comm,PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly "); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = MatPartitioningApply(part, &is);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = ISView(is, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 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