xref: /petsc/src/mat/impls/aij/mpi/mpb_aij.c (revision 7b6bb2c608b6fc6714ef38fda02c2dbb91c82665)
1 #include <../src/mat/impls/aij/mpi/mpiaij.h>
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "MatGetMultiProcBlock_MPIAIJ"
5 PetscErrorCode  MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, Mat* subMat)
6 {
7   PetscErrorCode ierr;
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;
13 
14   PetscFunctionBegin;
15   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&commsize);CHKERRQ(ierr);
16   ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr);
17 
18   /* create subMat object with the relavent layout */
19   ierr = MatCreate(subComm,subMat);CHKERRQ(ierr);
20   ierr = MatSetType(*subMat,MATMPIAIJ);CHKERRQ(ierr);
21   ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
22   /* need to setup rmap and cmap before Preallocation */
23   ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr);
24   ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr);
25   ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr);
26   ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr);
27 
28   /* create a map of comm_rank from subComm to comm */
29   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&commRank);CHKERRQ(ierr);
30   ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr);
31   ierr = PetscMalloc(subCommSize*sizeof(PetscMPIInt),&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 = PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr);
38   ierr = PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
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   }
49 
50   /* Now compute preallocation for the offdiag mat */
51   ierr = PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
52   ierr = PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
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   ierr = MatMPIAIJSetPreallocation(*(subMat),PETSC_NULL,PETSC_NULL,PETSC_NULL,nnz);CHKERRQ(ierr);
59 
60   /* reuse diag block with the new submat */
61   ierr = MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr);
62   ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
63   ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
64 
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         ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr);
73       }
74     }
75   }
76 
77   /* assemble the submat */
78   ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79   ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80 
81   /* deallocate temporary data */
82   ierr = PetscFree(commRankMap);CHKERRQ(ierr);
83   ierr = PetscFree(garrayCMap);CHKERRQ(ierr);
84   ierr = PetscFree(nnz);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87