xref: /petsc/src/mat/impls/aij/mpi/mpb_aij.c (revision f692024ea6ceda132bc9ff30ca7a31e85cfbbcb2)
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, MatReuse scall,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   if (scall == MAT_INITIAL_MATRIX){
20     ierr = MatCreate(subComm,subMat);CHKERRQ(ierr);
21     ierr = MatSetType(*subMat,MATMPIAIJ);CHKERRQ(ierr);
22     ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
23 
24     /* need to setup rmap and cmap before Preallocation */
25     ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr);
26     ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr);
27     ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr);
28     ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr);
29   }
30 
31   /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
32   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&commRank);CHKERRQ(ierr);
33   ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr);
34   ierr = PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);CHKERRQ(ierr);
35   ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRQ(ierr);
36 
37   /* Traverse garray and identify column indices [of offdiag mat] that
38    should be discarded. For the ones not discarded, store the newCol+1
39    value in garrayCMap */
40   ierr = PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);CHKERRQ(ierr);
41   ierr = PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
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   }
52 
53   if (scall == MAT_INITIAL_MATRIX){
54     /* Now compute preallocation for the offdiag mat */
55     ierr = PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
56     ierr = PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
57     for (i=0; i<aij->B->rmap->n; i++) {
58       for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
59         if (garrayCMap[aijB->j[j]]) nnz[i]++;
60       }
61     }
62     ierr = MatMPIAIJSetPreallocation(*(subMat),0,PETSC_NULL,0,nnz);CHKERRQ(ierr);
63 
64     /* reuse diag block with the new submat */
65     ierr = MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);CHKERRQ(ierr);
66     ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
67     ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
68   } else if ( ((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A  ){
69     PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A;
70     ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr);
71     ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
72     ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
73   }
74 
75   /* Now traverse aij->B and insert values into subMat */
76   for (i=0; i<aij->B->rmap->n; i++) {
77     newRow = (*subMat)->rmap->range[subCommRank] + i;
78     for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
79       newCol = garrayCMap[aijB->j[j]];
80       if (newCol) {
81         newCol--; /* remove the increment */
82         ierr = MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);CHKERRQ(ierr);
83       }
84     }
85   }
86 
87   /* assemble the submat */
88   ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89   ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90 
91   /* deallocate temporary data */
92   ierr = PetscFree(commRankMap);CHKERRQ(ierr);
93   ierr = PetscFree(garrayCMap);CHKERRQ(ierr);
94   if (scall == MAT_INITIAL_MATRIX){
95     ierr = PetscFree(nnz);CHKERRQ(ierr);
96   }
97   PetscFunctionReturn(0);
98 }
99