1 #include <../src/mat/impls/baij/mpi/mpibaij.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "MatGetMultiProcBlock_MPIBAIJ" 5 PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat) 6 { 7 PetscErrorCode ierr; 8 Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)mat->data; 9 Mat_SeqBAIJ *aijB = (Mat_SeqBAIJ*)aij->B->data; 10 PetscMPIInt commRank,subCommSize,subCommRank; 11 PetscMPIInt *commRankMap,subRank,rank,commsize; 12 PetscInt *garrayCMap,col,i,j,*nnz,newRow,newCol; 13 PetscInt bs=mat->rmap->bs; 14 PetscScalar vals[bs*bs]; 15 PetscInt newbRow[bs],newbCol[bs],k,k1,k2; 16 17 PetscFunctionBegin; 18 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);CHKERRQ(ierr); 19 ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr); 20 21 /* create subMat object with the relavent layout */ 22 if (scall == MAT_INITIAL_MATRIX) { 23 ierr = MatCreate(subComm,subMat);CHKERRQ(ierr); 24 ierr = MatSetType(*subMat,MATMPIBAIJ);CHKERRQ(ierr); 25 ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 26 ierr = MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);CHKERRQ(ierr); 27 28 /* need to setup rmap and cmap before Preallocation */ 29 ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr); 30 ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr); 31 ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr); 32 ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr); 33 } 34 35 /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */ 36 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);CHKERRQ(ierr); 37 ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr); 38 ierr = PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);CHKERRQ(ierr); 39 ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRQ(ierr); 40 41 /* Traverse garray and identify blocked column indices [of offdiag mat] that 42 should be discarded. For the ones not discarded, store the newCol+1 43 value in garrayCMap */ 44 ierr = PetscMalloc(aij->B->cmap->n/bs*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr); 45 ierr = PetscMemzero(garrayCMap,aij->B->cmap->n/bs*sizeof(PetscInt));CHKERRQ(ierr); 46 for (i=0; i<aij->B->cmap->n/bs; i++) { 47 col = aij->garray[i]; /* blocked column index */ 48 for (subRank=0; subRank<subCommSize; subRank++) { 49 rank = commRankMap[subRank]; 50 if ((col >= mat->cmap->range[rank]/bs) && (col < mat->cmap->range[rank+1]/bs)) { 51 garrayCMap[i] = ((*subMat)->cmap->range[subRank]/bs + col - mat->cmap->range[rank]/bs+1); 52 break; 53 } 54 } 55 } 56 57 if (scall == MAT_INITIAL_MATRIX) { 58 /* Now compute preallocation for the offdiag mat */ 59 ierr = PetscMalloc(aij->B->rmap->n/bs*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 60 ierr = PetscMemzero(nnz,aij->B->rmap->n/bs*sizeof(PetscInt));CHKERRQ(ierr); 61 for (i=0; i<aij->B->rmap->n/bs; i++) { 62 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 63 if (garrayCMap[aijB->j[j]]) nnz[i]++; 64 } 65 } 66 ierr = MatMPIBAIJSetPreallocation(*(subMat),bs,0,NULL,0,nnz);CHKERRQ(ierr); 67 68 /* reuse diag block with the new submat */ 69 ierr = MatDestroy(&((Mat_MPIBAIJ*)((*subMat)->data))->A);CHKERRQ(ierr); 70 71 ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A; 72 73 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 74 } else if (((Mat_MPIBAIJ*)(*subMat)->data)->A != aij->A) { 75 PetscObject obj = (PetscObject)((Mat_MPIBAIJ*)((*subMat)->data))->A; 76 77 ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr); 78 79 ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A; 80 81 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 82 } 83 84 /* Now traverse aij->B and insert values into subMat */ 85 for (i=0; i<aij->B->rmap->n/bs; i++) { 86 newRow = (*subMat)->rmap->range[subCommRank] + i*bs; 87 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 88 newCol = garrayCMap[aijB->j[j]]; 89 if (newCol) { 90 newCol--; /* remove the increment */ 91 newCol *= bs; 92 for (k=0; k<bs; k++) { 93 newbRow[k] = newRow + k; 94 newbCol[k] = newCol + k; 95 } 96 /* copy column-oriented aijB->a into row-oriented vals */ 97 k=0; 98 for (k1=0; k1<bs; k1++) { 99 for (k2=0; k2<bs; k2++) { 100 vals[k1+k2*bs] = *(aijB->a+j*bs*bs + k); k++; 101 } 102 } 103 ierr = MatSetValues(*subMat,bs,newbRow,bs,newbCol,vals,INSERT_VALUES);CHKERRQ(ierr); 104 } 105 } 106 } 107 108 /* assemble the submat */ 109 ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 110 ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 111 112 /* deallocate temporary data */ 113 ierr = PetscFree(commRankMap);CHKERRQ(ierr); 114 ierr = PetscFree(garrayCMap);CHKERRQ(ierr); 115 if (scall == MAT_INITIAL_MATRIX) { 116 ierr = PetscFree(nnz);CHKERRQ(ierr); 117 } 118 PetscFunctionReturn(0); 119 } 120