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