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(PetscObjectComm((PetscObject)mat),&commsize);CHKERRQ(ierr); 16 ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr); 17 18 /* create subMat object with the relavent layout */ 19 if (scall == MAT_INITIAL_MATRIX) { 20 ierr = MatCreate(subComm,subMat);CHKERRQ(ierr); 21 ierr = MatSetType(*subMat,MATMPIAIJ);CHKERRQ(ierr); 22 ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 23 ierr = MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);CHKERRQ(ierr); 24 25 /* need to setup rmap and cmap before Preallocation */ 26 ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr); 27 ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr); 28 ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr); 29 ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr); 30 } 31 32 /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */ 33 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);CHKERRQ(ierr); 34 ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr); 35 ierr = PetscMalloc1(subCommSize,&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 = PetscCalloc1(aij->B->cmap->n,&garrayCMap);CHKERRQ(ierr); 42 for (i=0; i<aij->B->cmap->n; i++) { 43 col = aij->garray[i]; 44 for (subRank=0; subRank<subCommSize; subRank++) { 45 rank = commRankMap[subRank]; 46 if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) { 47 garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+1; 48 break; 49 } 50 } 51 } 52 53 if (scall == MAT_INITIAL_MATRIX) { 54 /* Now compute preallocation for the offdiag mat */ 55 ierr = PetscCalloc1(aij->B->rmap->n,&nnz);CHKERRQ(ierr); 56 for (i=0; i<aij->B->rmap->n; i++) { 57 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 58 if (garrayCMap[aijB->j[j]]) nnz[i]++; 59 } 60 } 61 ierr = MatMPIAIJSetPreallocation(*(subMat),0,NULL,0,nnz);CHKERRQ(ierr); 62 63 /* reuse diag block with the new submat */ 64 ierr = MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr); 65 66 ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A; 67 68 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 69 } else if (((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A) { 70 PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A; 71 72 ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr); 73 74 ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A; 75 76 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 77 } 78 79 /* Now traverse aij->B and insert values into subMat */ 80 for (i=0; i<aij->B->rmap->n; i++) { 81 newRow = (*subMat)->rmap->range[subCommRank] + i; 82 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 83 newCol = garrayCMap[aijB->j[j]]; 84 if (newCol) { 85 newCol--; /* remove the increment */ 86 ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr); 87 } 88 } 89 } 90 91 /* assemble the submat */ 92 ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93 ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94 95 /* deallocate temporary data */ 96 ierr = PetscFree(commRankMap);CHKERRQ(ierr); 97 ierr = PetscFree(garrayCMap);CHKERRQ(ierr); 98 if (scall == MAT_INITIAL_MATRIX) { 99 ierr = PetscFree(nnz);CHKERRQ(ierr); 100 } 101 PetscFunctionReturn(0); 102 } 103