xref: /petsc/src/mat/impls/aij/mpi/mpb_aij.c (revision 5a7e1895ba1221176dc63d0bbec6b671c7ff4efe)
1 #include "../src/mat/impls/aij/mpi/mpiaij.h"
2 
3 /*
4 
5   This routine creates multiple [bjacobi] 'parallel submatrices' from
6   a given 'mat' object. Each submatrix can span multiple procs.
7 
8   The submatrix partition across processors is dicated by 'subComm' a
9   communicator obtained by com_split(comm). Note: the comm_split
10   is not restriced to be grouped with consequitive original ranks.
11 
12   Due the comm_split() usage, the parallel layout of the submatrices
13   map directly to the layout of the original matrix [wrt the local
14   row,col partitioning]. So the original 'DiagonalMat' naturally maps
15   into the 'DiagonalMat' of the subMat, hence it is used directly from
16   the subMat. However the offDiagMat looses some columns - and this is
17   reconstructed with MatSetValues()
18 
19  */
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "MatGetMultiProcBlock_MPIAIJ"
23 PetscErrorCode  MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, Mat* subMat)
24 {
25   PetscErrorCode ierr;
26   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
27   Mat_SeqAIJ*    aijB = (Mat_SeqAIJ*)aij->B->data;
28   PetscMPIInt    commRank,subCommSize,subCommRank;
29   PetscMPIInt    *commRankMap,subRank,rank,commsize;
30   PetscInt       *garrayCMap,col,i,j,*nnz,newRow,newCol;
31 
32   PetscFunctionBegin;
33   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&commsize);CHKERRQ(ierr);
34   ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr);
35   if (subCommSize > commsize) SETERRQ2(((PetscObject)mat)->comm,PETSC_ERR_ARG_OUTOFRANGE,"CommSize %D < SubCommZize %D",commsize,subCommSize);
36   if (commsize == 1){
37     ierr = MatDuplicate(mat,MAT_COPY_VALUES,subMat);CHKERRQ(ierr);
38     PetscFunctionReturn(0);
39   }
40 
41   /* create subMat object with the relavent layout */
42   ierr = MatCreate(subComm,subMat);CHKERRQ(ierr);
43   ierr = MatSetType(*subMat,MATMPIAIJ);CHKERRQ(ierr);
44   ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
45   /* need to setup rmap and cmap before Preallocation */
46   ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr);
47   ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr);
48   ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr);
49   ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr);
50 
51   /* create a map of comm_rank from subComm to comm */
52   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&commRank);CHKERRQ(ierr);
53   ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr);
54   ierr = PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);CHKERRQ(ierr);
55   ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRQ(ierr);
56 
57   /* Traverse garray and identify column indices [of offdiag mat] that
58    should be discarded. For the ones not discarded, store the newCol+1
59    value in garrayCMap */
60   ierr = PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr);
61   ierr = PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
62   for (i=0; i<aij->B->cmap->n; i++) {
63     col = aij->garray[i];
64     for (subRank=0; subRank<subCommSize; subRank++) {
65       rank = commRankMap[subRank];
66       if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) {
67         garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+1;
68         break;
69       }
70     }
71   }
72 
73   /* Now compute preallocation for the offdiag mat */
74   ierr = PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
75   ierr = PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
76   for (i=0; i<aij->B->rmap->n; i++) {
77     for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
78       if (garrayCMap[aijB->j[j]]) nnz[i]++;
79     }
80   }
81   ierr = MatMPIAIJSetPreallocation(*(subMat),PETSC_NULL,PETSC_NULL,PETSC_NULL,nnz);CHKERRQ(ierr);
82 
83   /* reuse diag block with the new submat */
84   ierr = MatDestroy(((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr);
85   ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
86   ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
87 
88   /* Now traverse aij->B and insert values into subMat */
89   for (i=0; i<aij->B->rmap->n; i++) {
90     newRow = (*subMat)->rmap->range[subCommRank] + i;
91     for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
92       newCol = garrayCMap[aijB->j[j]];
93       if (newCol) {
94         newCol--; /* remove the increment */
95         ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr);
96       }
97     }
98   }
99 
100   /* assemble the submat */
101   ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
102   ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103 
104   /* deallocate temporary data */
105   ierr = PetscFree(commRankMap);CHKERRQ(ierr);
106   ierr = PetscFree(garrayCMap);CHKERRQ(ierr);
107   ierr = PetscFree(nnz);CHKERRQ(ierr);
108   PetscFunctionReturn(0);
109 }
110