xref: /petsc/src/dm/partitioner/tests/ex33.c (revision 56991f02f941884f2dfd92c62b471da0713cee35)
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;
9abe9303eSLisandro Dalcin   PetscSection     partSection;
10abe9303eSLisandro Dalcin   IS               partition,is;
11abe9303eSLisandro Dalcin   PetscMPIInt      size,rank;
12abe9303eSLisandro Dalcin   PetscInt         npar;
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};
16*56991f02SLisandro Dalcin   PetscBool        sequential;
17abe9303eSLisandro Dalcin 
18abe9303eSLisandro Dalcin   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
19abe9303eSLisandro Dalcin   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20abe9303eSLisandro Dalcin   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
21abe9303eSLisandro Dalcin   npar = size;
22abe9303eSLisandro Dalcin   ierr = PetscOptionsGetInt(NULL,NULL,"-nparts",&npar,NULL);CHKERRQ(ierr);
23abe9303eSLisandro Dalcin 
24abe9303eSLisandro Dalcin   /* create PetscPartitioner */
25abe9303eSLisandro Dalcin   ierr = PetscSectionCreate(PETSC_COMM_WORLD,&partSection);CHKERRQ(ierr);
26abe9303eSLisandro Dalcin   ierr = PetscPartitionerCreate(PETSC_COMM_WORLD,&p);CHKERRQ(ierr);
27abe9303eSLisandro Dalcin   ierr = PetscPartitionerSetType(p,PETSCPARTITIONERSIMPLE);CHKERRQ(ierr);
28abe9303eSLisandro Dalcin   ierr = PetscPartitionerSetFromOptions(p);CHKERRQ(ierr);
29abe9303eSLisandro Dalcin 
30abe9303eSLisandro Dalcin   /* test partitioning an empty graph */
31abe9303eSLisandro Dalcin   ierr = PetscPartitionerPartition(p,npar,0,NULL,NULL,NULL,NULL,partSection,&partition);CHKERRQ(ierr);
32abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"NULL SECTION");
33abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
34abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
35abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"NULL PARTITION");
36abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
37abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
38abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
39abe9303eSLisandro Dalcin 
40abe9303eSLisandro Dalcin   /* test partitioning a graph on one process only (not master) */
41abe9303eSLisandro Dalcin   if (rank == size - 1) {
42abe9303eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,npar,nv,vv,vadj,NULL,NULL,partSection,&partition);CHKERRQ(ierr);
43abe9303eSLisandro Dalcin   } else {
44abe9303eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,npar,0,NULL,NULL,NULL,NULL,partSection,&partition);CHKERRQ(ierr);
45abe9303eSLisandro Dalcin   }
46abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"SEQ SECTION");
47abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
48abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
49abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"SEQ PARTITION");
50abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
51abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
52abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
53abe9303eSLisandro Dalcin 
54*56991f02SLisandro Dalcin   ierr = PetscObjectTypeCompareAny((PetscObject)p,&sequential,PETSCPARTITIONERCHACO,NULL);CHKERRQ(ierr);
55*56991f02SLisandro Dalcin   if (sequential) goto finally;
56*56991f02SLisandro Dalcin 
57abe9303eSLisandro Dalcin   /* test partitioning a graph on a subset of the processess only */
58abe9303eSLisandro Dalcin   if (rank%2) {
59abe9303eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,npar,0,NULL,NULL,NULL,NULL,partSection,&partition);CHKERRQ(ierr);
60abe9303eSLisandro Dalcin   } else {
61abe9303eSLisandro Dalcin     PetscInt i,totv = nv*((size+1)/2),*pvadj;
62abe9303eSLisandro Dalcin 
63abe9303eSLisandro Dalcin     ierr = PetscMalloc1(2*nv,&pvadj);CHKERRQ(ierr);
64abe9303eSLisandro Dalcin     for (i = 0; i < nv; i++) {
65abe9303eSLisandro Dalcin       pvadj[2*i]   = (nv*(rank/2) + totv + i - 1)%totv;
66abe9303eSLisandro Dalcin       pvadj[2*i+1] = (nv*(rank/2) + totv + i + 1)%totv;
67abe9303eSLisandro Dalcin     }
68abe9303eSLisandro Dalcin     ierr = PetscPartitionerPartition(p,npar,nv,vv,pvadj,NULL,NULL,partSection,&partition);CHKERRQ(ierr);
69abe9303eSLisandro Dalcin     ierr = PetscFree(pvadj);CHKERRQ(ierr);
70abe9303eSLisandro Dalcin   }
71abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)partSection,"PARVOID SECTION");
72abe9303eSLisandro Dalcin   ierr = PetscSectionView(partSection,NULL);CHKERRQ(ierr);
73abe9303eSLisandro Dalcin   ierr = ISOnComm(partition,PETSC_COMM_WORLD,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
74abe9303eSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)is,"PARVOID PARTITION");
75abe9303eSLisandro Dalcin   ierr = ISView(is,NULL);CHKERRQ(ierr);
76abe9303eSLisandro Dalcin   ierr = ISDestroy(&is);CHKERRQ(ierr);
77abe9303eSLisandro Dalcin   ierr = ISDestroy(&partition);CHKERRQ(ierr);
78abe9303eSLisandro Dalcin 
79*56991f02SLisandro Dalcin finally:
80abe9303eSLisandro Dalcin   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
81abe9303eSLisandro Dalcin   ierr = PetscPartitionerDestroy(&p);CHKERRQ(ierr);
82abe9303eSLisandro Dalcin   ierr = PetscFinalize();
83abe9303eSLisandro Dalcin   return ierr;
84abe9303eSLisandro Dalcin }
85abe9303eSLisandro Dalcin 
86abe9303eSLisandro Dalcin /*TEST
87abe9303eSLisandro Dalcin 
88abe9303eSLisandro Dalcin   test:
89abe9303eSLisandro Dalcin     suffix: simple
90abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
91abe9303eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type simple -petscpartitioner_view -petscpartitioner_view_graph
92abe9303eSLisandro Dalcin 
93abe9303eSLisandro Dalcin   test:
94abe9303eSLisandro Dalcin     requires: parmetis
95abe9303eSLisandro Dalcin     suffix: parmetis
96abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
97abe9303eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type parmetis -petscpartitioner_view -petscpartitioner_view_graph
98abe9303eSLisandro Dalcin 
99abe9303eSLisandro Dalcin   test:
100abe9303eSLisandro Dalcin     requires: ptscotch
101abe9303eSLisandro Dalcin     suffix: ptscotch
102abe9303eSLisandro Dalcin     nsize: {{1 2 3}separate output}
103abe9303eSLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_view_graph
104abe9303eSLisandro Dalcin 
105*56991f02SLisandro Dalcin   test:
106*56991f02SLisandro Dalcin     TODO: broken
107*56991f02SLisandro Dalcin     requires: chaco
108*56991f02SLisandro Dalcin     suffix: chaco
109*56991f02SLisandro Dalcin     nsize: {{1 2 3}separate output}
110*56991f02SLisandro Dalcin     args: -nparts {{1 2 3}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph
111*56991f02SLisandro Dalcin 
112abe9303eSLisandro Dalcin TEST*/
113