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: #include <petscmat.h>
5: int main(int argc, char **args)
6: {
7: Mat A;
8: MatPartitioning part;
9: IS is;
10: PetscInt i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks;
11: PetscMPIInt rank,size;
12: PetscErrorCode ierr;
14: PetscInitialize(&argc,&args,(char*)0,help);
15: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
16: MPI_Comm_size(PETSC_COMM_WORLD,&size);
18: nemptyranks = 10;
19: nbigranks = 10;
20: PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks);
22: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL);
23: PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL);
24: PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL);
25: PetscOptionsEnd();
27: m = 1;
28: for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0;
29: for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5;
31: MatCreate(PETSC_COMM_WORLD,&A);
32: MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
33: MatSetFromOptions(A);
34: MatSeqAIJSetPreallocation(A,3,NULL);
35: MatMPIAIJSetPreallocation(A,3,NULL,2,NULL);
36: MatSeqBAIJSetPreallocation(A,1,3,NULL);
37: MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL);
38: MatSeqSBAIJSetPreallocation(A,1,2,NULL);
39: MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL);
41: MatGetSize(A,NULL,&N);
42: MatGetOwnershipRange(A,&rstart,&rend);
43: for (i=rstart; i<rend; i++) {
44: const PetscInt cols[] = {(i+N-1)%N,i,(i+1)%N};
45: const PetscScalar vals[] = {1,1,1};
46: MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES);
47: }
48: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
49: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
50: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
52: MatPartitioningCreate(PETSC_COMM_WORLD,&part);
53: MatPartitioningSetAdjacency(part,A);
54: MatPartitioningSetFromOptions(part);
55: MatPartitioningApply(part,&is);
56: ISView(is,PETSC_VIEWER_STDOUT_WORLD);
57: ISDestroy(&is);
58: MatPartitioningDestroy(&part);
59: MatDestroy(&A);
60: PetscFree2(emptyranks,bigranks);
61: PetscFinalize();
62: return 0;
63: }
65: /*TEST
67: test:
68: nsize: 8
69: args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
70: # cannot test with external package partitioners since they produce different results on different systems
72: TEST*/