Actual source code: mpb_aij.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <../src/mat/impls/aij/mpi/mpiaij.h>

  5: PetscErrorCode  MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
  6: {
  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;

 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,MATMPIAIJ);
 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 column indices [of offdiag mat] that
 39:    should be discarded. For the ones not discarded, store the newCol+1
 40:    value in garrayCMap */
 41:   PetscCalloc1(aij->B->cmap->n,&garrayCMap);
 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:   }

 53:   if (scall == MAT_INITIAL_MATRIX) {
 54:     /* Now compute preallocation for the offdiag mat */
 55:     PetscCalloc1(aij->B->rmap->n,&nnz);
 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:     MatMPIAIJSetPreallocation(*(subMat),0,NULL,0,nnz);

 63:     /* reuse diag block with the new submat */
 64:     MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);

 66:     ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;

 68:     PetscObjectReference((PetscObject)aij->A);
 69:   } else if (((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A) {
 70:     PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A;

 72:     PetscObjectReference((PetscObject)obj);

 74:     ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;

 76:     PetscObjectReference((PetscObject)aij->A);
 77:   }

 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:         MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);
 87:       }
 88:     }
 89:   }

 91:   /* assemble the submat */
 92:   MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);

 95:   /* deallocate temporary data */
 96:   PetscFree(commRankMap);
 97:   PetscFree(garrayCMap);
 98:   if (scall == MAT_INITIAL_MATRIX) {
 99:     PetscFree(nnz);
100:   }
101:   return(0);
102: }