1 #include <../src/mat/impls/baij/mpi/mpibaij.h> 2 3 PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat) 4 { 5 PetscErrorCode ierr; 6 Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)mat->data; 7 Mat_SeqBAIJ *aijB = (Mat_SeqBAIJ*)aij->B->data; 8 PetscMPIInt commRank,subCommSize,subCommRank; 9 PetscMPIInt *commRankMap,subRank,rank,commsize; 10 PetscInt *garrayCMap,col,i,j,*nnz,newRow,newCol,*newbRow,*newbCol,k,k1; 11 PetscInt bs=mat->rmap->bs; 12 PetscScalar *vals,*aijBvals; 13 14 PetscFunctionBegin; 15 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);CHKERRMPI(ierr); 16 ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRMPI(ierr); 17 18 /* create subMat object with the relevant layout */ 19 if (scall == MAT_INITIAL_MATRIX) { 20 ierr = MatCreate(subComm,subMat);CHKERRQ(ierr); 21 ierr = MatSetType(*subMat,MATMPIBAIJ);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);CHKERRMPI(ierr); 34 ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRMPI(ierr); 35 ierr = PetscMalloc1(subCommSize,&commRankMap);CHKERRQ(ierr); 36 ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRMPI(ierr); 37 38 /* Traverse garray and identify blocked 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/bs,&garrayCMap);CHKERRQ(ierr); 42 for (i=0; i<aij->B->cmap->n/bs; i++) { 43 col = aij->garray[i]; /* blocked column index */ 44 for (subRank=0; subRank<subCommSize; subRank++) { 45 rank = commRankMap[subRank]; 46 if ((col >= mat->cmap->range[rank]/bs) && (col < mat->cmap->range[rank+1]/bs)) { 47 garrayCMap[i] = (((*subMat)->cmap->range[subRank]- mat->cmap->range[rank])/bs + col + 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/bs,&nnz);CHKERRQ(ierr); 56 for (i=0; i<aij->B->rmap->n/bs; 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 = MatMPIBAIJSetPreallocation(*(subMat),bs,0,NULL,0,nnz);CHKERRQ(ierr); 62 63 /* reuse diag block with the new submat */ 64 ierr = MatDestroy(&((Mat_MPIBAIJ*)((*subMat)->data))->A);CHKERRQ(ierr); 65 66 ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A; 67 68 ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr); 69 } else if (((Mat_MPIBAIJ*)(*subMat)->data)->A != aij->A) { 70 PetscObject obj = (PetscObject)((Mat_MPIBAIJ*)((*subMat)->data))->A; 71 72 ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr); 73 74 ((Mat_MPIBAIJ*)((*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 ierr = PetscMalloc3(bs,&newbRow,bs,&newbCol,bs*bs,&vals);CHKERRQ(ierr); 81 for (i=0; i<aij->B->rmap->n/bs; i++) { 82 newRow = (*subMat)->rmap->range[subCommRank] + i*bs; 83 for (j=aijB->i[i]; j<aijB->i[i+1]; j++) { 84 newCol = garrayCMap[aijB->j[j]]; 85 if (newCol) { 86 newCol--; /* remove the increment */ 87 newCol *= bs; 88 for (k=0; k<bs; k++) { 89 newbRow[k] = newRow + k; 90 newbCol[k] = newCol + k; 91 } 92 /* copy column-oriented aijB->a into row-oriented vals */ 93 aijBvals = aijB->a + j*bs*bs; 94 for (k1=0; k1<bs; k1++) { 95 for (k=0; k<bs; k++) { 96 vals[k1+k*bs] = *aijBvals++; 97 } 98 } 99 ierr = MatSetValues(*subMat,bs,newbRow,bs,newbCol,vals,INSERT_VALUES);CHKERRQ(ierr); 100 } 101 } 102 } 103 ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 104 ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105 106 /* deallocate temporary data */ 107 ierr = PetscFree3(newbRow,newbCol,vals);CHKERRQ(ierr); 108 ierr = PetscFree(commRankMap);CHKERRQ(ierr); 109 ierr = PetscFree(garrayCMap);CHKERRQ(ierr); 110 if (scall == MAT_INITIAL_MATRIX) { 111 ierr = PetscFree(nnz);CHKERRQ(ierr); 112 } 113 PetscFunctionReturn(0); 114 } 115