Actual source code: mpb_baij.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1:  #include <../src/mat/impls/baij/mpi/mpibaij.h>

  3: PetscErrorCode  MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
  4: {
  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;

 15:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);
 16:   MPI_Comm_size(subComm,&subCommSize);

 18:   /* create subMat object with the relavent layout */
 19:   if (scall == MAT_INITIAL_MATRIX) {
 20:     MatCreate(subComm,subMat);
 21:     MatSetType(*subMat,MATMPIBAIJ);
 22:     MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
 23:     MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);

 25:     /* need to setup rmap and cmap before Preallocation */
 26:     PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);
 27:     PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);
 28:     PetscLayoutSetUp((*subMat)->rmap);
 29:     PetscLayoutSetUp((*subMat)->cmap);
 30:   }

 32:   /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
 33:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);
 34:   MPI_Comm_rank(subComm,&subCommRank);
 35:   PetscMalloc1(subCommSize,&commRankMap);
 36:   MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);

 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:   PetscMalloc1(aij->B->cmap->n/bs,&garrayCMap);
 42:   PetscMemzero(garrayCMap,aij->B->cmap->n/bs*sizeof(PetscInt));
 43:   for (i=0; i<aij->B->cmap->n/bs; i++) {
 44:     col = aij->garray[i]; /* blocked column index */
 45:     for (subRank=0; subRank<subCommSize; subRank++) {
 46:       rank = commRankMap[subRank];
 47:       if ((col >= mat->cmap->range[rank]/bs) && (col < mat->cmap->range[rank+1]/bs)) {
 48:         garrayCMap[i] = (((*subMat)->cmap->range[subRank]- mat->cmap->range[rank])/bs + col + 1);
 49:         break;
 50:       }
 51:     }
 52:   }

 54:   if (scall == MAT_INITIAL_MATRIX) {
 55:     /* Now compute preallocation for the offdiag mat */
 56:     PetscMalloc1(aij->B->rmap->n/bs,&nnz);
 57:     PetscMemzero(nnz,aij->B->rmap->n/bs*sizeof(PetscInt));
 58:     for (i=0; i<aij->B->rmap->n/bs; i++) {
 59:       for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
 60:         if (garrayCMap[aijB->j[j]]) nnz[i]++;
 61:       }
 62:     }
 63:     MatMPIBAIJSetPreallocation(*(subMat),bs,0,NULL,0,nnz);

 65:     /* reuse diag block with the new submat */
 66:     MatDestroy(&((Mat_MPIBAIJ*)((*subMat)->data))->A);

 68:     ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;

 70:     PetscObjectReference((PetscObject)aij->A);
 71:   } else if (((Mat_MPIBAIJ*)(*subMat)->data)->A != aij->A) {
 72:     PetscObject obj = (PetscObject)((Mat_MPIBAIJ*)((*subMat)->data))->A;

 74:     PetscObjectReference((PetscObject)obj);

 76:     ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;

 78:     PetscObjectReference((PetscObject)aij->A);
 79:   }

 81:   /* Now traverse aij->B and insert values into subMat */
 82:   PetscMalloc3(bs,&newbRow,bs,&newbCol,bs*bs,&vals);
 83:   for (i=0; i<aij->B->rmap->n/bs; i++) {
 84:     newRow = (*subMat)->rmap->range[subCommRank] + i*bs;
 85:     for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
 86:       newCol = garrayCMap[aijB->j[j]];
 87:       if (newCol) {
 88:         newCol--; /* remove the increment */
 89:         newCol *= bs;
 90:         for (k=0; k<bs; k++) {
 91:           newbRow[k] = newRow + k;
 92:           newbCol[k] = newCol + k;
 93:         }
 94:         /* copy column-oriented aijB->a into row-oriented vals */
 95:         aijBvals = aijB->a + j*bs*bs;
 96:         for (k1=0; k1<bs; k1++) {
 97:           for (k=0; k<bs; k++) {
 98:             vals[k1+k*bs] = *aijBvals++;
 99:           }
100:         }
101:         MatSetValues(*subMat,bs,newbRow,bs,newbCol,vals,INSERT_VALUES);
102:       }
103:     }
104:   }
105:   MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
106:   MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);

108:   /* deallocate temporary data */
109:   PetscFree3(newbRow,newbCol,vals);
110:   PetscFree(commRankMap);
111:   PetscFree(garrayCMap);
112:   if (scall == MAT_INITIAL_MATRIX) {
113:     PetscFree(nnz);
114:   }
115:   return(0);
116: }