xref: /petsc/src/mat/tutorials/ex15.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
21*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &args, (char*) 0, help));
22c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL));
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm, &A));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(A, 3, NULL));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_vertex_weights",&set_vweights,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_use_edge_weights",&use_edge_weights,NULL));
32c4762a1bSJed Brown   /* Create a linear mesh */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A, &start, &end));
34c4762a1bSJed Brown   if (set_vweights) {
355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(end-start,&vweights));
36c4762a1bSJed Brown     for (r = start; r < end; ++r)
37c4762a1bSJed Brown       vweights[r-start] = rank+1;
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   for (r = start; r < end; ++r) {
40c4762a1bSJed Brown     if (r == 0) {
41c4762a1bSJed Brown       PetscInt    cols[2];
42c4762a1bSJed Brown       PetscScalar vals[2];
43c4762a1bSJed Brown 
44c4762a1bSJed Brown       cols[0] = r;   cols[1] = r+1;
45c4762a1bSJed Brown       vals[0] = 1.0; vals[1] = use_edge_weights? 2.0: 1.0;
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
48c4762a1bSJed Brown     } else if (r == N-1) {
49c4762a1bSJed Brown       PetscInt    cols[2];
50c4762a1bSJed Brown       PetscScalar vals[2];
51c4762a1bSJed Brown 
52c4762a1bSJed Brown       cols[0] = r-1; cols[1] = r;
53c4762a1bSJed Brown       vals[0] = use_edge_weights? 3.0:1.0; vals[1] = 1.0;
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
56c4762a1bSJed Brown     } else {
57c4762a1bSJed Brown       PetscInt    cols[3];
58c4762a1bSJed Brown       PetscScalar vals[3];
59c4762a1bSJed Brown 
60c4762a1bSJed Brown       cols[0] = r-1; cols[1] = r;   cols[2] = r+1;
61c4762a1bSJed Brown       /* ADJ matrix needs to be symmetric */
62c4762a1bSJed Brown       vals[0] = use_edge_weights? (cols[0]==0? 2.0:5.0):1.0;
63c4762a1bSJed Brown       vals[1] = 1.0;
64c4762a1bSJed Brown       vals[2] = use_edge_weights? (cols[2]==N-1? 3.0:5.0):1.0;
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES));
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown   }
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningCreate(comm, &part));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetAdjacency(part, A));
74c4762a1bSJed Brown   if (set_vweights) {
755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningSetVertexWeights(part,vweights));
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown   if (use_edge_weights) {
785f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningSetUseEdgeWeights(part,use_edge_weights));
79c4762a1bSJed Brown 
805f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningGetUseEdgeWeights(part,&use_edge_weights));
8128b400f6SJacob Faibussowitsch     PetscCheck(use_edge_weights,comm,PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly ");
82c4762a1bSJed Brown   }
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningSetFromOptions(part));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningApply(part, &is));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPartitioningDestroy(&part));
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
90*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
91*b122ec5aSJacob Faibussowitsch   return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*TEST
95c4762a1bSJed Brown 
96c4762a1bSJed Brown    test:
97c4762a1bSJed Brown       nsize: 3
98c4762a1bSJed Brown       requires: parmetis
99c4762a1bSJed Brown       args: -mat_partitioning_type parmetis
100c4762a1bSJed Brown 
101c4762a1bSJed Brown    test:
102c4762a1bSJed Brown       suffix: 2
103c4762a1bSJed Brown       nsize: 3
104c4762a1bSJed Brown       requires: ptscotch
105c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch
106c4762a1bSJed Brown 
107c4762a1bSJed Brown    test:
108c4762a1bSJed Brown       suffix: 3
109c4762a1bSJed Brown       nsize: 4
110c4762a1bSJed Brown       requires: party
111c4762a1bSJed Brown       args: -mat_partitioning_type party
112c4762a1bSJed Brown 
113c4762a1bSJed Brown    test:
114c4762a1bSJed Brown       suffix: 4
115c4762a1bSJed Brown       nsize: 3
116c4762a1bSJed Brown       requires: chaco
117c4762a1bSJed Brown       args: -mat_partitioning_type chaco
118c4762a1bSJed Brown 
119c4762a1bSJed Brown    test:
120c4762a1bSJed Brown       suffix: 5
121c4762a1bSJed Brown       nsize: 3
122c4762a1bSJed Brown       requires: parmetis
123c4762a1bSJed Brown       args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100
124c4762a1bSJed Brown 
125c4762a1bSJed Brown    test:
126c4762a1bSJed Brown       suffix: 6
127c4762a1bSJed Brown       nsize: 3
128c4762a1bSJed Brown       requires: parmetis
129c4762a1bSJed 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
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    test:
132c4762a1bSJed Brown       suffix: 7
133c4762a1bSJed Brown       nsize: 2
134c4762a1bSJed Brown       requires: parmetis
135c4762a1bSJed 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
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: 8
139c4762a1bSJed Brown       nsize: 2
140c4762a1bSJed Brown       requires: parmetis
141c4762a1bSJed Brown       args: -mat_partitioning_type parmetis -mat_partitioning_nparts 3 -test_use_edge_weights 1
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       suffix: 9
145c4762a1bSJed Brown       nsize: 2
146c4762a1bSJed Brown       requires: ptscotch
147c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch -mat_partitioning_nparts 3 -test_use_edge_weights 1 -mat_partitioning_ptscotch_proc_weight 0
148c4762a1bSJed Brown 
149c4762a1bSJed Brown TEST*/
150