xref: /petsc/src/dm/partitioner/tests/ex33.c (revision d7cc930e14e615e9907267aaa472dd0ccceeab82)
1abe9303eSLisandro Dalcin static char help[] = "Tests PetscPartitioner.\n\n";
2abe9303eSLisandro Dalcin 
3abe9303eSLisandro Dalcin #include <petscpartitioner.h>
4abe9303eSLisandro Dalcin 
5abe9303eSLisandro Dalcin int main(int argc, char **argv)
6abe9303eSLisandro Dalcin {
7abe9303eSLisandro Dalcin   PetscErrorCode   ierr;
8abe9303eSLisandro Dalcin   PetscPartitioner p;
9*d7cc930eSLisandro Dalcin   PetscSection     partSection, vertexSection = NULL, targetSection = NULL;
10abe9303eSLisandro Dalcin   IS               partition,is;
11abe9303eSLisandro Dalcin   PetscMPIInt      size,rank;
12*d7cc930eSLisandro Dalcin   PetscInt         nparts,i;
13abe9303eSLisandro Dalcin   PetscInt         nv = 4;
14abe9303eSLisandro Dalcin   PetscInt         vv[5] = {0,2,4,6,8};
15abe9303eSLisandro Dalcin   PetscInt         vadj[8] = {3,1,0,2,1,3,2,0};
1656991f02SLisandro Dalcin   PetscBool        sequential;
17*d7cc930eSLisandro Dalcin   PetscBool        vwgts = PETSC_FALSE;
18*d7cc930eSLisandro Dalcin   PetscBool        pwgts = PETSC_FALSE;
19abe9303eSLisandro Dalcin 
20abe9303eSLisandro Dalcin   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
21abe9303eSLisandro Dalcin   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
22abe9303eSLisandro Dalcin   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
23*d7cc930eSLisandro Dalcin   nparts = size;
24*d7cc930eSLisandro Dalcin   ierr = PetscOptionsGetInt(NULL,NULL,"-nparts",&nparts,NULL);CHKERRQ(ierr);
25*d7cc930eSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-vwgts",&vwgts,NULL);CHKERRQ(ierr);
26*d7cc930eSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-pwgts",&pwgts,NULL);CHKERRQ(ierr);
27abe9303eSLisandro Dalcin 
28abe9303eSLisandro Dalcin   /* create PetscPartitioner */
29abe9303eSLisandro Dalcin   ierr = PetscPartitionerCreate(PETSC_COMM_WORLD,&p);CHKERRQ(ierr);
30abe9303eSLisandro Dalcin   ierr = PetscPartitionerSetType(p,PETSCPARTITIONERSIMPLE);CHKERRQ(ierr);
31abe9303eSLisandro Dalcin   ierr = PetscPartitionerSetFromOptions(p);CHKERRQ(ierr);
32abe9303eSLisandro Dalcin 
33*d7cc930eSLisandro Dalcin   /* create partition section */
34*d7cc930eSLisandro Dalcin   ierr = PetscSectionCreate(PETSC_COMM_WORLD,&partSection);CHKERRQ(ierr);
35*d7cc930eSLisandro Dalcin 
36*d7cc930eSLisandro Dalcin   if (vwgts) { /* create vertex weights section */
37*d7cc930eSLisandro Dalcin     ierr = PetscSectionCreate(PETSC_COMM_WORLD,&vertexSection);CHKERRQ(ierr);
38*d7cc930eSLisandro Dalcin     ierr = PetscSectionSetChart(vertexSection,0,nv);CHKERRQ(ierr);
39*d7cc930eSLisandro Dalcin     for (i = 0; i< nv; i++) {ierr = PetscSectionSetDof(vertexSection,i,1);CHKERRQ(ierr);}
40*d7cc930eSLisandro Dalcin     ierr = PetscSectionSetUp(vertexSection);CHKERRQ(ierr);
41*d7cc930eSLisandro Dalcin   }
42*d7cc930eSLisandro Dalcin 
43*d7cc930eSLisandro Dalcin   if (pwgts) { /* create partition weights section */
44*d7cc930eSLisandro Dalcin     ierr = PetscSectionCreate(PETSC_COMM_WORLD,&targetSection);CHKERRQ(ierr);
45*d7cc930eSLisandro Dalcin     ierr = PetscSectionSetChart(targetSection,0,nparts);CHKERRQ(ierr);
46*d7cc930eSLisandro Dalcin     for (i = 0; i< nparts; i++) {ierr = PetscSectionSetDof(targetSection,i,1);CHKERRQ(ierr);}
47*d7cc930eSLisandro Dalcin     ierr = PetscSectionSetUp(targetSection);CHKERRQ(ierr);
48*d7cc930eSLisandro Dalcin   }
49*d7cc930eSLisandro Dalcin 
50*d7cc930eSLisandro Dalcin #if defined(PETSC_USE_LOG)
51*d7cc930eSLisandro Dalcin   { /* Test logging */
52*d7cc930eSLisandro Dalcin     PetscLogEvent event;
53*d7cc930eSLisandro Dalcin 
54*d7cc930eSLisandro Dalcin     ierr = PetscLogEventRegister("MyPartitionerEvent",PETSCPARTITIONER_CLASSID,&event);CHKERRQ(ierr);
55*d7cc930eSLisandro Dalcin     { /* PetscLogEventExcludeClass is broken, new events are not deactivated */
56*d7cc930eSLisandro Dalcin       char      logList[256];
57*d7cc930eSLisandro Dalcin       PetscBool opt,pkg;
58*d7cc930eSLisandro Dalcin 
59*d7cc930eSLisandro Dalcin       ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr);
60*d7cc930eSLisandro Dalcin       if (opt) {
61*d7cc930eSLisandro Dalcin         ierr = PetscStrInList("partitioner",logList,',',&pkg);CHKERRQ(ierr);
62*d7cc930eSLisandro Dalcin         if (pkg) {ierr = PetscLogEventExcludeClass(PETSCPARTITIONER_CLASSID);CHKERRQ(ierr);}
63*d7cc930eSLisandro Dalcin       }
64*d7cc930eSLisandro Dalcin     }
65*d7cc930eSLisandro Dalcin     ierr = PetscLogEventBegin(event,p,NULL,NULL,NULL);CHKERRQ(ierr);
66*d7cc930eSLisandro Dalcin     ierr = PetscLogEventEnd(event,p,NULL,NULL,NULL);CHKERRQ(ierr);
67*d7cc930eSLisandro Dalcin   }
68*d7cc930eSLisandro Dalcin #endif
69*d7cc930eSLisandro Dalcin 
70*d7cc930eSLisandro Dalcin   /* test setup and reset */
71*d7cc930eSLisandro Dalcin   ierr = PetscPartitionerSetUp(p);CHKERRQ(ierr);
72*d7cc930eSLisandro Dalcin   ierr = PetscPartitionerReset(p);CHKERRQ(ierr);
73*d7cc930eSLisandro Dalcin 
74abe9303eSLisandro Dalcin   /* test partitioning an empty graph */
75*d7cc930eSLisandro Dalcin   ierr = PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition);CHKERRQ(ierr);
76abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"NULL SECTION");
77abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
78abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
79abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"NULL PARTITION");
80abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
81abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
82abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
83abe9303eSLisandro Dalcin 
84*d7cc930eSLisandro Dalcin   /* test view from options */
85*d7cc930eSLisandro Dalcin   ierr = PetscPartitionerViewFromOptions(p,NULL,"-part_view");CHKERRQ(ierr);
86*d7cc930eSLisandro Dalcin 
87abe9303eSLisandro Dalcin   /* test partitioning a graph on one process only (not master) */
88abe9303eSLisandro Dalcin   if (rank == size - 1) {
89*d7cc930eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,nparts,nv,vv,vadj,vertexSection,targetSection,partSection,&partition);CHKERRQ(ierr);
90abe9303eSLisandro Dalcin   } else {
91*d7cc930eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition);CHKERRQ(ierr);
92abe9303eSLisandro Dalcin   }
93abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"SEQ SECTION");
94abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
95abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
96abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"SEQ PARTITION");
97abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
98abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
99abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
100abe9303eSLisandro Dalcin 
10156991f02SLisandro Dalcin   ierr = PetscObjectTypeCompareAny((PetscObject)p,&sequential,PETSCPARTITIONERCHACO,NULL);CHKERRQ(ierr);
10256991f02SLisandro Dalcin   if (sequential) goto finally;
10356991f02SLisandro Dalcin 
104abe9303eSLisandro Dalcin   /* test partitioning a graph on a subset of the processess only */
105abe9303eSLisandro Dalcin   if (rank%2) {
106*d7cc930eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,nparts,0,NULL,NULL,NULL,targetSection,partSection,&partition);CHKERRQ(ierr);
107abe9303eSLisandro Dalcin   } else {
108abe9303eSLisandro Dalcin     PetscInt i,totv = nv*((size+1)/2),*pvadj;
109abe9303eSLisandro Dalcin 
110abe9303eSLisandro Dalcin     ierr = PetscMalloc1(2*nv,&pvadj);CHKERRQ(ierr);
111abe9303eSLisandro Dalcin     for (i = 0; i < nv; i++) {
112abe9303eSLisandro Dalcin       pvadj[2*i]   = (nv*(rank/2) + totv + i - 1)%totv;
113abe9303eSLisandro Dalcin       pvadj[2*i+1] = (nv*(rank/2) + totv + i + 1)%totv;
114abe9303eSLisandro Dalcin     }
115*d7cc930eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,nparts,nv,vv,pvadj,NULL,targetSection,partSection,&partition);CHKERRQ(ierr);
116abe9303eSLisandro Dalcin     ierr = PetscFree(pvadj);CHKERRQ(ierr);
117abe9303eSLisandro Dalcin   }
118abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"PARVOID SECTION");
119abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
120abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
121abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"PARVOID PARTITION");
122abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
123abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
124abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
125abe9303eSLisandro Dalcin 
12656991f02SLisandro Dalcin finally:
127abe9303eSLisandro Dalcin   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
128*d7cc930eSLisandro Dalcin   ierr = PetscSectionDestroy(&vertexSection);CHKERRQ(ierr);
129*d7cc930eSLisandro Dalcin   ierr = PetscSectionDestroy(&targetSection);CHKERRQ(ierr);
130abe9303eSLisandro Dalcin   ierr = PetscPartitionerDestroy(&p);CHKERRQ(ierr);
131abe9303eSLisandro Dalcin   ierr = PetscFinalize();
132abe9303eSLisandro Dalcin   return ierr;
133abe9303eSLisandro Dalcin }
134abe9303eSLisandro Dalcin 
135abe9303eSLisandro Dalcin /*TEST
136abe9303eSLisandro Dalcin 
137abe9303eSLisandro Dalcin   test:
138*d7cc930eSLisandro Dalcin     suffix: default
139*d7cc930eSLisandro Dalcin 
140*d7cc930eSLisandro Dalcin   testset:
141*d7cc930eSLisandro Dalcin     requires: defined(PETSC_USE_LOG)
142*d7cc930eSLisandro Dalcin     args: -petscpartitioner_type simple -log_summary
143*d7cc930eSLisandro Dalcin     filter: grep MyPartitionerEvent | cut -d ' ' -f 1
144*d7cc930eSLisandro Dalcin     test:
145*d7cc930eSLisandro Dalcin        suffix: log_include
146*d7cc930eSLisandro Dalcin     test:
147*d7cc930eSLisandro Dalcin       suffix: log_exclude
148*d7cc930eSLisandro Dalcin       args: -log_exclude partitioner
149*d7cc930eSLisandro Dalcin 
150*d7cc930eSLisandro Dalcin   test:
151abe9303eSLisandro Dalcin     suffix: simple
152abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
153*d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -petscpartitioner_type simple -petscpartitioner_view
154*d7cc930eSLisandro Dalcin 
155*d7cc930eSLisandro Dalcin   test:
156*d7cc930eSLisandro Dalcin     suffix: shell
157*d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
158*d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type shell -petscpartitioner_shell_random -petscpartitioner_view
159abe9303eSLisandro Dalcin 
160abe9303eSLisandro Dalcin   test:
1615130b6c4SLisandro Dalcin     suffix: gather
1625130b6c4SLisandro Dalcin     nsize: {{1 2 3}separate output}
1635130b6c4SLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type gather -petscpartitioner_view -petscpartitioner_view_graph
1645130b6c4SLisandro Dalcin 
1655130b6c4SLisandro Dalcin   test:
166abe9303eSLisandro Dalcin     requires: parmetis
167abe9303eSLisandro Dalcin     suffix: parmetis
168abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
169*d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}} -vwgts {{false true}}
170*d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -petscpartitioner_view -petscpartitioner_view_graph
171*d7cc930eSLisandro Dalcin 
172*d7cc930eSLisandro Dalcin   test:
173*d7cc930eSLisandro Dalcin     requires: parmetis
174*d7cc930eSLisandro Dalcin     suffix: parmetis_type
175*d7cc930eSLisandro Dalcin     nsize: {{1 2}}
176*d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -part_view
177*d7cc930eSLisandro Dalcin     args: -petscpartitioner_parmetis_type {{kway rb}separate output}
178*d7cc930eSLisandro Dalcin     filter: grep "ParMetis type"
179abe9303eSLisandro Dalcin 
180abe9303eSLisandro Dalcin   test:
181abe9303eSLisandro Dalcin     requires: ptscotch
182abe9303eSLisandro Dalcin     suffix: ptscotch
183abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
184*d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -vwgts {{false true}}
185*d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_view_graph
186abe9303eSLisandro Dalcin 
18756991f02SLisandro Dalcin   test:
188*d7cc930eSLisandro Dalcin     requires: ptscotch
189*d7cc930eSLisandro Dalcin     suffix: ptscotch_strategy
190*d7cc930eSLisandro Dalcin     nsize: {{1 2}}
191*d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -part_view
192*d7cc930eSLisandro Dalcin     args: -petscpartitioner_ptscotch_strategy {{DEFAULT QUALITY SPEED BALANCE SAFETY SCALABILITY RECURSIVE REMAP}separate output}
193*d7cc930eSLisandro Dalcin     filter: grep "partitioning strategy"
194*d7cc930eSLisandro Dalcin 
195*d7cc930eSLisandro Dalcin   test:
19656991f02SLisandro Dalcin     requires: chaco
19756991f02SLisandro Dalcin     suffix: chaco
19856991f02SLisandro Dalcin     nsize: {{1 2 3}separate output}
199*d7cc930eSLisandro Dalcin     args: -nparts {{1}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
200*d7cc930eSLisandro Dalcin 
201*d7cc930eSLisandro Dalcin   test:
202*d7cc930eSLisandro Dalcin     TODO: non reproducible (uses C stdlib rand())
203*d7cc930eSLisandro Dalcin     requires: chaco
204*d7cc930eSLisandro Dalcin     suffix: chaco
205*d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
206*d7cc930eSLisandro Dalcin     args: -nparts {{2 3}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
20756991f02SLisandro Dalcin 
208abe9303eSLisandro Dalcin TEST*/
209