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
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat A;
8c4762a1bSJed Brown MatPartitioning part;
9c4762a1bSJed Brown IS is;
10c4762a1bSJed Brown PetscInt r, N = 10, start, end, *vweights;
11c4762a1bSJed Brown PetscBool set_vweights = PETSC_FALSE, use_edge_weights = PETSC_FALSE;
12c4762a1bSJed Brown PetscMPIInt rank;
13c4762a1bSJed Brown MPI_Comm comm;
14c4762a1bSJed Brown
15327415f7SBarry Smith PetscFunctionBeginUser;
16*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
17c4762a1bSJed Brown comm = PETSC_COMM_WORLD;
189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
209566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A));
219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N));
229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
239566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 3, NULL));
249566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL));
259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_vertex_weights", &set_vweights, NULL));
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_use_edge_weights", &use_edge_weights, NULL));
27c4762a1bSJed Brown /* Create a linear mesh */
289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &start, &end));
29c4762a1bSJed Brown if (set_vweights) {
309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start, &vweights));
319371c9d4SSatish Balay for (r = start; r < end; ++r) vweights[r - start] = rank + 1;
32c4762a1bSJed Brown }
33c4762a1bSJed Brown for (r = start; r < end; ++r) {
34c4762a1bSJed Brown if (r == 0) {
35c4762a1bSJed Brown PetscInt cols[2];
36c4762a1bSJed Brown PetscScalar vals[2];
37c4762a1bSJed Brown
389371c9d4SSatish Balay cols[0] = r;
399371c9d4SSatish Balay cols[1] = r + 1;
409371c9d4SSatish Balay vals[0] = 1.0;
419371c9d4SSatish Balay vals[1] = use_edge_weights ? 2.0 : 1.0;
42c4762a1bSJed Brown
439566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
44c4762a1bSJed Brown } else if (r == N - 1) {
45c4762a1bSJed Brown PetscInt cols[2];
46c4762a1bSJed Brown PetscScalar vals[2];
47c4762a1bSJed Brown
489371c9d4SSatish Balay cols[0] = r - 1;
499371c9d4SSatish Balay cols[1] = r;
509371c9d4SSatish Balay vals[0] = use_edge_weights ? 3.0 : 1.0;
519371c9d4SSatish Balay vals[1] = 1.0;
52c4762a1bSJed Brown
539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES));
54c4762a1bSJed Brown } else {
55c4762a1bSJed Brown PetscInt cols[3];
56c4762a1bSJed Brown PetscScalar vals[3];
57c4762a1bSJed Brown
589371c9d4SSatish Balay cols[0] = r - 1;
599371c9d4SSatish Balay cols[1] = r;
609371c9d4SSatish Balay 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
669566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES));
67c4762a1bSJed Brown }
68c4762a1bSJed Brown }
699566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown
729566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &part));
739566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A));
741baa6e33SBarry Smith if (set_vweights) PetscCall(MatPartitioningSetVertexWeights(part, vweights));
75c4762a1bSJed Brown if (use_edge_weights) {
769566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetUseEdgeWeights(part, use_edge_weights));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(MatPartitioningGetUseEdgeWeights(part, &use_edge_weights));
7928b400f6SJacob Faibussowitsch PetscCheck(use_edge_weights, comm, PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly ");
80c4762a1bSJed Brown }
819566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part));
829566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part, &is));
839566063dSJacob Faibussowitsch PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
859566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part));
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch return 0;
90c4762a1bSJed Brown }
91c4762a1bSJed Brown
92c4762a1bSJed Brown /*TEST
93c4762a1bSJed Brown
94c4762a1bSJed Brown test:
95c4762a1bSJed Brown nsize: 3
96c4762a1bSJed Brown requires: parmetis
97c4762a1bSJed Brown args: -mat_partitioning_type parmetis
98c4762a1bSJed Brown
99c4762a1bSJed Brown test:
100c4762a1bSJed Brown suffix: 2
101c4762a1bSJed Brown nsize: 3
102c4762a1bSJed Brown requires: ptscotch
103c4762a1bSJed Brown args: -mat_partitioning_type ptscotch
104c4762a1bSJed Brown
105c4762a1bSJed Brown test:
106c4762a1bSJed Brown suffix: 3
107c4762a1bSJed Brown nsize: 4
108c4762a1bSJed Brown requires: party
109c4762a1bSJed Brown args: -mat_partitioning_type party
110c4762a1bSJed Brown
111c4762a1bSJed Brown test:
112c4762a1bSJed Brown suffix: 4
113c4762a1bSJed Brown nsize: 3
114c4762a1bSJed Brown requires: chaco
115c4762a1bSJed Brown args: -mat_partitioning_type chaco
116c4762a1bSJed Brown
117c4762a1bSJed Brown test:
118c4762a1bSJed Brown suffix: 5
119c4762a1bSJed Brown nsize: 3
120c4762a1bSJed Brown requires: parmetis
121c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100
122c4762a1bSJed Brown
123c4762a1bSJed Brown test:
124c4762a1bSJed Brown suffix: 6
125c4762a1bSJed Brown nsize: 3
126c4762a1bSJed Brown requires: parmetis
127c4762a1bSJed 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
128c4762a1bSJed Brown
129c4762a1bSJed Brown test:
130c4762a1bSJed Brown suffix: 7
131c4762a1bSJed Brown nsize: 2
132c4762a1bSJed Brown requires: parmetis
133c4762a1bSJed 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
134c4762a1bSJed Brown
135c4762a1bSJed Brown test:
136c4762a1bSJed Brown suffix: 8
137c4762a1bSJed Brown nsize: 2
138c4762a1bSJed Brown requires: parmetis
139c4762a1bSJed Brown args: -mat_partitioning_type parmetis -mat_partitioning_nparts 3 -test_use_edge_weights 1
140c4762a1bSJed Brown
141c4762a1bSJed Brown test:
142c4762a1bSJed Brown suffix: 9
143c4762a1bSJed Brown nsize: 2
144c4762a1bSJed Brown requires: ptscotch
145c4762a1bSJed Brown args: -mat_partitioning_type ptscotch -mat_partitioning_nparts 3 -test_use_edge_weights 1 -mat_partitioning_ptscotch_proc_weight 0
146c4762a1bSJed Brown
147c4762a1bSJed Brown TEST*/
148