Actual source code: mpb_aij.c
1: #include <../src/mat/impls/aij/mpi/mpiaij.h>
5: PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, 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(((PetscObject)mat)->comm,&commsize);
16: MPI_Comm_size(subComm,&subCommSize);
18: /* create subMat object with the relavent layout */
19: MatCreate(subComm,subMat);
20: MatSetType(*subMat,MATMPIAIJ);
21: MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
22: /* need to setup rmap and cmap before Preallocation */
23: PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);
24: PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);
25: PetscLayoutSetUp((*subMat)->rmap);
26: PetscLayoutSetUp((*subMat)->cmap);
28: /* create a map of comm_rank from subComm to comm */
29: MPI_Comm_rank(((PetscObject)mat)->comm,&commRank);
30: MPI_Comm_rank(subComm,&subCommRank);
31: PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);
32: MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);
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: PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);
38: PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));
39: for (i=0; i<aij->B->cmap->n; i++) {
40: col = aij->garray[i];
41: for (subRank=0; subRank<subCommSize; subRank++) {
42: rank = commRankMap[subRank];
43: if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) {
44: garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+1;
45: break;
46: }
47: }
48: }
50: /* Now compute preallocation for the offdiag mat */
51: PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);
52: PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));
53: for (i=0; i<aij->B->rmap->n; i++) {
54: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
55: if (garrayCMap[aijB->j[j]]) nnz[i]++;
56: }
57: }
58: MatMPIAIJSetPreallocation(*(subMat),PETSC_NULL,PETSC_NULL,PETSC_NULL,nnz);
60: /* reuse diag block with the new submat */
61: MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);
62: ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
63: PetscObjectReference((PetscObject)aij->A);
65: /* Now traverse aij->B and insert values into subMat */
66: for (i=0; i<aij->B->rmap->n; i++) {
67: newRow = (*subMat)->rmap->range[subCommRank] + i;
68: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
69: newCol = garrayCMap[aijB->j[j]];
70: if (newCol) {
71: newCol--; /* remove the increment */
72: MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);
73: }
74: }
75: }
77: /* assemble the submat */
78: MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);
80:
81: /* deallocate temporary data */
82: PetscFree(commRankMap);
83: PetscFree(garrayCMap);
84: PetscFree(nnz);
85: return(0);
86: }