1 #include <../src/mat/impls/aij/mpi/mpiaij.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "MatGetMultiProcBlock_MPIAIJ" 5 /* 6 Developers Note: This is used directly by some preconditioners, hence is PETSC_EXTERN 7 */ 8 PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat) 9 { 10 PetscErrorCode ierr; 11 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 12 Mat_SeqAIJ *aijB = (Mat_SeqAIJ*)aij->B->data; 13 PetscMPIInt commRank,subCommSize,subCommRank; 14 PetscMPIInt *commRankMap,subRank,rank,commsize; 15 PetscInt *garrayCMap,col,i,j,*nnz,newRow,newCol; 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,MATMPIAIJ);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 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*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr); 45 ierr = PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 46 for (i=0; i<aij->B->cmap->n; i++) { 47 col = aij->garray[i]; 48 for (subRank=0; subRank<subCommSize; subRank++) { 49 rank = commRankMap[subRank]; 50 if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) { 51 garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+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*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 60 ierr = PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 61 for (i=0; i<aij->B->rmap->n; 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 = MatMPIAIJSetPreallocation(*(subMat),0,NULL,0,nnz);CHKERRQ(ierr); 67 68 /* reuse diag block with the new submat */ 69 ierr = MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr); 70 71 ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A; 72 73 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 74 } else if (((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A) { 75 PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A; 76 77 ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr); 78 79 ((Mat_MPIAIJ*)((*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; i++) { 86 newRow = (*subMat)->rmap->range[subCommRank] + i; 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 ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr); 92 } 93 } 94 } 95 96 /* assemble the submat */ 97 ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98 ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99 100 /* deallocate temporary data */ 101 ierr = PetscFree(commRankMap);CHKERRQ(ierr); 102 ierr = PetscFree(garrayCMap);CHKERRQ(ierr); 103 if (scall == MAT_INITIAL_MATRIX) { 104 ierr = PetscFree(nnz);CHKERRQ(ierr); 105 } 106 PetscFunctionReturn(0); 107 } 108