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 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex183","Mat");PetscCall(ierr); 31c4762a1bSJed Brown m = 5; 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m","Local matrix size","MatSetSizes",m,&m,&flg)); 33c4762a1bSJed Brown total_subdomains = size-1; 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-total_subdomains","Number of submatrices where 0 < n < comm size","MatCreateSubMatricesMPI",total_subdomains,&total_subdomains,&flg)); 35c4762a1bSJed Brown permute_indices = PETSC_FALSE; 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-permute_indices","Whether to permute indices before breaking them into subdomains","ISCreateGeneral",permute_indices,&permute_indices,&flg)); 37c4762a1bSJed Brown hash = 7; 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hash","Permutation factor, which has to be relatively prime to M = size*m (total matrix size)","ISCreateGeneral",hash,&hash,&flg)); 39c4762a1bSJed Brown rep = 2; 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-rep","Number of times to carry out submatrix extractions; currently only 1 & 2 are supported",NULL,rep,&rep,&flg)); 419566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 42c4762a1bSJed Brown 43*08401ef6SPierre Jolivet PetscCheck(total_subdomains <= size,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of subdomains %" PetscInt_FMT " must not exceed comm size %d",total_subdomains,size); 442c71b3e2SJacob 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); 452c71b3e2SJacob 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. */ 499566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE)); 519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 529566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 539566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 549566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,NULL,&n)); 559566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 569566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,n,NULL)); 579566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,n,NULL,N-n,NULL)); 589566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A,bs,n/bs,NULL)); 599566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A,bs,n/bs,NULL,(N-n)/bs,NULL)); 609566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A,bs,n/bs,NULL)); 619566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A,bs,n/bs,NULL,(N-n)/bs,NULL)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N,&cols,N,&vals)); 649566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 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 } 709566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,N,cols,vals,INSERT_VALUES)); 71c4762a1bSJed Brown } 729566063dSJacob Faibussowitsch PetscCall(PetscFree2(cols,vals)); 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Initial matrix:\n")); 779566063dSJacob Faibussowitsch PetscCall(MatView(A,viewer)); 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; 869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD,color,rank,&subcomm)); 879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm,&subsize)); 889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(subcomm,&subrank)); 899566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 90c4762a1bSJed Brown nis = 1; 919566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(rend-rstart,&rowindices,rend-rstart,&colindices)); 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 } 1029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm,rend-rstart,rowindices,PETSC_COPY_VALUES,&rowis[0])); 1039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm,rend-rstart,colindices,PETSC_COPY_VALUES,&colis[0])); 1049566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0])); 1059566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0])); 1069566063dSJacob Faibussowitsch PetscCall(PetscFree2(rowindices,colindices)); 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 */ 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Subdomains")); 113c4762a1bSJed Brown if (permute_indices) { 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," (hash=%" PetscInt_FMT ")",hash)); 115c4762a1bSJed Brown } 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,":\n")); 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown nsubdomains = 1; 120c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1219566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)rowis,&gnsubdomains,gsubdomainnums)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm)); 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; 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)rowis[ss]),&subviewer)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Row IS %" PetscInt_FMT "\n",gs)); 1319566063dSJacob Faibussowitsch PetscCall(ISView(rowis[ss],subviewer)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(subviewer)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Col IS %" PetscInt_FMT "\n",gs)); 1349566063dSJacob Faibussowitsch PetscCall(ISView(colis[ss],subviewer)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)rowis[ss]),&subviewer)); 136c4762a1bSJed Brown ++s; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 1399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 140c4762a1bSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1429566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0])); 1439566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0])); 144c4762a1bSJed Brown nsubdomains = 1; 1459566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A,nsubdomains,rowis,colis,MAT_INITIAL_MATRIX,&submats)); 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 */ 1519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Submatrices (repetition 1):\n")); 152c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)submats,&gnsubdomains,gsubdomainnums)); 1549566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm)); 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; 1619566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer)); 1629566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss],subviewer)); 1639566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer)); 164c4762a1bSJed Brown ++s; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 168c4762a1bSJed Brown } 1699566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 170c4762a1bSJed Brown if (rep == 1) goto cleanup; 171c4762a1bSJed Brown nsubdomains = 1; 1729566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A,nsubdomains,rowis,colis,MAT_REUSE_MATRIX,&submats)); 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 */ 1789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Submatrices (repetition 2):\n")); 179c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1809566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD,1,(PetscObject*)submats,&gnsubdomains,gsubdomainnums)); 1819566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains,gsubdomainnums,gsubdomainperm)); 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; 1889566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer)); 1899566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss],subviewer)); 1909566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)submats[ss]),&subviewer)); 191c4762a1bSJed Brown ++s; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 1949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 195c4762a1bSJed Brown } 1969566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 197c4762a1bSJed Brown cleanup: 198c4762a1bSJed Brown for (k=0;k<nsubdomains;++k) { 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(submats+k)); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(submats)); 202c4762a1bSJed Brown for (k=0;k<nis;++k) { 2039566063dSJacob Faibussowitsch PetscCall(ISDestroy(rowis+k)); 2049566063dSJacob Faibussowitsch PetscCall(ISDestroy(colis+k)); 205c4762a1bSJed Brown } 2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm)); 2089566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 209b122ec5aSJacob Faibussowitsch return 0; 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