xref: /petsc/src/mat/tutorials/ex15.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Example of using graph partitioning to partition a graph\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5*9371c9d4SSatish Balay int main(int argc, char **args) {
6c4762a1bSJed Brown   Mat             A;
7c4762a1bSJed Brown   MatPartitioning part;
8c4762a1bSJed Brown   IS              is;
9c4762a1bSJed Brown   PetscInt        r, N = 10, start, end, *vweights;
10c4762a1bSJed Brown   PetscBool       set_vweights = PETSC_FALSE, use_edge_weights = PETSC_FALSE;
11c4762a1bSJed Brown   PetscMPIInt     rank;
12c4762a1bSJed Brown   MPI_Comm        comm;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
16c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
199566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
209566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N));
219566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
229566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, 3, NULL));
239566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_vertex_weights", &set_vweights, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_use_edge_weights", &use_edge_weights, NULL));
26c4762a1bSJed Brown   /* Create a linear mesh */
279566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &start, &end));
28c4762a1bSJed Brown   if (set_vweights) {
299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(end - start, &vweights));
30*9371c9d4SSatish Balay     for (r = start; r < end; ++r) vweights[r - start] = rank + 1;
31c4762a1bSJed Brown   }
32c4762a1bSJed Brown   for (r = start; r < end; ++r) {
33c4762a1bSJed Brown     if (r == 0) {
34c4762a1bSJed Brown       PetscInt    cols[2];
35c4762a1bSJed Brown       PetscScalar vals[2];
36c4762a1bSJed Brown 
37*9371c9d4SSatish Balay       cols[0] = r;
38*9371c9d4SSatish Balay       cols[1] = r + 1;
39*9371c9d4SSatish Balay       vals[0] = 1.0;
40*9371c9d4SSatish Balay       vals[1] = use_edge_weights ? 2.0 : 1.0;
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
43c4762a1bSJed Brown     } else if (r == N - 1) {
44c4762a1bSJed Brown       PetscInt    cols[2];
45c4762a1bSJed Brown       PetscScalar vals[2];
46c4762a1bSJed Brown 
47*9371c9d4SSatish Balay       cols[0] = r - 1;
48*9371c9d4SSatish Balay       cols[1] = r;
49*9371c9d4SSatish Balay       vals[0] = use_edge_weights ? 3.0 : 1.0;
50*9371c9d4SSatish Balay       vals[1] = 1.0;
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
53c4762a1bSJed Brown     } else {
54c4762a1bSJed Brown       PetscInt    cols[3];
55c4762a1bSJed Brown       PetscScalar vals[3];
56c4762a1bSJed Brown 
57*9371c9d4SSatish Balay       cols[0] = r - 1;
58*9371c9d4SSatish Balay       cols[1] = r;
59*9371c9d4SSatish Balay       cols[2] = r + 1;
60c4762a1bSJed Brown       /* ADJ matrix needs to be symmetric */
61c4762a1bSJed Brown       vals[0] = use_edge_weights ? (cols[0] == 0 ? 2.0 : 5.0) : 1.0;
62c4762a1bSJed Brown       vals[1] = 1.0;
63c4762a1bSJed Brown       vals[2] = use_edge_weights ? (cols[2] == N - 1 ? 3.0 : 5.0) : 1.0;
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES));
66c4762a1bSJed Brown     }
67c4762a1bSJed Brown   }
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm, &part));
729566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, A));
731baa6e33SBarry Smith   if (set_vweights) PetscCall(MatPartitioningSetVertexWeights(part, vweights));
74c4762a1bSJed Brown   if (use_edge_weights) {
759566063dSJacob Faibussowitsch     PetscCall(MatPartitioningSetUseEdgeWeights(part, use_edge_weights));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch     PetscCall(MatPartitioningGetUseEdgeWeights(part, &use_edge_weights));
7828b400f6SJacob Faibussowitsch     PetscCheck(use_edge_weights, comm, PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly ");
79c4762a1bSJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
819566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &is));
829566063dSJacob Faibussowitsch   PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
849566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
88b122ec5aSJacob Faibussowitsch   return 0;
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /*TEST
92c4762a1bSJed Brown 
93c4762a1bSJed Brown    test:
94c4762a1bSJed Brown       nsize: 3
95c4762a1bSJed Brown       requires: parmetis
96c4762a1bSJed Brown       args: -mat_partitioning_type parmetis
97c4762a1bSJed Brown 
98c4762a1bSJed Brown    test:
99c4762a1bSJed Brown       suffix: 2
100c4762a1bSJed Brown       nsize: 3
101c4762a1bSJed Brown       requires: ptscotch
102c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch
103c4762a1bSJed Brown 
104c4762a1bSJed Brown    test:
105c4762a1bSJed Brown       suffix: 3
106c4762a1bSJed Brown       nsize: 4
107c4762a1bSJed Brown       requires: party
108c4762a1bSJed Brown       args: -mat_partitioning_type party
109c4762a1bSJed Brown 
110c4762a1bSJed Brown    test:
111c4762a1bSJed Brown       suffix: 4
112c4762a1bSJed Brown       nsize: 3
113c4762a1bSJed Brown       requires: chaco
114c4762a1bSJed Brown       args: -mat_partitioning_type chaco
115c4762a1bSJed Brown 
116c4762a1bSJed Brown    test:
117c4762a1bSJed Brown       suffix: 5
118c4762a1bSJed Brown       nsize: 3
119c4762a1bSJed Brown       requires: parmetis
120c4762a1bSJed Brown       args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100
121c4762a1bSJed Brown 
122c4762a1bSJed Brown    test:
123c4762a1bSJed Brown       suffix: 6
124c4762a1bSJed Brown       nsize: 3
125c4762a1bSJed Brown       requires: parmetis
126c4762a1bSJed 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
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    test:
129c4762a1bSJed Brown       suffix: 7
130c4762a1bSJed Brown       nsize: 2
131c4762a1bSJed Brown       requires: parmetis
132c4762a1bSJed 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
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    test:
135c4762a1bSJed Brown       suffix: 8
136c4762a1bSJed Brown       nsize: 2
137c4762a1bSJed Brown       requires: parmetis
138c4762a1bSJed Brown       args: -mat_partitioning_type parmetis -mat_partitioning_nparts 3 -test_use_edge_weights 1
139c4762a1bSJed Brown 
140c4762a1bSJed Brown    test:
141c4762a1bSJed Brown       suffix: 9
142c4762a1bSJed Brown       nsize: 2
143c4762a1bSJed Brown       requires: ptscotch
144c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch -mat_partitioning_nparts 3 -test_use_edge_weights 1 -mat_partitioning_ptscotch_proc_weight 0
145c4762a1bSJed Brown 
146c4762a1bSJed Brown TEST*/
147