1c4762a1bSJed Brown static char help[] = "Example of extracting an array of MPI submatrices from a given MPI matrix.\n" 2c4762a1bSJed Brown "This test can only be run in parallel.\n" 3c4762a1bSJed Brown "\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /*T 6c4762a1bSJed Brown Concepts: Mat^mat submatrix, parallel 7c4762a1bSJed Brown Processors: n 8c4762a1bSJed Brown T*/ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown int main(int argc, char **args) 13c4762a1bSJed Brown { 14c4762a1bSJed Brown Mat A,*submats; 15c4762a1bSJed Brown MPI_Comm subcomm; 16c4762a1bSJed Brown PetscMPIInt rank,size,subrank,subsize,color; 17c4762a1bSJed Brown PetscInt m,n,N,bs,rstart,rend,i,j,k,total_subdomains,hash,nsubdomains=1; 18c4762a1bSJed Brown PetscInt nis,*cols,gnsubdomains,gsubdomainnums[1],gsubdomainperm[1],s,gs; 19c4762a1bSJed Brown PetscInt *rowindices,*colindices,idx,rep; 20c4762a1bSJed Brown PetscScalar *vals; 21c4762a1bSJed Brown IS rowis[1],colis[1]; 22c4762a1bSJed Brown PetscViewer viewer; 23c4762a1bSJed Brown PetscBool permute_indices,flg; 24c4762a1bSJed Brown PetscErrorCode ierr; 25c4762a1bSJed Brown 26c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 27ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 28ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex183","Mat");CHKERRQ(ierr); 31c4762a1bSJed Brown m = 5; 32c4762a1bSJed Brown ierr = PetscOptionsInt("-m","Local matrix size","MatSetSizes",m,&m,&flg);CHKERRQ(ierr); 33c4762a1bSJed Brown total_subdomains = size-1; 34c4762a1bSJed Brown ierr = PetscOptionsInt("-total_subdomains","Number of submatrices where 0 < n < comm size","MatCreateSubMatricesMPI",total_subdomains,&total_subdomains,&flg);CHKERRQ(ierr); 35c4762a1bSJed Brown permute_indices = PETSC_FALSE; 36c4762a1bSJed Brown ierr = PetscOptionsBool("-permute_indices","Whether to permute indices before breaking them into subdomains","ISCreateGeneral",permute_indices,&permute_indices,&flg);CHKERRQ(ierr); 37c4762a1bSJed Brown hash = 7; 38c4762a1bSJed Brown ierr = PetscOptionsInt("-hash","Permutation factor, which has to be relatively prime to M = size*m (total matrix size)","ISCreateGeneral",hash,&hash,&flg);CHKERRQ(ierr); 39c4762a1bSJed Brown rep = 2; 40c4762a1bSJed Brown ierr = PetscOptionsInt("-rep","Number of times to carry out submatrix extractions; currently only 1 & 2 are supported",NULL,rep,&rep,&flg);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 42c4762a1bSJed Brown 43*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(total_subdomains > size,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of subdomains %" PetscInt_FMT " must not exceed comm size %d",total_subdomains,size); 44*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(total_subdomains < 1 || total_subdomains > size,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subdomains must be > 0 and <= %d (comm size), got total_subdomains = %" PetscInt_FMT,size,total_subdomains); 45*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(rep != 1 && rep != 2,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid number of test repetitions: %" PetscInt_FMT "; must be 1 or 2",rep); 46c4762a1bSJed Brown 47c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_WORLD; 48c4762a1bSJed Brown /* Create logically sparse, but effectively dense matrix for easy verification of submatrix extraction correctness. */ 49c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatGetLocalSize(A,NULL,&n);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,n,NULL);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,n,NULL,N-n,NULL);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,n/bs,NULL);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,n/bs,NULL,(N-n)/bs,NULL);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = MatSeqSBAIJSetPreallocation(A,bs,n/bs,NULL);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = MatMPISBAIJSetPreallocation(A,bs,n/bs,NULL,(N-n)/bs,NULL);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown ierr = PetscMalloc2(N,&cols,N,&vals);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 65c4762a1bSJed Brown for (j = 0; j < N; ++j) cols[j] = j; 66c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 67c4762a1bSJed Brown for (j=0;j<N;++j) { 68c4762a1bSJed Brown vals[j] = i*10000+j; 69c4762a1bSJed Brown } 70c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,N,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown ierr = PetscFree2(cols,vals);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Initial matrix:\n");CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = MatView(A,viewer);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Create subcomms and ISs so that each rank participates in one IS. 81c4762a1bSJed Brown The IS either coalesces adjacent rank indices (contiguous), 82c4762a1bSJed Brown or selects indices by scrambling them using a hash. 83c4762a1bSJed Brown */ 84c4762a1bSJed Brown k = size/total_subdomains + (size%total_subdomains>0); /* There are up to k ranks to a color */ 85c4762a1bSJed Brown color = rank/k; 86ffc4695bSBarry Smith ierr = MPI_Comm_split(PETSC_COMM_WORLD,color,rank,&subcomm);CHKERRMPI(ierr); 87ffc4695bSBarry Smith ierr = MPI_Comm_size(subcomm,&subsize);CHKERRMPI(ierr); 88ffc4695bSBarry Smith ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRMPI(ierr); 89c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 90c4762a1bSJed Brown nis = 1; 91c4762a1bSJed Brown ierr = PetscMalloc2(rend-rstart,&rowindices,rend-rstart,&colindices);CHKERRQ(ierr); 92c4762a1bSJed Brown 93c4762a1bSJed Brown for (j = rstart; j < rend; ++j) { 94c4762a1bSJed Brown if (permute_indices) { 95c4762a1bSJed Brown idx = (j*hash); 96c4762a1bSJed Brown } else { 97c4762a1bSJed Brown idx = j; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown rowindices[j-rstart] = idx%N; 100c4762a1bSJed Brown colindices[j-rstart] = (idx+m)%N; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown ierr = ISCreateGeneral(subcomm,rend-rstart,rowindices,PETSC_COPY_VALUES,&rowis[0]);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = ISCreateGeneral(subcomm,rend-rstart,colindices,PETSC_COPY_VALUES,&colis[0]);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = ISSort(rowis[0]);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = ISSort(colis[0]);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscFree2(rowindices,colindices);CHKERRQ(ierr); 107c4762a1bSJed Brown /* 108c4762a1bSJed Brown Now view the ISs. To avoid deadlock when viewing a list of objects on different subcomms, 109c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 110c4762a1bSJed Brown number to be viewed. 111c4762a1bSJed Brown */ 112c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Subdomains");CHKERRQ(ierr); 113c4762a1bSJed Brown if (permute_indices) { 1148fffc762SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(viewer," (hash=%" PetscInt_FMT ")",hash);CHKERRQ(ierr); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,":\n");CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown nsubdomains = 1; 120c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 121c4762a1bSJed Brown ierr = PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)rowis,&gnsubdomains,gsubdomainnums);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm);CHKERRQ(ierr); 123c4762a1bSJed Brown for (gs=0,s=0; gs < gnsubdomains;++gs) { 124c4762a1bSJed Brown if (s < nsubdomains) { 125c4762a1bSJed Brown PetscInt ss; 126c4762a1bSJed Brown ss = gsubdomainperm[s]; 127c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 128c4762a1bSJed Brown PetscViewer subviewer = NULL; 129c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)rowis[ss]),&subviewer);CHKERRQ(ierr); 1308fffc762SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(subviewer,"Row IS %" PetscInt_FMT "\n",gs);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = ISView(rowis[ss],subviewer);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = PetscViewerFlush(subviewer);CHKERRQ(ierr); 1338fffc762SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(subviewer,"Col IS %" PetscInt_FMT "\n",gs);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = ISView(colis[ss],subviewer);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)rowis[ss]),&subviewer);CHKERRQ(ierr); 136c4762a1bSJed Brown ++s; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 139ffc4695bSBarry Smith ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRMPI(ierr); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = ISSort(rowis[0]);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = ISSort(colis[0]);CHKERRQ(ierr); 144c4762a1bSJed Brown nsubdomains = 1; 145c4762a1bSJed Brown ierr = MatCreateSubMatricesMPI(A,nsubdomains,rowis,colis,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr); 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms, 148c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 149c4762a1bSJed Brown number to be viewed. 150c4762a1bSJed Brown */ 151c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Submatrices (repetition 1):\n");CHKERRQ(ierr); 152c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 153c4762a1bSJed Brown ierr = PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)submats,&gnsubdomains,gsubdomainnums);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm);CHKERRQ(ierr); 155c4762a1bSJed Brown for (gs=0,s=0; gs < gnsubdomains;++gs) { 156c4762a1bSJed Brown if (s < nsubdomains) { 157c4762a1bSJed Brown PetscInt ss; 158c4762a1bSJed Brown ss = gsubdomainperm[s]; 159c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 160c4762a1bSJed Brown PetscViewer subviewer = NULL; 161c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatView(submats[ss],subviewer);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer);CHKERRQ(ierr); 164c4762a1bSJed Brown ++s; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167ffc4695bSBarry Smith ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRMPI(ierr); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 170c4762a1bSJed Brown if (rep == 1) goto cleanup; 171c4762a1bSJed Brown nsubdomains = 1; 172c4762a1bSJed Brown ierr = MatCreateSubMatricesMPI(A,nsubdomains,rowis,colis,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms, 175c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 176c4762a1bSJed Brown number to be viewed. 177c4762a1bSJed Brown */ 178c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Submatrices (repetition 2):\n");CHKERRQ(ierr); 179c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 180c4762a1bSJed Brown ierr = PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)submats,&gnsubdomains,gsubdomainnums);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm);CHKERRQ(ierr); 182c4762a1bSJed Brown for (gs=0,s=0; gs < gnsubdomains;++gs) { 183c4762a1bSJed Brown if (s < nsubdomains) { 184c4762a1bSJed Brown PetscInt ss; 185c4762a1bSJed Brown ss = gsubdomainperm[s]; 186c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 187c4762a1bSJed Brown PetscViewer subviewer = NULL; 188c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = MatView(submats[ss],subviewer);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer);CHKERRQ(ierr); 191c4762a1bSJed Brown ++s; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 194ffc4695bSBarry Smith ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRMPI(ierr); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 197c4762a1bSJed Brown cleanup: 198c4762a1bSJed Brown for (k=0;k<nsubdomains;++k) { 199c4762a1bSJed Brown ierr = MatDestroy(submats+k);CHKERRQ(ierr); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown ierr = PetscFree(submats);CHKERRQ(ierr); 202c4762a1bSJed Brown for (k=0;k<nis;++k) { 203c4762a1bSJed Brown ierr = ISDestroy(rowis+k);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = ISDestroy(colis+k);CHKERRQ(ierr); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 207ffc4695bSBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRMPI(ierr); 208c4762a1bSJed Brown ierr = PetscFinalize(); 209c4762a1bSJed Brown return ierr; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown /*TEST 213c4762a1bSJed Brown 214c4762a1bSJed Brown test: 215c4762a1bSJed Brown nsize: 2 216c4762a1bSJed Brown args: -total_subdomains 1 217c4762a1bSJed Brown output_file: output/ex183_2_1.out 218c4762a1bSJed Brown 219c4762a1bSJed Brown test: 220c4762a1bSJed Brown suffix: 2 221c4762a1bSJed Brown nsize: 3 222c4762a1bSJed Brown args: -total_subdomains 2 223c4762a1bSJed Brown output_file: output/ex183_3_2.out 224c4762a1bSJed Brown 225c4762a1bSJed Brown test: 226c4762a1bSJed Brown suffix: 3 227c4762a1bSJed Brown nsize: 4 228c4762a1bSJed Brown args: -total_subdomains 2 229c4762a1bSJed Brown output_file: output/ex183_4_2.out 230c4762a1bSJed Brown 231c4762a1bSJed Brown test: 232c4762a1bSJed Brown suffix: 4 233c4762a1bSJed Brown nsize: 6 234c4762a1bSJed Brown args: -total_subdomains 2 235c4762a1bSJed Brown output_file: output/ex183_6_2.out 236c4762a1bSJed Brown 237c4762a1bSJed Brown TEST*/ 238