1c4762a1bSJed Brown static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc, char **args) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown Mat A; 8c4762a1bSJed Brown MatPartitioning part; 9c4762a1bSJed Brown IS is; 10c4762a1bSJed Brown PetscInt i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks; 11c4762a1bSJed Brown PetscMPIInt rank,size; 12c4762a1bSJed Brown 13*327415f7SBarry Smith PetscFunctionBeginUser; 149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 17c4762a1bSJed Brown 18c4762a1bSJed Brown nemptyranks = 10; 19c4762a1bSJed Brown nbigranks = 10; 209566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks)); 21c4762a1bSJed Brown 22d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL)); 25d0609cedSBarry Smith PetscOptionsEnd(); 26c4762a1bSJed Brown 27c4762a1bSJed Brown m = 1; 28c4762a1bSJed Brown for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0; 29c4762a1bSJed Brown for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5; 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE)); 339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 349566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,3,NULL)); 359566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,3,NULL,2,NULL)); 369566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A,1,3,NULL)); 379566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL)); 389566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A,1,2,NULL)); 399566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL)); 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 43c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 44c4762a1bSJed Brown const PetscInt cols[] = {(i+N-1)%N,i,(i+1)%N}; 45c4762a1bSJed Brown const PetscScalar vals[] = {1,1,1}; 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES)); 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 509566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD,&part)); 539566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,A)); 549566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 559566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&is)); 569566063dSJacob Faibussowitsch PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD)); 579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 589566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree2(emptyranks,bigranks)); 619566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 62b122ec5aSJacob Faibussowitsch return 0; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /*TEST 66c4762a1bSJed Brown 67c4762a1bSJed Brown test: 68c4762a1bSJed Brown nsize: 8 69c4762a1bSJed Brown args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average 70c4762a1bSJed Brown # cannot test with external package partitioners since they produce different results on different systems 71c4762a1bSJed Brown 72c4762a1bSJed Brown TEST*/ 73