xref: /petsc/src/dm/partitioner/tests/ex33.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscPartitioner p;
8d7cc930eSLisandro Dalcin   PetscSection     partSection, vertexSection = NULL, targetSection = NULL;
9abe9303eSLisandro Dalcin   IS               partition,is;
10abe9303eSLisandro Dalcin   PetscMPIInt      size,rank;
11d7cc930eSLisandro Dalcin   PetscInt         nparts,i;
12abe9303eSLisandro Dalcin   PetscInt         nv = 4;
13abe9303eSLisandro Dalcin   PetscInt         vv[5] = {0,2,4,6,8};
14abe9303eSLisandro Dalcin   PetscInt         vadj[8] = {3,1,0,2,1,3,2,0};
1556991f02SLisandro Dalcin   PetscBool        sequential;
16d7cc930eSLisandro Dalcin   PetscBool        vwgts = PETSC_FALSE;
17d7cc930eSLisandro Dalcin   PetscBool        pwgts = PETSC_FALSE;
18abe9303eSLisandro Dalcin 
19*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23d7cc930eSLisandro Dalcin   nparts = size;
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nparts",&nparts,NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-vwgts",&vwgts,NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-pwgts",&pwgts,NULL));
27abe9303eSLisandro Dalcin 
28abe9303eSLisandro Dalcin   /* create PetscPartitioner */
299566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerCreate(PETSC_COMM_WORLD,&p));
309566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(p,PETSCPARTITIONERSIMPLE));
319566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(p));
32abe9303eSLisandro Dalcin 
33d7cc930eSLisandro Dalcin   /* create partition section */
349566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PETSC_COMM_WORLD,&partSection));
35d7cc930eSLisandro Dalcin 
36d7cc930eSLisandro Dalcin   if (vwgts) { /* create vertex weights section */
379566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PETSC_COMM_WORLD,&vertexSection));
389566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(vertexSection,0,nv));
399566063dSJacob Faibussowitsch     for (i = 0; i< nv; i++) PetscCall(PetscSectionSetDof(vertexSection,i,1));
409566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(vertexSection));
41d7cc930eSLisandro Dalcin   }
42d7cc930eSLisandro Dalcin 
43d7cc930eSLisandro Dalcin   if (pwgts) { /* create partition weights section */
449566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PETSC_COMM_WORLD,&targetSection));
459566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(targetSection,0,nparts));
469566063dSJacob Faibussowitsch     for (i = 0; i< nparts; i++) PetscCall(PetscSectionSetDof(targetSection,i,1));
479566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(targetSection));
48d7cc930eSLisandro Dalcin   }
49d7cc930eSLisandro Dalcin 
50d7cc930eSLisandro Dalcin #if defined(PETSC_USE_LOG)
51d7cc930eSLisandro Dalcin   { /* Test logging */
52d7cc930eSLisandro Dalcin     PetscLogEvent event;
53d7cc930eSLisandro Dalcin 
549566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("MyPartitionerEvent",PETSCPARTITIONER_CLASSID,&event));
55d7cc930eSLisandro Dalcin     { /* PetscLogEventExcludeClass is broken, new events are not deactivated */
56d7cc930eSLisandro Dalcin       char      logList[256];
57d7cc930eSLisandro Dalcin       PetscBool opt,pkg;
58d7cc930eSLisandro Dalcin 
599566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
60d7cc930eSLisandro Dalcin       if (opt) {
619566063dSJacob Faibussowitsch         PetscCall(PetscStrInList("partitioner",logList,',',&pkg));
629566063dSJacob Faibussowitsch         if (pkg) PetscCall(PetscLogEventExcludeClass(PETSCPARTITIONER_CLASSID));
63d7cc930eSLisandro Dalcin       }
64d7cc930eSLisandro Dalcin     }
659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(event,p,NULL,NULL,NULL));
669566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(event,p,NULL,NULL,NULL));
67d7cc930eSLisandro Dalcin   }
68d7cc930eSLisandro Dalcin #endif
69d7cc930eSLisandro Dalcin 
70d7cc930eSLisandro Dalcin   /* test setup and reset */
719566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetUp(p));
729566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerReset(p));
73d7cc930eSLisandro Dalcin 
74abe9303eSLisandro Dalcin   /* test partitioning an empty graph */
759566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition));
769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)partSection,"NULL SECTION"));
779566063dSJacob Faibussowitsch   PetscCall(PetscSectionView(partSection,NULL));
789566063dSJacob Faibussowitsch   PetscCall(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)is,"NULL PARTITION"));
809566063dSJacob Faibussowitsch   PetscCall(ISView(is,NULL));
819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
829566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&partition));
83abe9303eSLisandro Dalcin 
84d7cc930eSLisandro Dalcin   /* test view from options */
859566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerViewFromOptions(p,NULL,"-part_view"));
86d7cc930eSLisandro Dalcin 
879dddd249SSatish Balay   /* test partitioning a graph on one process only (not main) */
88abe9303eSLisandro Dalcin   if (rank == size - 1) {
899566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerPartition(p,nparts,nv,vv,vadj,vertexSection,targetSection,partSection,&partition));
90abe9303eSLisandro Dalcin   } else {
919566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition));
92abe9303eSLisandro Dalcin   }
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)partSection,"SEQ SECTION"));
949566063dSJacob Faibussowitsch   PetscCall(PetscSectionView(partSection,NULL));
959566063dSJacob Faibussowitsch   PetscCall(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
969566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)is,"SEQ PARTITION"));
979566063dSJacob Faibussowitsch   PetscCall(ISView(is,NULL));
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&partition));
100abe9303eSLisandro Dalcin 
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)p,&sequential,PETSCPARTITIONERCHACO,NULL));
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) {
1069566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerPartition(p,nparts,0,NULL,NULL,NULL,targetSection,partSection,&partition));
107abe9303eSLisandro Dalcin   } else {
108abe9303eSLisandro Dalcin     PetscInt i,totv = nv*((size+1)/2),*pvadj;
109abe9303eSLisandro Dalcin 
1109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2*nv,&pvadj));
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     }
1159566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerPartition(p,nparts,nv,vv,pvadj,NULL,targetSection,partSection,&partition));
1169566063dSJacob Faibussowitsch     PetscCall(PetscFree(pvadj));
117abe9303eSLisandro Dalcin   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)partSection,"PARVOID SECTION"));
1199566063dSJacob Faibussowitsch   PetscCall(PetscSectionView(partSection,NULL));
1209566063dSJacob Faibussowitsch   PetscCall(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
1219566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)is,"PARVOID PARTITION"));
1229566063dSJacob Faibussowitsch   PetscCall(ISView(is,NULL));
1239566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
1249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&partition));
125abe9303eSLisandro Dalcin 
12656991f02SLisandro Dalcin finally:
1279566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&partSection));
1289566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&vertexSection));
1299566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&targetSection));
1309566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDestroy(&p));
1319566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
132b122ec5aSJacob Faibussowitsch   return 0;
133abe9303eSLisandro Dalcin }
134abe9303eSLisandro Dalcin 
135abe9303eSLisandro Dalcin /*TEST
136abe9303eSLisandro Dalcin 
137abe9303eSLisandro Dalcin   test:
138d7cc930eSLisandro Dalcin     suffix: default
139d7cc930eSLisandro Dalcin 
140d7cc930eSLisandro Dalcin   testset:
141dfd57a17SPierre Jolivet     requires: defined(PETSC_USE_LOG)
142d7cc930eSLisandro Dalcin     args: -petscpartitioner_type simple -log_summary
14379229e12SLisandro Dalcin     filter: grep MyPartitionerEvent | cut -d " " -f 1
144d7cc930eSLisandro Dalcin     test:
145d7cc930eSLisandro Dalcin        suffix: log_include
146d7cc930eSLisandro Dalcin     test:
147d7cc930eSLisandro Dalcin       suffix: log_exclude
148d7cc930eSLisandro Dalcin       args: -log_exclude partitioner
149d7cc930eSLisandro Dalcin 
150d7cc930eSLisandro Dalcin   test:
151abe9303eSLisandro Dalcin     suffix: simple
152abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
153d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -petscpartitioner_type simple -petscpartitioner_view
154d7cc930eSLisandro Dalcin 
155d7cc930eSLisandro Dalcin   test:
156d7cc930eSLisandro Dalcin     suffix: shell
157d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
158d7cc930eSLisandro 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}
169d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}} -vwgts {{false true}}
170d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -petscpartitioner_view -petscpartitioner_view_graph
171d7cc930eSLisandro Dalcin 
172d7cc930eSLisandro Dalcin   test:
173d7cc930eSLisandro Dalcin     requires: parmetis
174d7cc930eSLisandro Dalcin     suffix: parmetis_type
175d7cc930eSLisandro Dalcin     nsize: {{1 2}}
176d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -part_view
177d7cc930eSLisandro Dalcin     args: -petscpartitioner_parmetis_type {{kway rb}separate output}
178d7cc930eSLisandro 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}
184d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -vwgts {{false true}}
185d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_view_graph
186abe9303eSLisandro Dalcin 
18756991f02SLisandro Dalcin   test:
188d7cc930eSLisandro Dalcin     requires: ptscotch
189d7cc930eSLisandro Dalcin     suffix: ptscotch_strategy
190d7cc930eSLisandro Dalcin     nsize: {{1 2}}
191d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -part_view
192d7cc930eSLisandro Dalcin     args: -petscpartitioner_ptscotch_strategy {{DEFAULT QUALITY SPEED BALANCE SAFETY SCALABILITY RECURSIVE REMAP}separate output}
193d7cc930eSLisandro Dalcin     filter: grep "partitioning strategy"
194d7cc930eSLisandro Dalcin 
195d7cc930eSLisandro Dalcin   test:
19656991f02SLisandro Dalcin     requires: chaco
19756991f02SLisandro Dalcin     suffix: chaco
19856991f02SLisandro Dalcin     nsize: {{1 2 3}separate output}
199d7cc930eSLisandro Dalcin     args: -nparts {{1}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
200d7cc930eSLisandro Dalcin 
201d7cc930eSLisandro Dalcin   test:
202d7cc930eSLisandro Dalcin     TODO: non reproducible (uses C stdlib rand())
203d7cc930eSLisandro Dalcin     requires: chaco
204d7cc930eSLisandro Dalcin     suffix: chaco
205d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
206d7cc930eSLisandro Dalcin     args: -nparts {{2 3}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
20756991f02SLisandro Dalcin 
208abe9303eSLisandro Dalcin TEST*/
209