Actual source code: mpb_aij.c
petsc-3.4.2 2013-07-02
1: #include <../src/mat/impls/aij/mpi/mpiaij.h>
5: PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
6: {
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;
15: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);
16: MPI_Comm_size(subComm,&subCommSize);
18: /* create subMat object with the relavent layout */
19: if (scall == MAT_INITIAL_MATRIX) {
20: MatCreate(subComm,subMat);
21: MatSetType(*subMat,MATMPIAIJ);
22: MatSetSizes(*subMat,mat->rmap->n,mat->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
23: MatSetBlockSizes(*subMat,mat->rmap->bs,mat->cmap->bs);
25: /* need to setup rmap and cmap before Preallocation */
26: PetscLayoutSetBlockSize((*subMat)->rmap,mat->rmap->bs);
27: PetscLayoutSetBlockSize((*subMat)->cmap,mat->cmap->bs);
28: PetscLayoutSetUp((*subMat)->rmap);
29: PetscLayoutSetUp((*subMat)->cmap);
30: }
32: /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
33: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&commRank);
34: MPI_Comm_rank(subComm,&subCommRank);
35: PetscMalloc(subCommSize*sizeof(PetscMPIInt),&commRankMap);
36: MPI_Allgather(&commRank,1,MPI_INT,commRankMap,1,MPI_INT,subComm);
38: /* Traverse garray and identify column indices [of offdiag mat] that
39: should be discarded. For the ones not discarded, store the newCol+1
40: value in garrayCMap */
41: PetscMalloc(aij->B->cmap->n*sizeof(PetscInt),&garrayCMap);
42: PetscMemzero(garrayCMap,aij->B->cmap->n*sizeof(PetscInt));
43: for (i=0; i<aij->B->cmap->n; i++) {
44: col = aij->garray[i];
45: for (subRank=0; subRank<subCommSize; subRank++) {
46: rank = commRankMap[subRank];
47: if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) {
48: garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+1;
49: break;
50: }
51: }
52: }
54: if (scall == MAT_INITIAL_MATRIX) {
55: /* Now compute preallocation for the offdiag mat */
56: PetscMalloc(aij->B->rmap->n*sizeof(PetscInt),&nnz);
57: PetscMemzero(nnz,aij->B->rmap->n*sizeof(PetscInt));
58: for (i=0; i<aij->B->rmap->n; i++) {
59: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
60: if (garrayCMap[aijB->j[j]]) nnz[i]++;
61: }
62: }
63: MatMPIAIJSetPreallocation(*(subMat),0,NULL,0,nnz);
65: /* reuse diag block with the new submat */
66: MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);
68: ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
70: PetscObjectReference((PetscObject)aij->A);
71: } else if (((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A) {
72: PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A;
74: PetscObjectReference((PetscObject)obj);
76: ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
78: PetscObjectReference((PetscObject)aij->A);
79: }
81: /* Now traverse aij->B and insert values into subMat */
82: for (i=0; i<aij->B->rmap->n; i++) {
83: newRow = (*subMat)->rmap->range[subCommRank] + i;
84: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
85: newCol = garrayCMap[aijB->j[j]];
86: if (newCol) {
87: newCol--; /* remove the increment */
88: MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);
89: }
90: }
91: }
93: /* assemble the submat */
94: MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
95: MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);
97: /* deallocate temporary data */
98: PetscFree(commRankMap);
99: PetscFree(garrayCMap);
100: if (scall == MAT_INITIAL_MATRIX) {
101: PetscFree(nnz);
102: }
103: return(0);
104: }