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