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