Actual source code: mpb_aij.c
petsc-3.6.1 2015-08-06
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: MatSetBlockSizesFromMats(*subMat,mat,mat);
25: /* need to setup rmap and cmap before Preallocation */
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 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,&garrayCMap);
40: for (i=0; i<aij->B->cmap->n; i++) {
41: col = aij->garray[i];
42: for (subRank=0; subRank<subCommSize; subRank++) {
43: rank = commRankMap[subRank];
44: if ((col >= mat->cmap->range[rank]) && (col < mat->cmap->range[rank+1])) {
45: garrayCMap[i] = (*subMat)->cmap->range[subRank] + col - mat->cmap->range[rank]+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,&nnz);
54: for (i=0; i<aij->B->rmap->n; i++) {
55: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
56: if (garrayCMap[aijB->j[j]]) nnz[i]++;
57: }
58: }
59: MatMPIAIJSetPreallocation(*(subMat),0,NULL,0,nnz);
61: /* reuse diag block with the new submat */
62: MatDestroy(&((Mat_MPIAIJ*)((*subMat)->data))->A);
64: ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
66: PetscObjectReference((PetscObject)aij->A);
67: } else if (((Mat_MPIAIJ*)(*subMat)->data)->A != aij->A) {
68: PetscObject obj = (PetscObject)((Mat_MPIAIJ*)((*subMat)->data))->A;
70: PetscObjectReference((PetscObject)obj);
72: ((Mat_MPIAIJ*)((*subMat)->data))->A = aij->A;
74: PetscObjectReference((PetscObject)aij->A);
75: }
77: /* Now traverse aij->B and insert values into subMat */
78: for (i=0; i<aij->B->rmap->n; i++) {
79: newRow = (*subMat)->rmap->range[subCommRank] + i;
80: for (j=aijB->i[i]; j<aijB->i[i+1]; j++) {
81: newCol = garrayCMap[aijB->j[j]];
82: if (newCol) {
83: newCol--; /* remove the increment */
84: MatSetValues(*subMat,1,&newRow,1,&newCol,(aijB->a+j),INSERT_VALUES);
85: }
86: }
87: }
89: /* assemble the submat */
90: MatAssemblyBegin(*subMat,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(*subMat,MAT_FINAL_ASSEMBLY);
93: /* deallocate temporary data */
94: PetscFree(commRankMap);
95: PetscFree(garrayCMap);
96: if (scall == MAT_INITIAL_MATRIX) {
97: PetscFree(nnz);
98: }
99: return(0);
100: }