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