Actual source code: mpb_baij.c
petsc-3.7.3 2016-08-01
1: #include <../src/mat/impls/baij/mpi/mpibaij.h>
5: PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
6: {
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;
17: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);
18: MPI_Comm_size(subComm,&subCommSize);
20: /* create subMat object with the relavent layout */
21: if (scall == MAT_INITIAL_MATRIX) {
22: MatCreate(subComm,subMat);
23: MatSetType(*subMat,MATMPIBAIJ);
24: MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
25: MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);
27: /* need to setup rmap and cmap before Preallocation */
28: PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);
29: PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);
30: PetscLayoutSetUp((*subMat)->rmap);
31: PetscLayoutSetUp((*subMat)->cmap);
32: }
34: /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
35: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);
36: MPI_Comm_rank(subComm,&subCommRank);
37: PetscMalloc1(subCommSize,&commRankMap);
38: MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);
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: PetscMalloc1(aij->B->cmap->n/bs,&garrayCMap);
44: PetscMemzero(garrayCMap,aij->B->cmap->n/bs*sizeof(PetscInt));
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: }
56: if (scall == MAT_INITIAL_MATRIX) {
57: /* Now compute preallocation for the offdiag mat */
58: PetscMalloc1(aij->B->rmap->n/bs,&nnz);
59: PetscMemzero(nnz,aij->B->rmap->n/bs*sizeof(PetscInt));
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: MatMPIBAIJSetPreallocation(*(subMat),bs,0,NULL,0,nnz);
67: /* reuse diag block with the new submat */
68: MatDestroy(&((Mat_MPIBAIJ*)((*subMat)->data))->A);
70: ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;
72: PetscObjectReference((PetscObject)aij->A);
73: } else if (((Mat_MPIBAIJ*)(*subMat)->data)->A != aij->A) {
74: PetscObject obj = (PetscObject)((Mat_MPIBAIJ*)((*subMat)->data))->A;
76: PetscObjectReference((PetscObject)obj);
78: ((Mat_MPIBAIJ*)((*subMat)->data))->A = aij->A;
80: PetscObjectReference((PetscObject)aij->A);
81: }
83: /* Now traverse aij->B and insert values into subMat */
84: PetscMalloc3(bs,&newbRow,bs,&newbCol,bs*bs,&vals);
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: MatSetValues(*subMat,bs,newbRow,bs,newbCol,vals,INSERT_VALUES);
104: }
105: }
106: }
107: MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);
110: /* deallocate temporary data */
111: PetscFree3(newbRow,newbCol,vals);
112: PetscFree(commRankMap);
113: PetscFree(garrayCMap);
114: if (scall == MAT_INITIAL_MATRIX) {
115: PetscFree(nnz);
116: }
117: return(0);
118: }