Actual source code: ex17.c
petsc-3.7.7 2017-09-25
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>
13: int main(int argc, char **args)
14: {
15: Mat A;
16: MatPartitioning part;
17: IS is;
18: PetscInt i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks;
19: PetscMPIInt rank,size;
20: PetscErrorCode ierr;
22: PetscInitialize(&argc,&args,(char*)0,help);
23: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
24: MPI_Comm_size(PETSC_COMM_WORLD,&size);
26: nemptyranks = 10;
27: nbigranks = 10;
28: PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks);
30: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL);
31: PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL);
32: PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL);
33: PetscOptionsEnd();
35: m = 1;
36: for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0;
37: for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5;
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
41: MatSetFromOptions(A);
42: MatSeqAIJSetPreallocation(A,3,NULL);
43: MatMPIAIJSetPreallocation(A,3,NULL,2,NULL);
44: MatSeqBAIJSetPreallocation(A,1,3,NULL);
45: MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL);
46: MatSeqSBAIJSetPreallocation(A,1,2,NULL);
47: MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL);
49: MatGetSize(A,NULL,&N);
50: MatGetOwnershipRange(A,&rstart,&rend);
51: for (i=rstart; i<rend; i++) {
52: const PetscInt cols[] = {(i+N-1)%N,i,(i+1)%N};
53: const PetscScalar vals[] = {1,1,1};
54: MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES);
55: }
56: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
58: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
60: MatPartitioningCreate(PETSC_COMM_WORLD,&part);
61: MatPartitioningSetAdjacency(part,A);
62: MatPartitioningSetFromOptions(part);
63: MatPartitioningApply(part,&is);
64: ISView(is,PETSC_VIEWER_STDOUT_WORLD);
65: ISDestroy(&is);
66: MatPartitioningDestroy(&part);
67: MatDestroy(&A);
68: PetscFree2(emptyranks,bigranks);
69: PetscFinalize();
70: return 0;
71: }