xref: /petsc/src/mat/tutorials/ex15.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Example of using graph partitioning to partition a graph\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*T
4*c4762a1bSJed Brown    Concepts: Mat^mat partitioning
5*c4762a1bSJed Brown    Concepts: Mat^image segmentation
6*c4762a1bSJed Brown    Processors: n
7*c4762a1bSJed Brown T*/
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown #include <petscmat.h>
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown int main(int argc, char **args)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   Mat             A;
14*c4762a1bSJed Brown   MatPartitioning part;
15*c4762a1bSJed Brown   IS              is;
16*c4762a1bSJed Brown   PetscInt        r,N = 10, start, end, *vweights;
17*c4762a1bSJed Brown   PetscBool       set_vweights=PETSC_FALSE,use_edge_weights=PETSC_FALSE;
18*c4762a1bSJed Brown   PetscMPIInt     rank;
19*c4762a1bSJed Brown   MPI_Comm        comm;
20*c4762a1bSJed Brown   PetscErrorCode  ierr;
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &args, (char*) 0, help);if (ierr) return ierr;
23*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
24*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatCreate(comm, &A);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(A, 3, NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(A, 3, NULL, 2, NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_vertex_weights",&set_vweights,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_use_edge_weights",&use_edge_weights,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown   /* Create a linear mesh */
34*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A, &start, &end);CHKERRQ(ierr);
35*c4762a1bSJed Brown   if (set_vweights) {
36*c4762a1bSJed Brown     ierr = PetscMalloc1(end-start,&vweights);CHKERRQ(ierr);
37*c4762a1bSJed Brown     for (r = start; r < end; ++r)
38*c4762a1bSJed Brown       vweights[r-start] = rank+1;
39*c4762a1bSJed Brown   }
40*c4762a1bSJed Brown   for (r = start; r < end; ++r) {
41*c4762a1bSJed Brown     if (r == 0) {
42*c4762a1bSJed Brown       PetscInt    cols[2];
43*c4762a1bSJed Brown       PetscScalar vals[2];
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown       cols[0] = r;   cols[1] = r+1;
46*c4762a1bSJed Brown       vals[0] = 1.0; vals[1] = use_edge_weights? 2.0: 1.0;
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown       ierr = MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
49*c4762a1bSJed Brown     } else if (r == N-1) {
50*c4762a1bSJed Brown       PetscInt    cols[2];
51*c4762a1bSJed Brown       PetscScalar vals[2];
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown       cols[0] = r-1; cols[1] = r;
54*c4762a1bSJed Brown       vals[0] = use_edge_weights? 3.0:1.0; vals[1] = 1.0;
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown       ierr = MatSetValues(A, 1, &r, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
57*c4762a1bSJed Brown     } else {
58*c4762a1bSJed Brown       PetscInt    cols[3];
59*c4762a1bSJed Brown       PetscScalar vals[3];
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown       cols[0] = r-1; cols[1] = r;   cols[2] = r+1;
62*c4762a1bSJed Brown       /* ADJ matrix needs to be symmetric */
63*c4762a1bSJed Brown       vals[0] = use_edge_weights? (cols[0]==0? 2.0:5.0):1.0;
64*c4762a1bSJed Brown       vals[1] = 1.0;
65*c4762a1bSJed Brown       vals[2] = use_edge_weights? (cols[2]==N-1? 3.0:5.0):1.0;
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown       ierr = MatSetValues(A, 1, &r, 3, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
68*c4762a1bSJed Brown     }
69*c4762a1bSJed Brown   }
70*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   ierr = MatPartitioningCreate(comm, &part);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = MatPartitioningSetAdjacency(part, A);CHKERRQ(ierr);
75*c4762a1bSJed Brown   if (set_vweights) {
76*c4762a1bSJed Brown     ierr = MatPartitioningSetVertexWeights(part,vweights);CHKERRQ(ierr);
77*c4762a1bSJed Brown   }
78*c4762a1bSJed Brown   if (use_edge_weights) {
79*c4762a1bSJed Brown     ierr = MatPartitioningSetUseEdgeWeights(part,use_edge_weights);CHKERRQ(ierr);
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown     ierr = MatPartitioningGetUseEdgeWeights(part,&use_edge_weights);CHKERRQ(ierr);
82*c4762a1bSJed Brown     if (!use_edge_weights) SETERRQ(comm,PETSC_ERR_ARG_INCOMP, "use_edge_weights flag does not setup correctly \n");
83*c4762a1bSJed Brown   }
84*c4762a1bSJed Brown   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = MatPartitioningApply(part, &is);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = ISView(is, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = PetscFinalize();
92*c4762a1bSJed Brown   return ierr;
93*c4762a1bSJed Brown }
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown /*TEST
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown    test:
99*c4762a1bSJed Brown       nsize: 3
100*c4762a1bSJed Brown       requires: parmetis
101*c4762a1bSJed Brown       args: -mat_partitioning_type parmetis
102*c4762a1bSJed Brown 
103*c4762a1bSJed Brown    test:
104*c4762a1bSJed Brown       suffix: 2
105*c4762a1bSJed Brown       nsize: 3
106*c4762a1bSJed Brown       requires: ptscotch
107*c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown    test:
110*c4762a1bSJed Brown       suffix: 3
111*c4762a1bSJed Brown       nsize: 4
112*c4762a1bSJed Brown       requires: party
113*c4762a1bSJed Brown       args: -mat_partitioning_type party
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown    test:
116*c4762a1bSJed Brown       suffix: 4
117*c4762a1bSJed Brown       nsize: 3
118*c4762a1bSJed Brown       requires: chaco
119*c4762a1bSJed Brown       args: -mat_partitioning_type chaco
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown    test:
122*c4762a1bSJed Brown       suffix: 5
123*c4762a1bSJed Brown       nsize: 3
124*c4762a1bSJed Brown       requires: parmetis
125*c4762a1bSJed Brown       args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown    test:
128*c4762a1bSJed Brown       suffix: 6
129*c4762a1bSJed Brown       nsize: 3
130*c4762a1bSJed Brown       requires: parmetis
131*c4762a1bSJed 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
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown    test:
134*c4762a1bSJed Brown       suffix: 7
135*c4762a1bSJed Brown       nsize: 2
136*c4762a1bSJed Brown       requires: parmetis
137*c4762a1bSJed 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
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown    test:
140*c4762a1bSJed Brown       suffix: 8
141*c4762a1bSJed Brown       nsize: 2
142*c4762a1bSJed Brown       requires: parmetis
143*c4762a1bSJed Brown       args: -mat_partitioning_type parmetis -mat_partitioning_nparts 3 -test_use_edge_weights 1
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown    test:
146*c4762a1bSJed Brown       suffix: 9
147*c4762a1bSJed Brown       nsize: 2
148*c4762a1bSJed Brown       requires: ptscotch
149*c4762a1bSJed Brown       args: -mat_partitioning_type ptscotch -mat_partitioning_nparts 3 -test_use_edge_weights 1 -mat_partitioning_ptscotch_proc_weight 0
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown TEST*/
152