1 #include <../src/mat/impls/aij/mpi/mpiaij.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "MatGetMultiProcBlock_MPIAIJ" 5 PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat* subMat) 6 { 7 PetscErrorCode ierr; 8 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 9 Mat_SeqAIJ* aijB = (Mat_SeqAIJ*)aij->B->data; 10 PetscMPIInt commRank,subCommSize,subCommRank; 11 PetscMPIInt *commRankMap,subRank,rank,commsize; 12 PetscInt *garrayCMap,col,i,j,*nnz,newRow,newCol; 13 14 PetscFunctionBegin; 15 ierr = MPI_Comm_size(((PetscObject)mat)->comm,&commsize);CHKERRQ(ierr); 16 ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr); 17 18 if (scall == MAT_REUSE_MATRIX && subMat){ 19 ierr = MatDestroy(subMat);CHKERRQ(ierr); /* fix this case later! */ 20 } 21 /* create subMat object with the relavent layout */ 22 ierr = MatCreate(subComm,subMat);CHKERRQ(ierr); 23 ierr = MatSetType(*subMat,MATMPIAIJ);CHKERRQ(ierr); 24 ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 25 26 /* need to setup rmap and cmap before Preallocation */ 27 ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr); 28 ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr); 29 ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr); 30 ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr); 31 32 /* create a map of comm_rank from subComm to comm */ 33 ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&commRank);CHKERRQ(ierr); 34 ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr); 35 ierr = PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);CHKERRQ(ierr); 36 ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRQ(ierr); 37 38 /* Traverse garray and identify column indices [of offdiag mat] that 39 should be discarded. For the ones not discarded, store the newCol+1 40 value in garrayCMap */ 41 ierr = PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr); 42 ierr = PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 43 for (i=0; i<aij->B->cmap->n; i++) { 44 col = aij->garray[i]; 45 for (subRank=0; subRank<subCommSize; subRank++) { 46 rank = commRankMap[subRank]; 47 if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) { 48 garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+1; 49 break; 50 } 51 } 52 } 53 54 /* Now compute preallocation for the offdiag mat */ 55 ierr = PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 56 ierr = PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 57 for (i=0; i<aij->B->rmap->n; i++) { 58 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 59 if (garrayCMap[aijB->j[j]]) nnz[i]++; 60 } 61 } 62 ierr = MatMPIAIJSetPreallocation(*(subMat),0,PETSC_NULL,0,nnz);CHKERRQ(ierr); 63 64 65 /* reuse diag block with the new submat */ 66 ierr = MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr); 67 ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A; 68 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 69 70 /* Now traverse aij->B and insert values into subMat */ 71 for (i=0; i<aij->B->rmap->n; i++) { 72 newRow = (*subMat)->rmap->range[subCommRank] + i; 73 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 74 newCol = garrayCMap[aijB->j[j]]; 75 if (newCol) { 76 newCol--; /* remove the increment */ 77 ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr); 78 } 79 } 80 } 81 82 /* assemble the submat */ 83 ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 84 ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 85 86 /* deallocate temporary data */ 87 ierr = PetscFree(commRankMap);CHKERRQ(ierr); 88 ierr = PetscFree(garrayCMap);CHKERRQ(ierr); 89 ierr = PetscFree(nnz);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92