xref: /petsc/src/mat/impls/baij/mpi/mpb_baij.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <../src/mat/impls/baij/mpi/mpibaij.h>
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "MatGetMultiProcBlock_MPIBAIJ"
5 PetscErrorCode  MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
6 {
7   PetscErrorCode ierr;
8   Mat_MPIBAIJ    *aij  = (Mat_MPIBAIJ*)mat->data;
9   Mat_SeqBAIJ    *aijB = (Mat_SeqBAIJ*)aij->B->data;
10   PetscMPIInt    commRank,subCommSize,subCommRank;
11   PetscMPIInt    *commRankMap,subRank,rank,commsize;
12   PetscInt       *garrayCMap,col,i,j,*nnz,newRow,newCol,*newbRow,*newbCol,k,k1;
13   PetscInt       bs=mat->rmap->bs;
14   PetscScalar    *vals,*aijBvals;
15 
16   PetscFunctionBegin;
17   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);CHKERRQ(ierr);
18   ierr = MPI_Comm_size(subComm,&subCommSize);CHKERRQ(ierr);
19 
20   /* create subMat object with the relavent layout */
21   if (scall == MAT_INITIAL_MATRIX) {
22     ierr = MatCreate(subComm,subMat);CHKERRQ(ierr);
23     ierr = MatSetType(*subMat,MATMPIBAIJ);CHKERRQ(ierr);
24     ierr = MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
25     ierr = MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);CHKERRQ(ierr);
26 
27     /* need to setup rmap and cmap before Preallocation */
28     ierr = PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);CHKERRQ(ierr);
29     ierr = PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);CHKERRQ(ierr);
30     ierr = PetscLayoutSetUp((*subMat)->rmap);CHKERRQ(ierr);
31     ierr = PetscLayoutSetUp((*subMat)->cmap);CHKERRQ(ierr);
32   }
33 
34   /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
35   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);CHKERRQ(ierr);
36   ierr = MPI_Comm_rank(subComm,&subCommRank);CHKERRQ(ierr);
37   ierr = PetscMalloc1(subCommSize,&commRankMap);CHKERRQ(ierr);
38   ierr = MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);CHKERRQ(ierr);
39 
40   /* Traverse garray and identify blocked column indices [of offdiag mat] that
41    should be discarded. For the ones not discarded, store the newCol+1
42    value in garrayCMap */
43   ierr = PetscMalloc1(aij->B->cmap->n/bs,&garrayCMap);CHKERRQ(ierr);
44   ierr = PetscMemzero(garrayCMap,aij->B->cmap->n/bs*sizeof(PetscInt));CHKERRQ(ierr);
45   for (i=0; i<aij->B->cmap->n/bs; i++) {
46     col = aij->garray[i]; /* blocked column index */
47     for (subRank=0; subRank<subCommSize; subRank++) {
48       rank = commRankMap[subRank];
49       if ((col >= mat->cmap->range[rank]/bs) && (col < mat->cmap->range[rank+1]/bs)) {
50         garrayCMap[i] = (((*subMat)->cmap->range[subRank]- mat->cmap->range[rank])/bs + col + 1);
51         break;
52       }
53     }
54   }
55 
56   if (scall == MAT_INITIAL_MATRIX) {
57     /* Now compute preallocation for the offdiag mat */
58     ierr = PetscMalloc1(aij->B->rmap->n/bs,&nnz);CHKERRQ(ierr);
59     ierr = PetscMemzero(nnz,aij->B->rmap->n/bs*sizeof(PetscInt));CHKERRQ(ierr);
60     for (i=0; i<aij->B->rmap->n/bs; i++) {
61       for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
62         if (garrayCMap[aijB->j[j]]) nnz[i]++;
63       }
64     }
65     ierr = MatMPIBAIJSetPreallocation(*(subMat),bs,0,NULL,0,nnz);CHKERRQ(ierr);
66 
67     /* reuse diag block with the new submat */
68     ierr = MatDestroy(&((Mat_MPIBAIJ*)((*subMat)->data))->A);CHKERRQ(ierr);
69 
70     ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;
71 
72     ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
73   } else if (((Mat_MPIBAIJ*)(*subMat)->data)->A != aij->A) {
74     PetscObject obj = (PetscObject)((Mat_MPIBAIJ*)((*subMat)->data))->A;
75 
76     ierr = PetscObjectReference((PetscObject)obj);CHKERRQ(ierr);
77 
78     ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;
79 
80     ierr = PetscObjectReference((PetscObject)aij->A);CHKERRQ(ierr);
81   }
82 
83   /* Now traverse aij->B and insert values into subMat */
84   ierr = PetscMalloc3(bs,&newbRow,bs,&newbCol,bs*bs,&vals);CHKERRQ(ierr);
85   for (i=0; i<aij->B->rmap->n/bs; i++) {
86     newRow = (*subMat)->rmap->range[subCommRank] + i*bs;
87     for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
88       newCol = garrayCMap[aijB->j[j]];
89       if (newCol) {
90         newCol--; /* remove the increment */
91         newCol *= bs;
92         for (k=0; k<bs; k++) {
93           newbRow[k] = newRow + k;
94           newbCol[k] = newCol + k;
95         }
96         /* copy column-oriented aijB->a into row-oriented vals */
97         aijBvals = aijB->a + j*bs*bs;
98         for (k1=0; k1<bs; k1++) {
99           for (k=0; k<bs; k++) {
100             vals[k1+k*bs] = *aijBvals++;
101           }
102         }
103         ierr = MatSetValues(*subMat,bs,newbRow,bs,newbCol,vals,INSERT_VALUES);CHKERRQ(ierr);
104       }
105     }
106   }
107   ierr = MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108   ierr = MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109 
110   /* deallocate temporary data */
111   ierr = PetscFree3(newbRow,newbCol,vals);CHKERRQ(ierr);
112   ierr = PetscFree(commRankMap);CHKERRQ(ierr);
113   ierr = PetscFree(garrayCMap);CHKERRQ(ierr);
114   if (scall == MAT_INITIAL_MATRIX) {
115     ierr = PetscFree(nnz);CHKERRQ(ierr);
116   }
117   PetscFunctionReturn(0);
118 }
119