Actual source code: ex17.c
1: static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n";
3: /*T
4: Concepts: Mat^mat partitioning
5: Concepts: Mat^image segmentation
6: Processors: n
7: T*/
9: #include <petscmat.h>
11: int main(int argc, char **args)
12: {
13: Mat A;
14: MatPartitioning part;
15: IS is;
16: PetscInt i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks;
17: PetscMPIInt rank,size;
18: PetscErrorCode ierr;
20: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: nemptyranks = 10;
25: nbigranks = 10;
26: PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks);
28: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL);
29: PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL);
30: PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL);
31: PetscOptionsEnd();
33: m = 1;
34: for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0;
35: for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5;
37: MatCreate(PETSC_COMM_WORLD,&A);
38: MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
39: MatSetFromOptions(A);
40: MatSeqAIJSetPreallocation(A,3,NULL);
41: MatMPIAIJSetPreallocation(A,3,NULL,2,NULL);
42: MatSeqBAIJSetPreallocation(A,1,3,NULL);
43: MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL);
44: MatSeqSBAIJSetPreallocation(A,1,2,NULL);
45: MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL);
47: MatGetSize(A,NULL,&N);
48: MatGetOwnershipRange(A,&rstart,&rend);
49: for (i=rstart; i<rend; i++) {
50: const PetscInt cols[] = {(i+N-1)%N,i,(i+1)%N};
51: const PetscScalar vals[] = {1,1,1};
52: MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES);
53: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
58: MatPartitioningCreate(PETSC_COMM_WORLD,&part);
59: MatPartitioningSetAdjacency(part,A);
60: MatPartitioningSetFromOptions(part);
61: MatPartitioningApply(part,&is);
62: ISView(is,PETSC_VIEWER_STDOUT_WORLD);
63: ISDestroy(&is);
64: MatPartitioningDestroy(&part);
65: MatDestroy(&A);
66: PetscFree2(emptyranks,bigranks);
67: PetscFinalize();
68: return ierr;
69: }
71: /*TEST
73: test:
74: nsize: 8
75: args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
76: # cannot test with external package partitioners since they produce different results on different systems
78: TEST*/