xref: /petsc/src/dm/partitioner/tests/ex33.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22d7cc930eSLisandro Dalcin   nparts = size;
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nparts",&nparts,NULL));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-vwgts",&vwgts,NULL));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-pwgts",&pwgts,NULL));
26abe9303eSLisandro Dalcin 
27abe9303eSLisandro Dalcin   /* create PetscPartitioner */
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerCreate(PETSC_COMM_WORLD,&p));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetType(p,PETSCPARTITIONERSIMPLE));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(p));
31abe9303eSLisandro Dalcin 
32d7cc930eSLisandro Dalcin   /* create partition section */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(PETSC_COMM_WORLD,&partSection));
34d7cc930eSLisandro Dalcin 
35d7cc930eSLisandro Dalcin   if (vwgts) { /* create vertex weights section */
365f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(PETSC_COMM_WORLD,&vertexSection));
375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(vertexSection,0,nv));
385f80ce2aSJacob Faibussowitsch     for (i = 0; i< nv; i++) CHKERRQ(PetscSectionSetDof(vertexSection,i,1));
395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(vertexSection));
40d7cc930eSLisandro Dalcin   }
41d7cc930eSLisandro Dalcin 
42d7cc930eSLisandro Dalcin   if (pwgts) { /* create partition weights section */
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(PETSC_COMM_WORLD,&targetSection));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(targetSection,0,nparts));
455f80ce2aSJacob Faibussowitsch     for (i = 0; i< nparts; i++) CHKERRQ(PetscSectionSetDof(targetSection,i,1));
465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(targetSection));
47d7cc930eSLisandro Dalcin   }
48d7cc930eSLisandro Dalcin 
49d7cc930eSLisandro Dalcin #if defined(PETSC_USE_LOG)
50d7cc930eSLisandro Dalcin   { /* Test logging */
51d7cc930eSLisandro Dalcin     PetscLogEvent event;
52d7cc930eSLisandro Dalcin 
535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventRegister("MyPartitionerEvent",PETSCPARTITIONER_CLASSID,&event));
54d7cc930eSLisandro Dalcin     { /* PetscLogEventExcludeClass is broken, new events are not deactivated */
55d7cc930eSLisandro Dalcin       char      logList[256];
56d7cc930eSLisandro Dalcin       PetscBool opt,pkg;
57d7cc930eSLisandro Dalcin 
585f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
59d7cc930eSLisandro Dalcin       if (opt) {
605f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscStrInList("partitioner",logList,',',&pkg));
615f80ce2aSJacob Faibussowitsch         if (pkg) CHKERRQ(PetscLogEventExcludeClass(PETSCPARTITIONER_CLASSID));
62d7cc930eSLisandro Dalcin       }
63d7cc930eSLisandro Dalcin     }
645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventBegin(event,p,NULL,NULL,NULL));
655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventEnd(event,p,NULL,NULL,NULL));
66d7cc930eSLisandro Dalcin   }
67d7cc930eSLisandro Dalcin #endif
68d7cc930eSLisandro Dalcin 
69d7cc930eSLisandro Dalcin   /* test setup and reset */
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetUp(p));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerReset(p));
72d7cc930eSLisandro Dalcin 
73abe9303eSLisandro Dalcin   /* test partitioning an empty graph */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)partSection,"NULL SECTION"));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionView(partSection,NULL));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)is,"NULL PARTITION"));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(is,NULL));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&partition));
82abe9303eSLisandro Dalcin 
83d7cc930eSLisandro Dalcin   /* test view from options */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerViewFromOptions(p,NULL,"-part_view"));
85d7cc930eSLisandro Dalcin 
869dddd249SSatish Balay   /* test partitioning a graph on one process only (not main) */
87abe9303eSLisandro Dalcin   if (rank == size - 1) {
885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerPartition(p,nparts,nv,vv,vadj,vertexSection,targetSection,partSection,&partition));
89abe9303eSLisandro Dalcin   } else {
905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerPartition(p,nparts,0,NULL,NULL,vertexSection,targetSection,partSection,&partition));
91abe9303eSLisandro Dalcin   }
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)partSection,"SEQ SECTION"));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionView(partSection,NULL));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)is,"SEQ PARTITION"));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(is,NULL));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&partition));
99abe9303eSLisandro Dalcin 
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompareAny((PetscObject)p,&sequential,PETSCPARTITIONERCHACO,NULL));
10156991f02SLisandro Dalcin   if (sequential) goto finally;
10256991f02SLisandro Dalcin 
103abe9303eSLisandro Dalcin   /* test partitioning a graph on a subset of the processess only */
104abe9303eSLisandro Dalcin   if (rank%2) {
1055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerPartition(p,nparts,0,NULL,NULL,NULL,targetSection,partSection,&partition));
106abe9303eSLisandro Dalcin   } else {
107abe9303eSLisandro Dalcin     PetscInt i,totv = nv*((size+1)/2),*pvadj;
108abe9303eSLisandro Dalcin 
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2*nv,&pvadj));
110abe9303eSLisandro Dalcin     for (i = 0; i < nv; i++) {
111abe9303eSLisandro Dalcin       pvadj[2*i]   = (nv*(rank/2) + totv + i - 1)%totv;
112abe9303eSLisandro Dalcin       pvadj[2*i+1] = (nv*(rank/2) + totv + i + 1)%totv;
113abe9303eSLisandro Dalcin     }
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerPartition(p,nparts,nv,vv,pvadj,NULL,targetSection,partSection,&partition));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(pvadj));
116abe9303eSLisandro Dalcin   }
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)partSection,"PARVOID SECTION"));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionView(partSection,NULL));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)is,"PARVOID PARTITION"));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(is,NULL));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&partition));
124abe9303eSLisandro Dalcin 
12556991f02SLisandro Dalcin finally:
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&partSection));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&vertexSection));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&targetSection));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerDestroy(&p));
130*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
131*b122ec5aSJacob Faibussowitsch   return 0;
132abe9303eSLisandro Dalcin }
133abe9303eSLisandro Dalcin 
134abe9303eSLisandro Dalcin /*TEST
135abe9303eSLisandro Dalcin 
136abe9303eSLisandro Dalcin   test:
137d7cc930eSLisandro Dalcin     suffix: default
138d7cc930eSLisandro Dalcin 
139d7cc930eSLisandro Dalcin   testset:
140dfd57a17SPierre Jolivet     requires: defined(PETSC_USE_LOG)
141d7cc930eSLisandro Dalcin     args: -petscpartitioner_type simple -log_summary
14279229e12SLisandro Dalcin     filter: grep MyPartitionerEvent | cut -d " " -f 1
143d7cc930eSLisandro Dalcin     test:
144d7cc930eSLisandro Dalcin        suffix: log_include
145d7cc930eSLisandro Dalcin     test:
146d7cc930eSLisandro Dalcin       suffix: log_exclude
147d7cc930eSLisandro Dalcin       args: -log_exclude partitioner
148d7cc930eSLisandro Dalcin 
149d7cc930eSLisandro Dalcin   test:
150abe9303eSLisandro Dalcin     suffix: simple
151abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
152d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -petscpartitioner_type simple -petscpartitioner_view
153d7cc930eSLisandro Dalcin 
154d7cc930eSLisandro Dalcin   test:
155d7cc930eSLisandro Dalcin     suffix: shell
156d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
157d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type shell -petscpartitioner_shell_random -petscpartitioner_view
158abe9303eSLisandro Dalcin 
159abe9303eSLisandro Dalcin   test:
1605130b6c4SLisandro Dalcin     suffix: gather
1615130b6c4SLisandro Dalcin     nsize: {{1 2 3}separate output}
1625130b6c4SLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type gather -petscpartitioner_view -petscpartitioner_view_graph
1635130b6c4SLisandro Dalcin 
1645130b6c4SLisandro Dalcin   test:
165abe9303eSLisandro Dalcin     requires: parmetis
166abe9303eSLisandro Dalcin     suffix: parmetis
167abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
168d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}} -vwgts {{false true}}
169d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -petscpartitioner_view -petscpartitioner_view_graph
170d7cc930eSLisandro Dalcin 
171d7cc930eSLisandro Dalcin   test:
172d7cc930eSLisandro Dalcin     requires: parmetis
173d7cc930eSLisandro Dalcin     suffix: parmetis_type
174d7cc930eSLisandro Dalcin     nsize: {{1 2}}
175d7cc930eSLisandro Dalcin     args: -petscpartitioner_type parmetis -part_view
176d7cc930eSLisandro Dalcin     args: -petscpartitioner_parmetis_type {{kway rb}separate output}
177d7cc930eSLisandro Dalcin     filter: grep "ParMetis type"
178abe9303eSLisandro Dalcin 
179abe9303eSLisandro Dalcin   test:
180abe9303eSLisandro Dalcin     requires: ptscotch
181abe9303eSLisandro Dalcin     suffix: ptscotch
182abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
183d7cc930eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -vwgts {{false true}}
184d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_view_graph
185abe9303eSLisandro Dalcin 
18656991f02SLisandro Dalcin   test:
187d7cc930eSLisandro Dalcin     requires: ptscotch
188d7cc930eSLisandro Dalcin     suffix: ptscotch_strategy
189d7cc930eSLisandro Dalcin     nsize: {{1 2}}
190d7cc930eSLisandro Dalcin     args: -petscpartitioner_type ptscotch -part_view
191d7cc930eSLisandro Dalcin     args: -petscpartitioner_ptscotch_strategy {{DEFAULT QUALITY SPEED BALANCE SAFETY SCALABILITY RECURSIVE REMAP}separate output}
192d7cc930eSLisandro Dalcin     filter: grep "partitioning strategy"
193d7cc930eSLisandro Dalcin 
194d7cc930eSLisandro Dalcin   test:
19556991f02SLisandro Dalcin     requires: chaco
19656991f02SLisandro Dalcin     suffix: chaco
19756991f02SLisandro Dalcin     nsize: {{1 2 3}separate output}
198d7cc930eSLisandro Dalcin     args: -nparts {{1}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
199d7cc930eSLisandro Dalcin 
200d7cc930eSLisandro Dalcin   test:
201d7cc930eSLisandro Dalcin     TODO: non reproducible (uses C stdlib rand())
202d7cc930eSLisandro Dalcin     requires: chaco
203d7cc930eSLisandro Dalcin     suffix: chaco
204d7cc930eSLisandro Dalcin     nsize: {{1 2 3}separate output}
205d7cc930eSLisandro Dalcin     args: -nparts {{2 3}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
20656991f02SLisandro Dalcin 
207abe9303eSLisandro Dalcin TEST*/
208