xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision e497d3c8e8f8bb645f0a9e4d8aea8ec4981f22b9)
1be1d678aSKris Buschelman 
22499ec78SHong Zhang /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
42499ec78SHong Zhang           C = A * B
52499ec78SHong Zhang */
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
11af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
121634b032SHong Zhang 
132499ec78SHong Zhang #undef __FUNCT__
142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
152499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
162499ec78SHong Zhang {
172499ec78SHong Zhang   PetscErrorCode ierr;
180fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
190d3441aeSHong Zhang   PetscInt       alg=1; /* set default algorithm */
202499ec78SHong Zhang 
212499ec78SHong Zhang   PetscFunctionBegin;
222499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
230fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
240fc8cf34SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
250fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
260fc8cf34SHong Zhang 
273ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
280fc8cf34SHong Zhang     switch (alg) {
290fc8cf34SHong Zhang     case 1:
300fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
310fc8cf34SHong Zhang       break;
320fc8cf34SHong Zhang     default:
33b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
340fc8cf34SHong Zhang       break;
350fc8cf34SHong Zhang     }
363ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
372499ec78SHong Zhang   }
383ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
39598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
403ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
412499ec78SHong Zhang   PetscFunctionReturn(0);
422499ec78SHong Zhang }
432499ec78SHong Zhang 
44f33d1a9aSHong Zhang #undef __FUNCT__
45a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
46a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
47598bc09dSHong Zhang {
48598bc09dSHong Zhang   PetscErrorCode ierr;
49598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
50598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
51598bc09dSHong Zhang 
52598bc09dSHong Zhang   PetscFunctionBegin;
53b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
54598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
55a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
56a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
57ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
58a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
60598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
61598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
62598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
63598bc09dSHong Zhang   PetscFunctionReturn(0);
64598bc09dSHong Zhang }
65598bc09dSHong Zhang 
662499ec78SHong Zhang #undef __FUNCT__
67a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
68a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
694ae0847bSHong Zhang {
704ae0847bSHong Zhang   PetscErrorCode ierr;
714ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
724ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
734ae0847bSHong Zhang 
744ae0847bSHong Zhang   PetscFunctionBegin;
754ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
7626fbe8dcSKarl Rupp 
774ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
784ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
794ae0847bSHong Zhang   PetscFunctionReturn(0);
804ae0847bSHong Zhang }
814ae0847bSHong Zhang 
82*e497d3c8SHong Zhang #if 0
83*e497d3c8SHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth */
84*e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
85c12557a1SHong Zhang {\
86c12557a1SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
87c12557a1SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
88c12557a1SHong Zhang   /* diagonal portion of A */\
89c12557a1SHong Zhang   _ai  = ad->i;\
90c12557a1SHong Zhang   _anz = _ai[i+1] - _ai[i];\
91c12557a1SHong Zhang   _aj  = ad->j + _ai[i];\
92c12557a1SHong Zhang   _aa  = ad->a + _ai[i];\
93c12557a1SHong Zhang   for (_j=0; _j<_anz; _j++) {\
94c12557a1SHong Zhang     _row = _aj[_j]; \
95c12557a1SHong Zhang     _pi  = p_loc->i;                                 \
96c12557a1SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
97c12557a1SHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
98c12557a1SHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
99c12557a1SHong Zhang     /* perform dense axpy */                    \
100c12557a1SHong Zhang     valtmp = _aa[_j];                           \
101c12557a1SHong Zhang     for (_k=0; _k<_pnz; _k++) {                    \
102c12557a1SHong Zhang       apa[_pj[_k]] += valtmp*_pa[_k];               \
103c12557a1SHong Zhang     }                                           \
104c12557a1SHong Zhang     PetscLogFlops(2.0*_pnz);                    \
105c12557a1SHong Zhang   }                                             \
106c12557a1SHong Zhang   /* off-diagonal portion of A */               \
107c12557a1SHong Zhang   _ai  = ao->i;\
108c12557a1SHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
109c12557a1SHong Zhang   _aj  = ao->j + _ai[i];                         \
110c12557a1SHong Zhang   _aa  = ao->a + _ai[i];                         \
111c12557a1SHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
112c12557a1SHong Zhang     _row = _aj[_j];    \
113c12557a1SHong Zhang     _pi  = p_oth->i;                         \
114c12557a1SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
115c12557a1SHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
116c12557a1SHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
117c12557a1SHong Zhang     /* perform dense axpy */                     \
118c12557a1SHong Zhang     _valtmp = _aa[_j];                             \
119c12557a1SHong Zhang     for (_k=0; _k<_pnz; _k++) {                     \
120c12557a1SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];                \
121c12557a1SHong Zhang     }                                            \
122c12557a1SHong Zhang     PetscLogFlops(2.0*_pnz);                     \
123c12557a1SHong Zhang   }                                              \
124c12557a1SHong Zhang }
125*e497d3c8SHong Zhang #endif
126c12557a1SHong Zhang 
1274ae0847bSHong Zhang #undef __FUNCT__
128f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
129f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
130598bc09dSHong Zhang {
131598bc09dSHong Zhang   PetscErrorCode ierr;
1324ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
133598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
134598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
135c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
136598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
137c12557a1SHong Zhang   PetscScalar    *apa,valtmp,*ca;
138c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
139598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
140c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
14129825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
142598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1439467f45fSHong Zhang   MPI_Comm       comm;
1449467f45fSHong Zhang   PetscMPIInt    size;
145598bc09dSHong Zhang 
146598bc09dSHong Zhang   PetscFunctionBegin;
1479467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1489467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1499467f45fSHong Zhang 
150a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
151598bc09dSHong Zhang   /*-----------------------------------------------------*/
152a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
153b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
154a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
155598bc09dSHong Zhang 
156598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
157598bc09dSHong Zhang   /*----------------------------------------------------------*/
158598bc09dSHong Zhang   /* get data from symbolic products */
159a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
160c12557a1SHong Zhang   p_oth = NULL;
1619467f45fSHong Zhang   if (size >1) {
1629467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1639467f45fSHong Zhang   }
164598bc09dSHong Zhang 
165598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
166598bc09dSHong Zhang   apa = ptap->apa;
167598bc09dSHong Zhang 
16829825010SHong Zhang   api = ptap->api;
16929825010SHong Zhang   apj = ptap->apj;
170598bc09dSHong Zhang   for (i=0; i<cm; i++) {
171c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
172*e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
173598bc09dSHong Zhang 
174598bc09dSHong Zhang     /* set values in C */
175598bc09dSHong Zhang     apJ  = apj + api[i];
176598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
177598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
178598bc09dSHong Zhang 
179598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
180598bc09dSHong Zhang     ca = coa + co->i[i];
181598bc09dSHong Zhang     k  = 0;
182598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
183598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
184598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
185598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
186598bc09dSHong Zhang       k++;
187598bc09dSHong Zhang     }
188598bc09dSHong Zhang 
189598bc09dSHong Zhang     /* diagonal part of C */
190598bc09dSHong Zhang     ca = cda + cd->i[i];
191598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
192598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
193598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
194598bc09dSHong Zhang       k++;
195598bc09dSHong Zhang     }
196598bc09dSHong Zhang 
197598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
198598bc09dSHong Zhang     ca = coa + co->i[i];
199598bc09dSHong Zhang     for (; k0<conz; k0++) {
200598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
201598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
202598bc09dSHong Zhang       k++;
203598bc09dSHong Zhang     }
204598bc09dSHong Zhang   }
205598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
206598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207598bc09dSHong Zhang   PetscFunctionReturn(0);
208598bc09dSHong Zhang }
209598bc09dSHong Zhang 
210598bc09dSHong Zhang #undef __FUNCT__
2110fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
2120fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
213598bc09dSHong Zhang {
214598bc09dSHong Zhang   PetscErrorCode     ierr;
215ce94432eSBarry Smith   MPI_Comm           comm;
2169467f45fSHong Zhang   PetscMPIInt        size;
217bfcd1a73SHong Zhang   Mat                Cmpi;
218598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2190298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2204ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
221bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2224ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2234ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
224d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
225bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
226598bc09dSHong Zhang   PetscBT            lnkbt;
227598bc09dSHong Zhang   PetscScalar        *apa;
228bfcd1a73SHong Zhang   PetscReal          afill;
229f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
230598bc09dSHong Zhang 
231598bc09dSHong Zhang   PetscFunctionBegin;
232ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2339467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2349467f45fSHong Zhang 
235a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
236cf1a0672SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
237a1a4e84aSHong Zhang   }
238152983bfSHong Zhang 
239598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
240b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
241598bc09dSHong Zhang 
242598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
243b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
24419f0e57aSHong Zhang 
245598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
246a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
247598bc09dSHong Zhang 
248a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
249598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2509467f45fSHong Zhang   if (size > 1) {
2519467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
252598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
25320e3dc75SHong Zhang   } else {
25420e3dc75SHong Zhang     p_oth = NULL;
25520e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2569467f45fSHong Zhang   }
257598bc09dSHong Zhang 
258598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
259598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
260854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
261a1a4e84aSHong Zhang   ptap->api = api;
262598bc09dSHong Zhang   api[0]    = 0;
263598bc09dSHong Zhang 
264598bc09dSHong Zhang   /* create and initialize a linked list */
265f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
2669467f45fSHong Zhang   if (size >1) {
267f84c536bSHong Zhang     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
2689467f45fSHong Zhang   } else {
2699467f45fSHong Zhang     prmax = p_loc->rmax;
2709467f45fSHong Zhang   }
271f84c536bSHong Zhang   nlnk_max = armax*prmax;
272f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2730d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
274598bc09dSHong Zhang 
275bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
276bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2772205254eSKarl Rupp 
278598bc09dSHong Zhang   current_space = free_space;
279598bc09dSHong Zhang 
280598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
281598bc09dSHong Zhang   for (i=0; i<am; i++) {
282598bc09dSHong Zhang     /* diagonal portion of A */
283598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
284598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
285598bc09dSHong Zhang       row  = *adj++;
286598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
287598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
288598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2891fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
290598bc09dSHong Zhang     }
291598bc09dSHong Zhang     /* off-diagonal portion of A */
292598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
293598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
294598bc09dSHong Zhang       row  = *aoj++;
295598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
296598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2971fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
298598bc09dSHong Zhang     }
299598bc09dSHong Zhang 
300d14fa243SHong Zhang     apnz     = lnk[0];
301598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
302598bc09dSHong Zhang 
303598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
304598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
305598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
306598bc09dSHong Zhang       nspacedouble++;
307598bc09dSHong Zhang     }
308598bc09dSHong Zhang 
309598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
3101fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
311598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
3122205254eSKarl Rupp 
313598bc09dSHong Zhang     current_space->array           += apnz;
314598bc09dSHong Zhang     current_space->local_used      += apnz;
315598bc09dSHong Zhang     current_space->local_remaining -= apnz;
316598bc09dSHong Zhang   }
317598bc09dSHong Zhang 
318598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
319598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
320854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
321a1a4e84aSHong Zhang   apj  = ptap->apj;
322a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
323598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
324598bc09dSHong Zhang 
325f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
3261795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
3272205254eSKarl Rupp 
328f84c536bSHong Zhang   ptap->apa = apa;
329f84c536bSHong Zhang 
330bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
331598bc09dSHong Zhang   /*----------------------------------------------------*/
332bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
333bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
33433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
335a2f3521dSMark F. Adams 
336bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
337bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
338598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
339598bc09dSHong Zhang   for (i=0; i<am; i++) {
340598bc09dSHong Zhang     row  = i + rstart;
341598bc09dSHong Zhang     apnz = api[i+1] - api[i];
342bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
343598bc09dSHong Zhang     apj += apnz;
344598bc09dSHong Zhang   }
345bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
346bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
347598bc09dSHong Zhang 
348bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
349bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
350f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
351bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
352bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
353598bc09dSHong Zhang 
354bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
355bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
356598bc09dSHong Zhang   c->ptap = ptap;
357598bc09dSHong Zhang 
358bfcd1a73SHong Zhang   *C = Cmpi;
359bfcd1a73SHong Zhang 
360bfcd1a73SHong Zhang   /* set MatInfo */
361118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
362bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
363bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
364bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
365bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
366bfcd1a73SHong Zhang 
367bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
368bfcd1a73SHong Zhang   if (api[am]) {
36957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
37057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
371bfcd1a73SHong Zhang   } else {
372bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
373bfcd1a73SHong Zhang   }
374bfcd1a73SHong Zhang #endif
375598bc09dSHong Zhang   PetscFunctionReturn(0);
376598bc09dSHong Zhang }
377598bc09dSHong Zhang 
3789123193aSHong Zhang #undef __FUNCT__
3799123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3809123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3819123193aSHong Zhang {
3829123193aSHong Zhang   PetscErrorCode ierr;
3839123193aSHong Zhang 
3849123193aSHong Zhang   PetscFunctionBegin;
3859123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3863ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3879123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3883ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3899123193aSHong Zhang   }
3903ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3919123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3923ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3939123193aSHong Zhang   PetscFunctionReturn(0);
3949123193aSHong Zhang }
3959123193aSHong Zhang 
3964324174eSBarry Smith typedef struct {
3974324174eSBarry Smith   Mat         workB;
3984324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3994324174eSBarry Smith   MPI_Request *rwaits,*swaits;
4004324174eSBarry Smith } MPIAIJ_MPIDense;
4014324174eSBarry Smith 
4027af9e4fdSHong Zhang #undef __FUNCT__
40396b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
40496b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
4054324174eSBarry Smith {
4064324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
4074324174eSBarry Smith   PetscErrorCode  ierr;
4084324174eSBarry Smith 
4094324174eSBarry Smith   PetscFunctionBegin;
4106bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
4114324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
4124324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
4134324174eSBarry Smith   PetscFunctionReturn(0);
4144324174eSBarry Smith }
4154324174eSBarry Smith 
4169123193aSHong Zhang #undef __FUNCT__
417fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense"
418fd4e9aacSBarry Smith /*
419fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
420fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
421fd4e9aacSBarry Smith 
422fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
423fd4e9aacSBarry Smith */
424fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
425fd4e9aacSBarry Smith {
426fd4e9aacSBarry Smith   PetscErrorCode         ierr;
427fd4e9aacSBarry Smith   PetscBool              flg;
428fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
429fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
430fd4e9aacSBarry Smith   PetscContainer         container;
431fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
432fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
433fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
434fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
435fd4e9aacSBarry Smith 
436fd4e9aacSBarry Smith   PetscFunctionBegin;
437fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
438fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
439fd4e9aacSBarry Smith 
440fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
441fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
442fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
443fd4e9aacSBarry Smith 
444fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
445fd4e9aacSBarry Smith 
446fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
447fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
448fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
449fd4e9aacSBarry Smith   /* Create work arrays needed */
450fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
451fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
452fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
453fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
454fd4e9aacSBarry Smith 
455fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
456fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
457fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
458fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
459fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
460fd4e9aacSBarry Smith 
461fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
462fd4e9aacSBarry Smith   PetscFunctionReturn(0);
463fd4e9aacSBarry Smith }
464fd4e9aacSBarry Smith 
465fd4e9aacSBarry Smith #undef __FUNCT__
4669123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4679123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4689123193aSHong Zhang {
4699123193aSHong Zhang   PetscErrorCode         ierr;
470f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
471d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
472bf0cc555SLisandro Dalcin   PetscContainer         container;
4734324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4744324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4754324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4764324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
477d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4789123193aSHong Zhang 
4799123193aSHong Zhang   PetscFunctionBegin;
480ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
481d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
48233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
483cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4840298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
485cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
486cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4872205254eSKarl Rupp 
488d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
489f32d5d43SBarry Smith 
490b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
491f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4920298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4934324174eSBarry Smith   /* Create work arrays needed */
494dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
495dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
496dcca6d9dSJed Brown                       from->n,&contents->rwaits,
497dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4984324174eSBarry Smith 
499ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
500bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
50196b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
502bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
503bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5049123193aSHong Zhang   PetscFunctionReturn(0);
5059123193aSHong Zhang }
5069123193aSHong Zhang 
5077af9e4fdSHong Zhang #undef __FUNCT__
5087af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
509f32d5d43SBarry Smith /*
5102636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
5112636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
512f32d5d43SBarry Smith */
5134324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
514f32d5d43SBarry Smith {
515f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
516f32d5d43SBarry Smith   PetscErrorCode         ierr;
517f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
518f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
519f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
520f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
521f32d5d43SBarry Smith   PetscInt               i,j,k;
522aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
523aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
524f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
525ce94432eSBarry Smith   MPI_Comm               comm;
526d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
527f32d5d43SBarry Smith   MPI_Status             status;
5284324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
529bf0cc555SLisandro Dalcin   PetscContainer         container;
5304324174eSBarry Smith   Mat                    workB;
531f32d5d43SBarry Smith 
532f32d5d43SBarry Smith   PetscFunctionBegin;
533ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
534bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
53529825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
536bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
5374324174eSBarry Smith 
5384324174eSBarry Smith   workB = *outworkB = contents->workB;
539cf1a0672SHong Zhang   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
540f32d5d43SBarry Smith   sindices = to->indices;
541f32d5d43SBarry Smith   sstarts  = to->starts;
542f32d5d43SBarry Smith   sprocs   = to->procs;
5434324174eSBarry Smith   swaits   = contents->swaits;
5444324174eSBarry Smith   svalues  = contents->svalues;
545f32d5d43SBarry Smith 
546f32d5d43SBarry Smith   rindices = from->indices;
547f32d5d43SBarry Smith   rstarts  = from->starts;
548f32d5d43SBarry Smith   rprocs   = from->procs;
5494324174eSBarry Smith   rwaits   = contents->rwaits;
5504324174eSBarry Smith   rvalues  = contents->rvalues;
551f32d5d43SBarry Smith 
5528c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5538c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
554f32d5d43SBarry Smith 
555f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
556f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
557f32d5d43SBarry Smith   }
558f32d5d43SBarry Smith 
559f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
560f32d5d43SBarry Smith     /* pack a message at a time */
561f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
562f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
563f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5642636ff26SBarry Smith       }
5652636ff26SBarry Smith     }
566f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
567f32d5d43SBarry Smith   }
5682636ff26SBarry Smith 
569f32d5d43SBarry Smith   nrecvs = from->n;
570f32d5d43SBarry Smith   while (nrecvs) {
571f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
572f32d5d43SBarry Smith     nrecvs--;
573f32d5d43SBarry Smith     /* unpack a message at a time */
574f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
575f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
576f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5772636ff26SBarry Smith       }
5782636ff26SBarry Smith     }
579f32d5d43SBarry Smith   }
580cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
581f32d5d43SBarry Smith 
5828c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5838c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
584f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
585f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
586f32d5d43SBarry Smith   PetscFunctionReturn(0);
587f32d5d43SBarry Smith }
588f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
589f32d5d43SBarry Smith 
5909123193aSHong Zhang #undef __FUNCT__
5919123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5929123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5939123193aSHong Zhang {
5949123193aSHong Zhang   PetscErrorCode ierr;
595f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
596f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
597f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
598f32d5d43SBarry Smith   Mat            workB;
5999123193aSHong Zhang 
6009123193aSHong Zhang   PetscFunctionBegin;
601f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
602f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
603f32d5d43SBarry Smith 
604f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
6054324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
606f32d5d43SBarry Smith 
607f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
608f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
6099123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6109123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6119123193aSHong Zhang   PetscFunctionReturn(0);
6129123193aSHong Zhang }
613cf1a0672SHong Zhang 
6141634b032SHong Zhang #undef __FUNCT__
615f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
616f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
6171634b032SHong Zhang {
6181634b032SHong Zhang   PetscErrorCode ierr;
619cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
620cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
621cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
622cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
623cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
624cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
625cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
62629825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
62729825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
628cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
62929825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
630cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
63129825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
63229825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
633a7c7454dSHong Zhang   MPI_Comm       comm;
634a7c7454dSHong Zhang   PetscMPIInt    size;
6351634b032SHong Zhang 
6361634b032SHong Zhang   PetscFunctionBegin;
637a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
638a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
639a7c7454dSHong Zhang 
640cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
641cf1a0672SHong Zhang   /*-----------------------------------------------------*/
642cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
643b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
644cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
645cf1a0672SHong Zhang 
646cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
647cf1a0672SHong Zhang   /*----------------------------------------------------------*/
648cf1a0672SHong Zhang   /* get data from symbolic products */
649cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
650cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
651a7c7454dSHong Zhang   if (size >1) {
652a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
653cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
65420e3dc75SHong Zhang   } else {
65520e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
656a7c7454dSHong Zhang   }
657cf1a0672SHong Zhang 
65829825010SHong Zhang   api = ptap->api;
65929825010SHong Zhang   apj = ptap->apj;
660cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
66129825010SHong Zhang     apJ = apj + api[i];
66229825010SHong Zhang 
663cf1a0672SHong Zhang     /* diagonal portion of A */
664cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
665cf1a0672SHong Zhang     adj = ad->j + adi[i];
666cf1a0672SHong Zhang     ada = ad->a + adi[i];
667cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
668cf1a0672SHong Zhang       row = adj[j];
669cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
670cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
671cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
67229825010SHong Zhang       /* perform sparse axpy */
673cf1a0672SHong Zhang       valtmp = ada[j];
67429825010SHong Zhang       nextp  = 0;
67529825010SHong Zhang       for (k=0; nextp<pnz; k++) {
67629825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
67729825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
67829825010SHong Zhang         }
6791634b032SHong Zhang       }
680cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
681cf1a0672SHong Zhang     }
6821634b032SHong Zhang 
683cf1a0672SHong Zhang     /* off-diagonal portion of A */
684cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
685cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
686cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
687cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
688cf1a0672SHong Zhang       row = aoj[j];
689cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
690cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
691cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
69229825010SHong Zhang       /* perform sparse axpy */
693cf1a0672SHong Zhang       valtmp = aoa[j];
69429825010SHong Zhang       nextp  = 0;
69529825010SHong Zhang       for (k=0; nextp<pnz; k++) {
69629825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
69729825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
69829825010SHong Zhang         }
699cf1a0672SHong Zhang       }
700cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
701cf1a0672SHong Zhang     }
7021634b032SHong Zhang 
703cf1a0672SHong Zhang     /* set values in C */
704cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
705cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
7061634b032SHong Zhang 
707cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
708cf1a0672SHong Zhang     ca = coa + co->i[i];
709cf1a0672SHong Zhang     k  = 0;
710cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
711cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
71229825010SHong Zhang       ca[k0]        = apa_sparse[k];
71329825010SHong Zhang       apa_sparse[k] = 0.0;
714cf1a0672SHong Zhang       k++;
715cf1a0672SHong Zhang     }
7161634b032SHong Zhang 
717cf1a0672SHong Zhang     /* diagonal part of C */
718cf1a0672SHong Zhang     ca = cda + cd->i[i];
719cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
72029825010SHong Zhang       ca[k1]        = apa_sparse[k];
72129825010SHong Zhang       apa_sparse[k] = 0.0;
722cf1a0672SHong Zhang       k++;
723cf1a0672SHong Zhang     }
724cf1a0672SHong Zhang 
725cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
726cf1a0672SHong Zhang     ca = coa + co->i[i];
727cf1a0672SHong Zhang     for (; k0<conz; k0++) {
72829825010SHong Zhang       ca[k0]        = apa_sparse[k];
72929825010SHong Zhang       apa_sparse[k] = 0.0;
730cf1a0672SHong Zhang       k++;
731cf1a0672SHong Zhang     }
732cf1a0672SHong Zhang   }
733cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
734cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
735cf1a0672SHong Zhang   PetscFunctionReturn(0);
736cf1a0672SHong Zhang }
737cf1a0672SHong Zhang 
7380fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
739cf1a0672SHong Zhang #undef __FUNCT__
740b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
741b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
742cf1a0672SHong Zhang {
743cf1a0672SHong Zhang   PetscErrorCode     ierr;
744ce94432eSBarry Smith   MPI_Comm           comm;
745a7c7454dSHong Zhang   PetscMPIInt        size;
746cf1a0672SHong Zhang   Mat                Cmpi;
747cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
7480298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
749cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
750cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
751cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
752cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
753f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
754cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
755cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
756cf1a0672SHong Zhang   PetscReal          afill;
757cf1a0672SHong Zhang   PetscScalar        *apa;
758cf1a0672SHong Zhang 
759cf1a0672SHong Zhang   PetscFunctionBegin;
760ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
761a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
762a7c7454dSHong Zhang 
763cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
764b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
765cf1a0672SHong Zhang 
766cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
767b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
76819f0e57aSHong Zhang 
769cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
770cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
771cf1a0672SHong Zhang 
772cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
773cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
774a7c7454dSHong Zhang   if (size > 1) {
775a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
77620e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
777968382fdSHong Zhang   } else {
77820e3dc75SHong Zhang     p_oth  = NULL;
77920e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
780a7c7454dSHong Zhang   }
781cf1a0672SHong Zhang 
782cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
783cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
784854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
785cf1a0672SHong Zhang   ptap->api = api;
786cf1a0672SHong Zhang   api[0]    = 0;
787cf1a0672SHong Zhang 
788cf1a0672SHong Zhang   /* create and initialize a linked list */
789cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
790a7c7454dSHong Zhang   if (size >1) {
791cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
792a7c7454dSHong Zhang   } else {
793a7c7454dSHong Zhang     prmax = p_loc->rmax;
794a7c7454dSHong Zhang   }
795cf1a0672SHong Zhang   nlnk_max = armax*prmax;
796cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
797f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
798cf1a0672SHong Zhang 
799cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
800cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
8012205254eSKarl Rupp 
802cf1a0672SHong Zhang   current_space = free_space;
803cf1a0672SHong Zhang 
804cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
805cf1a0672SHong Zhang   for (i=0; i<am; i++) {
806cf1a0672SHong Zhang     /* diagonal portion of A */
807cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
808cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
809cf1a0672SHong Zhang       row  = *adj++;
810cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
811cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
812cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
813f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
814cf1a0672SHong Zhang     }
815cf1a0672SHong Zhang     /* off-diagonal portion of A */
816cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
817cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
818cf1a0672SHong Zhang       row  = *aoj++;
819cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
820cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
821f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
822cf1a0672SHong Zhang     }
823cf1a0672SHong Zhang 
824f84c536bSHong Zhang     apnz     = *lnk;
825cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
826cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
827cf1a0672SHong Zhang 
828cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
829cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
830cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
831cf1a0672SHong Zhang       nspacedouble++;
832cf1a0672SHong Zhang     }
833cf1a0672SHong Zhang 
834cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
835f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
836cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8372205254eSKarl Rupp 
838cf1a0672SHong Zhang     current_space->array           += apnz;
839cf1a0672SHong Zhang     current_space->local_used      += apnz;
840cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
841cf1a0672SHong Zhang   }
842cf1a0672SHong Zhang 
843cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
844cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
845854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
846cf1a0672SHong Zhang   apj  = ptap->apj;
847cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
848f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
849cf1a0672SHong Zhang 
850cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
851cf1a0672SHong Zhang   /*----------------------------------------------------*/
852cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
853cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
85433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
855cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
856cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
857cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
858cf1a0672SHong Zhang 
85929825010SHong Zhang   /* malloc apa for assembly Cmpi */
8601795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8612205254eSKarl Rupp 
86229825010SHong Zhang   ptap->apa = apa;
863cf1a0672SHong Zhang   for (i=0; i<am; i++) {
864cf1a0672SHong Zhang     row  = i + rstart;
865cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
866cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
867cf1a0672SHong Zhang     apj += apnz;
868cf1a0672SHong Zhang   }
869cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
870cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
871cf1a0672SHong Zhang 
872cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
873cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
874f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
875cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
876cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
877cf1a0672SHong Zhang 
878cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
879cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
880cf1a0672SHong Zhang   c->ptap = ptap;
881cf1a0672SHong Zhang 
882cf1a0672SHong Zhang   *C = Cmpi;
883cf1a0672SHong Zhang 
884cf1a0672SHong Zhang   /* set MatInfo */
885118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
886cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
887cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
888cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
889cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
890cf1a0672SHong Zhang 
891cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
892cf1a0672SHong Zhang   if (api[am]) {
89357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
89457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
895cf1a0672SHong Zhang   } else {
896cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
897cf1a0672SHong Zhang   }
898cf1a0672SHong Zhang #endif
8991634b032SHong Zhang   PetscFunctionReturn(0);
9001634b032SHong Zhang }
901f96d37fcSHong Zhang 
902f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
903f96d37fcSHong Zhang #undef __FUNCT__
904f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
905c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
906f96d37fcSHong Zhang {
907f96d37fcSHong Zhang   PetscErrorCode ierr;
908c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
909c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
910f96d37fcSHong Zhang 
911f96d37fcSHong Zhang   PetscFunctionBegin;
912c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
913c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
914c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
915c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
916c216dbf3SHong Zhang 
917c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
918c216dbf3SHong Zhang     switch (alg) {
919c216dbf3SHong Zhang     case 1:
9202bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
921c216dbf3SHong Zhang       break;
922c216dbf3SHong Zhang     case 2:
923275476c6SMatthew G. Knepley     {
924ab1b013aSHong Zhang       Mat         Pt;
925ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
926ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
927ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
928ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
929ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
930ab1b013aSHong Zhang       ptap     = c->ptap;
931ab1b013aSHong Zhang       ptap->Pt = Pt;
932c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
933c216dbf3SHong Zhang       PetscFunctionReturn(0);
934275476c6SMatthew G. Knepley     }
935c216dbf3SHong Zhang       break;
936c216dbf3SHong Zhang     default:
9376da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
938c216dbf3SHong Zhang       break;
939c216dbf3SHong Zhang     }
940c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
941c216dbf3SHong Zhang   }
942c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
943c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
944c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
945ab1b013aSHong Zhang   PetscFunctionReturn(0);
946ab1b013aSHong Zhang }
947ab1b013aSHong Zhang 
948c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
949c216dbf3SHong Zhang #undef __FUNCT__
950c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
951c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
952c216dbf3SHong Zhang {
953c216dbf3SHong Zhang   PetscErrorCode ierr;
9542bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9552bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9562bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
957c216dbf3SHong Zhang 
958c216dbf3SHong Zhang   PetscFunctionBegin;
959c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
960c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
961f96d37fcSHong Zhang   PetscFunctionReturn(0);
962f96d37fcSHong Zhang }
963f96d37fcSHong Zhang 
9646da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
965f96d37fcSHong Zhang #undef __FUNCT__
9666da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9676da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
968f96d37fcSHong Zhang {
969c5af039cSHong Zhang   PetscErrorCode      ierr;
970c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
971497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
972c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
973c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
974d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
975497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
976e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
977c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
978ce94432eSBarry Smith   MPI_Comm            comm;
979c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
980c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
981c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
982c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
983c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
984c5af039cSHong Zhang   MPI_Status          *status;
985c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
986d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
987a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
988c5af039cSHong Zhang   Mat                 A_loc;
989c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
990f96d37fcSHong Zhang 
991f96d37fcSHong Zhang   PetscFunctionBegin;
992ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
993c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
994c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
995c5af039cSHong Zhang 
996c5af039cSHong Zhang   ptap  = c->ptap;
997c5af039cSHong Zhang   merge = ptap->merge;
998c5af039cSHong Zhang 
999c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1000c5af039cSHong Zhang   /*--------------------------------------------------------------*/
1001c5af039cSHong Zhang   /* get data from symbolic products */
1002c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
1003854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1004c5af039cSHong Zhang 
1005c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
1006c5af039cSHong Zhang   owners = merge->rowmap->range;
1007854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1008c5af039cSHong Zhang 
1009c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
1010c5af039cSHong Zhang   A_loc = ptap->A_loc;
1011c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1012c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1013d6ab1dc0SHong Zhang   ai    = a_loc->i;
1014d6ab1dc0SHong Zhang   aj    = a_loc->j;
1015c5af039cSHong Zhang 
1016854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
1017c5af039cSHong Zhang 
1018c5af039cSHong Zhang   for (i=0; i<am; i++) {
1019e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
1020d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1021d6ab1dc0SHong Zhang     adj = aj + ai[i];
1022d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1023c5af039cSHong Zhang     for (j=0; j<anz; j++) {
1024e745005bSHong Zhang       aval[adj[j]] = ada[j];
1025c5af039cSHong Zhang     }
1026c5af039cSHong Zhang 
1027c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1028c5af039cSHong Zhang     /*--------------------------------------------------------------*/
1029c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1030c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
1031c5af039cSHong Zhang     poJ = po->j + po->i[i];
1032c5af039cSHong Zhang     pA  = po->a + po->i[i];
1033c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1034c5af039cSHong Zhang       row = poJ[j];
1035c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
1036c5af039cSHong Zhang       cj  = coj + coi[row];
1037c5af039cSHong Zhang       ca  = coa + coi[row];
1038c5af039cSHong Zhang       /* perform dense axpy */
1039c5af039cSHong Zhang       valtmp = pA[j];
1040c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1041e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1042c5af039cSHong Zhang       }
1043c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1044c5af039cSHong Zhang     }
1045c5af039cSHong Zhang 
1046c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1047c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1048c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1049c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1050c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1051c5af039cSHong Zhang       row = pdJ[j];
1052c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1053c5af039cSHong Zhang       cj  = bj + bi[row];
1054c5af039cSHong Zhang       ca  = ba + bi[row];
1055c5af039cSHong Zhang       /* perform dense axpy */
1056c5af039cSHong Zhang       valtmp = pA[j];
1057c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1058e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1059c5af039cSHong Zhang       }
1060c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1061c5af039cSHong Zhang     }
1062c5af039cSHong Zhang 
1063d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1064d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1065e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1066c5af039cSHong Zhang   }
1067c5af039cSHong Zhang 
1068c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1069c5af039cSHong Zhang   /*------------------------------------*/
1070c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1071c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1072c5af039cSHong Zhang   len_s  = merge->len_s;
1073c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1074c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1075c5af039cSHong Zhang 
1076dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1077c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1078c5af039cSHong Zhang     if (!len_s[proc]) continue;
1079c5af039cSHong Zhang     i    = merge->owners_co[proc];
1080c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1081c5af039cSHong Zhang     k++;
1082c5af039cSHong Zhang   }
1083c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1084c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1085c5af039cSHong Zhang 
1086c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1087c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1088c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1089c5af039cSHong Zhang 
1090c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1091c5af039cSHong Zhang   /*----------------------------------------------------*/
1092dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1093c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1094c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1095c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1096c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1097c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1098c5af039cSHong Zhang   }
1099c5af039cSHong Zhang 
1100c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1101c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1102c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1103c5af039cSHong Zhang     ba_i = ba + bi[i];
1104c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1105c5af039cSHong Zhang     /* add received vals into ba */
1106c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1107c5af039cSHong Zhang       /* i-th row */
1108c5af039cSHong Zhang       if (i == *nextrow[k]) {
1109c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1110c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1111c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1112c5af039cSHong Zhang         nextcj = 0;
1113c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1114c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1115c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1116c5af039cSHong Zhang           }
1117c5af039cSHong Zhang         }
1118c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1119c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1120c5af039cSHong Zhang       }
1121c5af039cSHong Zhang     }
1122c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1123c5af039cSHong Zhang   }
1124c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1125c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1126c5af039cSHong Zhang 
1127c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1128c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1129c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1130c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1131e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1132f96d37fcSHong Zhang   PetscFunctionReturn(0);
1133f96d37fcSHong Zhang }
1134f96d37fcSHong Zhang 
1135c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1136c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1137f96d37fcSHong Zhang #undef __FUNCT__
11382bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
11392bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1140f96d37fcSHong Zhang {
1141f96d37fcSHong Zhang   PetscErrorCode      ierr;
11424e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1143f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
11440298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1145f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1146f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1147f96d37fcSHong Zhang   PetscInt            nnz;
11484e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1149497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1150f96d37fcSHong Zhang   PetscBT             lnkbt;
1151ce94432eSBarry Smith   MPI_Comm            comm;
1152f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1153f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1154f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1155f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1156f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1157f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1158f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1159f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1160d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1161f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1162c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1163f96d37fcSHong Zhang   PetscScalar         *vals;
11644e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1165f96d37fcSHong Zhang 
1166f96d37fcSHong Zhang   PetscFunctionBegin;
1167ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1168c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1169c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1170c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1171c5af039cSHong Zhang   }
1172c5af039cSHong Zhang 
1173f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1174f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1175f96d37fcSHong Zhang 
1176f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1177b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1178f96d37fcSHong Zhang 
1179f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1180f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11812205254eSKarl Rupp 
1182c5af039cSHong Zhang   ptap->A_loc = A_loc;
11832205254eSKarl Rupp 
11841c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1185d6ab1dc0SHong Zhang   ai    = a_loc->i;
1186d6ab1dc0SHong Zhang   aj    = a_loc->j;
1187f96d37fcSHong Zhang 
1188f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1189f96d37fcSHong Zhang   /*----------------------------------------------------*/
11904e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11914e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11924e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11934e938277SHong Zhang 
11944e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11954e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11964e938277SHong Zhang   poti = pot->i; potj = pot->j;
1197f96d37fcSHong Zhang 
1198f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11992205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1200854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1201f96d37fcSHong Zhang   coi[0] = 0;
1202f96d37fcSHong Zhang 
1203f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1204d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1205a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1206f96d37fcSHong Zhang   current_space = free_space;
120719f0e57aSHong Zhang 
1208f96d37fcSHong Zhang   /* create and initialize a linked list */
12094e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1210c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
12114e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
12124e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1213f96d37fcSHong Zhang 
1214f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1215f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1216f96d37fcSHong Zhang     ptJ = potj + poti[i];
1217f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1218f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1219d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1220d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1221f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1222d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1223f96d37fcSHong Zhang     }
12244e938277SHong Zhang     nnz = lnk[0];
1225f96d37fcSHong Zhang 
1226f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1227f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1228f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1229f96d37fcSHong Zhang       nspacedouble++;
1230f96d37fcSHong Zhang     }
1231f96d37fcSHong Zhang 
1232f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
12334e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
12342205254eSKarl Rupp 
1235f96d37fcSHong Zhang     current_space->array           += nnz;
1236f96d37fcSHong Zhang     current_space->local_used      += nnz;
1237f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
12382205254eSKarl Rupp 
1239f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1240f96d37fcSHong Zhang   }
1241f96d37fcSHong Zhang 
1242854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1243f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
12442205254eSKarl Rupp 
1245118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1246f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1247f96d37fcSHong Zhang 
1248f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1249f96d37fcSHong Zhang   /*----------------------------------------------*/
1250f96d37fcSHong Zhang   /* determine row ownership */
1251b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1252f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12532205254eSKarl Rupp 
1254f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1255f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12562205254eSKarl Rupp 
1257f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1258f96d37fcSHong Zhang   owners = merge->rowmap->range;
1259f96d37fcSHong Zhang 
1260f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12611795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1262785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12632205254eSKarl Rupp 
1264f96d37fcSHong Zhang   len_s        = merge->len_s;
1265f96d37fcSHong Zhang   merge->nsend = 0;
1266f96d37fcSHong Zhang 
1267854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1268f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1269f96d37fcSHong Zhang 
1270f96d37fcSHong Zhang   proc = 0;
1271f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1272f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1273f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1274f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1275f96d37fcSHong Zhang   }
1276f96d37fcSHong Zhang 
1277f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1278f96d37fcSHong Zhang   owners_co[0] = 0;
1279f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1280f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1281f96d37fcSHong Zhang     if (len_si[proc]) {
1282f96d37fcSHong Zhang       merge->nsend++;
1283f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1284f96d37fcSHong Zhang       len         += len_si[proc];
1285f96d37fcSHong Zhang     }
1286f96d37fcSHong Zhang   }
1287f96d37fcSHong Zhang 
1288f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12890298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1290f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1291f96d37fcSHong Zhang 
1292f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1293f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1294f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1295854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1297f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1298f96d37fcSHong Zhang     i    = owners_co[proc];
1299f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1300f96d37fcSHong Zhang     k++;
1301f96d37fcSHong Zhang   }
1302f96d37fcSHong Zhang 
1303f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1304785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1305f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1306c280ed6aSJed Brown     PetscMPIInt icompleted;
1307c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1308f96d37fcSHong Zhang   }
1309f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1310f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1311f96d37fcSHong Zhang 
1312f96d37fcSHong Zhang   /* send and recv coi */
1313f96d37fcSHong Zhang   /*-------------------*/
1314f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1315f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1316854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1317f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1318f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1319f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1320f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1321f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1322f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1323f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1324f96d37fcSHong Zhang     */
1325f96d37fcSHong Zhang     /*-------------------------------------------*/
1326f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1327f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1328f96d37fcSHong Zhang     buf_si[0]   = nrows;
1329f96d37fcSHong Zhang     buf_si_i[0] = 0;
1330f96d37fcSHong Zhang     nrows       = 0;
1331f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1332f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1333f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1334f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1335f96d37fcSHong Zhang       nrows++;
1336f96d37fcSHong Zhang     }
1337f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1338f96d37fcSHong Zhang     k++;
1339f96d37fcSHong Zhang     buf_si += len_si[proc];
1340f96d37fcSHong Zhang   }
1341f96d37fcSHong Zhang   i = merge->nrecv;
1342f96d37fcSHong Zhang   while (i--) {
1343c280ed6aSJed Brown     PetscMPIInt icompleted;
1344c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1345f96d37fcSHong Zhang   }
1346f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1347f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1348f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1349f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1350f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1351f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1352f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1353f96d37fcSHong Zhang 
1354f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1355f96d37fcSHong Zhang   /*------------------------------------------*/
1356f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1357854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1358f96d37fcSHong Zhang   bi[0] = 0;
1359f96d37fcSHong Zhang 
1360c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1361d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1362a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1363f96d37fcSHong Zhang   current_space = free_space;
1364f96d37fcSHong Zhang 
1365dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1366f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1367f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1368f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1369f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1370f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1371f96d37fcSHong Zhang   }
1372f96d37fcSHong Zhang 
13731c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1374f96d37fcSHong Zhang   rmax = 0;
1375f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1376f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1377f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1378f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1379f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1380f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1381d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1382d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1383f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1384d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1385f96d37fcSHong Zhang     }
13864e938277SHong Zhang 
1387f96d37fcSHong Zhang     /* add received col data into lnk */
1388f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1389f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1390f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1391f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13924e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1393f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1394f96d37fcSHong Zhang       }
1395f96d37fcSHong Zhang     }
13964e938277SHong Zhang     nnz = lnk[0];
1397f96d37fcSHong Zhang 
1398f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1399f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1400f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1401f96d37fcSHong Zhang       nspacedouble++;
1402f96d37fcSHong Zhang     }
1403f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
14044e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1405f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
14062205254eSKarl Rupp 
1407f96d37fcSHong Zhang     current_space->array           += nnz;
1408f96d37fcSHong Zhang     current_space->local_used      += nnz;
1409f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
14102205254eSKarl Rupp 
1411f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1412f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1413f96d37fcSHong Zhang   }
1414f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1415f96d37fcSHong Zhang 
1416854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1417f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
14182205254eSKarl Rupp 
1419118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1420f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1421d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
14224e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
14234e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1424f96d37fcSHong Zhang 
14251c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
14261c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
14271795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1428f96d37fcSHong Zhang 
1429f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1430f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
143133d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1432f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1433f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1434f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1435f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1436f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1437f96d37fcSHong Zhang     row  = i + rstart;
1438f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1439f96d37fcSHong Zhang     Jptr = bj + bi[i];
1440f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1441f96d37fcSHong Zhang   }
1442f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1443f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1444f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1445f96d37fcSHong Zhang 
1446f96d37fcSHong Zhang   merge->bi        = bi;
1447f96d37fcSHong Zhang   merge->bj        = bj;
1448f96d37fcSHong Zhang   merge->coi       = coi;
1449f96d37fcSHong Zhang   merge->coj       = coj;
1450f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1451f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1452f96d37fcSHong Zhang   merge->owners_co = owners_co;
1453f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1454f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1455f96d37fcSHong Zhang 
14566da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1457c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1458c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1459f96d37fcSHong Zhang 
1460f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1461f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1462f96d37fcSHong Zhang   c->ptap     = ptap;
14630298fd71SBarry Smith   ptap->api   = NULL;
14640298fd71SBarry Smith   ptap->apj   = NULL;
1465f96d37fcSHong Zhang   ptap->merge = merge;
1466d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1467d6ab1dc0SHong Zhang 
1468d6ab1dc0SHong Zhang   *C = Cmpi;
1469d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1470d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
147157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
147257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1473d6ab1dc0SHong Zhang   } else {
1474d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1475d6ab1dc0SHong Zhang   }
1476d6ab1dc0SHong Zhang #endif
1477d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1478d6ab1dc0SHong Zhang }
1479d6ab1dc0SHong Zhang 
1480d6ab1dc0SHong Zhang #undef __FUNCT__
14816da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14826da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1483d6ab1dc0SHong Zhang {
1484d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1485d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1486d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1487d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1488d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1489e745005bSHong Zhang   PetscInt            *adj;
1490e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1491e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1492d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1493ce94432eSBarry Smith   MPI_Comm            comm;
1494d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1495d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1496d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1497d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1498d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1499d6ab1dc0SHong Zhang   MPI_Status          *status;
1500d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1501d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1502a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1503d6ab1dc0SHong Zhang   Mat                 A_loc;
1504d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1505d6ab1dc0SHong Zhang 
1506d6ab1dc0SHong Zhang   PetscFunctionBegin;
1507ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1508d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1509d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1510d6ab1dc0SHong Zhang 
1511d6ab1dc0SHong Zhang   ptap  = c->ptap;
1512d6ab1dc0SHong Zhang   merge = ptap->merge;
1513d6ab1dc0SHong Zhang 
1514e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1515e745005bSHong Zhang   /*------------------------------------------*/
1516d6ab1dc0SHong Zhang   /* get data from symbolic products */
1517d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1518854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1519d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1520d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
15211795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1522d6ab1dc0SHong Zhang 
1523d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1524d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1525d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1526d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1527d6ab1dc0SHong Zhang   ai    = a_loc->i;
1528d6ab1dc0SHong Zhang   aj    = a_loc->j;
1529d6ab1dc0SHong Zhang 
1530d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1531d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1532d6ab1dc0SHong Zhang     adj = aj + ai[i];
1533d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1534d6ab1dc0SHong Zhang 
1535d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1536e745005bSHong Zhang     /*-------------------------------------------------------------*/
1537d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1538d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1539d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1540d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1541d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1542d6ab1dc0SHong Zhang       row = poJ[j];
1543d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1544d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1545e745005bSHong Zhang       /* perform sparse axpy */
1546e745005bSHong Zhang       nexta  = 0;
1547d6ab1dc0SHong Zhang       valtmp = pA[j];
1548e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1549e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1550e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1551e745005bSHong Zhang           nexta++;
1552d6ab1dc0SHong Zhang         }
1553e745005bSHong Zhang       }
1554e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1555d6ab1dc0SHong Zhang     }
1556d6ab1dc0SHong Zhang 
1557d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1558d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1559d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1560d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1561d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1562d6ab1dc0SHong Zhang       row = pdJ[j];
1563d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1564d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1565e745005bSHong Zhang       /* perform sparse axpy */
1566e745005bSHong Zhang       nexta  = 0;
1567d6ab1dc0SHong Zhang       valtmp = pA[j];
1568e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1569e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1570e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1571e745005bSHong Zhang           nexta++;
1572d6ab1dc0SHong Zhang         }
1573e745005bSHong Zhang       }
1574e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1575d6ab1dc0SHong Zhang     }
1576d6ab1dc0SHong Zhang   }
1577d6ab1dc0SHong Zhang 
1578d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1579d6ab1dc0SHong Zhang   /*------------------------------------*/
1580d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1581d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1582d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1583d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1584d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1585d6ab1dc0SHong Zhang 
1586dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1587d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1588d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1589d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1590d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1591d6ab1dc0SHong Zhang     k++;
1592d6ab1dc0SHong Zhang   }
1593d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1594d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1595d6ab1dc0SHong Zhang 
1596d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1597d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1598d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1599d6ab1dc0SHong Zhang 
1600d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1601d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1602dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1603d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1604e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1605d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1606d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1607d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1608d6ab1dc0SHong Zhang   }
1609d6ab1dc0SHong Zhang 
1610d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1611d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1612d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1613d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1614d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1615d6ab1dc0SHong Zhang     /* add received vals into ba */
1616d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1617d6ab1dc0SHong Zhang       /* i-th row */
1618d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1619d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1620d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1621d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1622d6ab1dc0SHong Zhang         nextcj = 0;
1623d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1624d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1625d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1626d6ab1dc0SHong Zhang           }
1627d6ab1dc0SHong Zhang         }
1628d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1629d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1630d6ab1dc0SHong Zhang       }
1631d6ab1dc0SHong Zhang     }
1632d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1633d6ab1dc0SHong Zhang   }
1634d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1635d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1636d6ab1dc0SHong Zhang 
1637d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1639d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1640d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1641d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1642d6ab1dc0SHong Zhang }
1643d6ab1dc0SHong Zhang 
1644c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
16452bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
16462bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1647d6ab1dc0SHong Zhang #undef __FUNCT__
16486da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
16496da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1650d6ab1dc0SHong Zhang {
1651d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1652d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1653d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16540298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1655d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1656d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1657d6ab1dc0SHong Zhang   PetscInt            nnz;
1658d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1659d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1660ce94432eSBarry Smith   MPI_Comm            comm;
1661d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1662d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1663d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1664d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1665d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1666d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1667d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1668d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1669d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1670d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1671c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1672d6ab1dc0SHong Zhang   PetscScalar         *vals;
1673d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1674d6ab1dc0SHong Zhang 
1675d6ab1dc0SHong Zhang   PetscFunctionBegin;
1676ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1677d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1678d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1679d6ab1dc0SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1680d6ab1dc0SHong Zhang   }
1681d6ab1dc0SHong Zhang 
1682d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1683d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1684d6ab1dc0SHong Zhang 
1685d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1686b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1687d6ab1dc0SHong Zhang 
1688d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1689d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16902205254eSKarl Rupp 
1691d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1692d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1693d6ab1dc0SHong Zhang   ai          = a_loc->i;
1694d6ab1dc0SHong Zhang   aj          = a_loc->j;
1695d6ab1dc0SHong Zhang 
1696d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1697d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1698d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1699d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1700d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1701d6ab1dc0SHong Zhang 
1702d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1703d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1704d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1705d6ab1dc0SHong Zhang 
1706d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1707d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1708d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1709854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1710d6ab1dc0SHong Zhang   coi[0] = 0;
1711d6ab1dc0SHong Zhang 
1712d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1713d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1714a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1715d6ab1dc0SHong Zhang   current_space = free_space;
171619f0e57aSHong Zhang 
1717d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1718d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1719c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1720d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1721d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1722d6ab1dc0SHong Zhang 
1723d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1724d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1725d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1726d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1727d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1728d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1729d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1730d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1731d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1732d6ab1dc0SHong Zhang     }
1733d6ab1dc0SHong Zhang     nnz = lnk[0];
1734d6ab1dc0SHong Zhang 
1735d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1736d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1737d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1738d6ab1dc0SHong Zhang       nspacedouble++;
1739d6ab1dc0SHong Zhang     }
1740d6ab1dc0SHong Zhang 
1741d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1742d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17432205254eSKarl Rupp 
1744d6ab1dc0SHong Zhang     current_space->array           += nnz;
1745d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1746d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17472205254eSKarl Rupp 
1748d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1749d6ab1dc0SHong Zhang   }
1750d6ab1dc0SHong Zhang 
1751854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1752d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17532205254eSKarl Rupp 
1754118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1755d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1756d6ab1dc0SHong Zhang 
1757d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1758d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1759d6ab1dc0SHong Zhang   /* determine row ownership */
1760b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17622205254eSKarl Rupp 
1763d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1764d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17652205254eSKarl Rupp 
1766d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1767d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1768d6ab1dc0SHong Zhang 
1769d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17701795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1771785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17722205254eSKarl Rupp 
1773d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1774d6ab1dc0SHong Zhang   merge->nsend = 0;
1775d6ab1dc0SHong Zhang 
1776854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1777d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1778d6ab1dc0SHong Zhang 
1779d6ab1dc0SHong Zhang   proc = 0;
1780d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1781d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1782d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1783d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1784d6ab1dc0SHong Zhang   }
1785d6ab1dc0SHong Zhang 
1786d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1787d6ab1dc0SHong Zhang   owners_co[0] = 0;
1788d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1789d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1790d6ab1dc0SHong Zhang     if (len_si[proc]) {
1791d6ab1dc0SHong Zhang       merge->nsend++;
1792d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1793d6ab1dc0SHong Zhang       len         += len_si[proc];
1794d6ab1dc0SHong Zhang     }
1795d6ab1dc0SHong Zhang   }
1796d6ab1dc0SHong Zhang 
1797d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17980298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1799d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1800d6ab1dc0SHong Zhang 
1801d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1802d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1803d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1804854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1805d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1806d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1807d6ab1dc0SHong Zhang     i    = owners_co[proc];
1808d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1809d6ab1dc0SHong Zhang     k++;
1810d6ab1dc0SHong Zhang   }
1811d6ab1dc0SHong Zhang 
1812d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1813785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1815c280ed6aSJed Brown     PetscMPIInt icompleted;
1816c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1817d6ab1dc0SHong Zhang   }
1818d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1820d6ab1dc0SHong Zhang 
1821d6ab1dc0SHong Zhang   /* send and recv coi */
1822d6ab1dc0SHong Zhang   /*-------------------*/
1823d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1825854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1826d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1827d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1828d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1829d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1830d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1831d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1832d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1833d6ab1dc0SHong Zhang     */
1834d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1835d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1836d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1837d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1838d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1839d6ab1dc0SHong Zhang     nrows       = 0;
1840d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1841d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1842d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1843d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1844d6ab1dc0SHong Zhang       nrows++;
1845d6ab1dc0SHong Zhang     }
1846d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang     k++;
1848d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1849d6ab1dc0SHong Zhang   }
1850d6ab1dc0SHong Zhang   i = merge->nrecv;
1851d6ab1dc0SHong Zhang   while (i--) {
1852c280ed6aSJed Brown     PetscMPIInt icompleted;
1853c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1854d6ab1dc0SHong Zhang   }
1855d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1856d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1857d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1859d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1860d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1862d6ab1dc0SHong Zhang 
1863d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1864d6ab1dc0SHong Zhang   /*------------------------------------------*/
1865d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1866854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1867d6ab1dc0SHong Zhang   bi[0] = 0;
1868d6ab1dc0SHong Zhang 
1869d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1870d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1871a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1872d6ab1dc0SHong Zhang   current_space = free_space;
1873d6ab1dc0SHong Zhang 
1874dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1876d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1877d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1878d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18792205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1880d6ab1dc0SHong Zhang   }
1881d6ab1dc0SHong Zhang 
1882d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1883d6ab1dc0SHong Zhang   rmax = 0;
1884d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1885d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1886d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1887d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1888d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1889d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1890d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1891d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1892d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1893d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1894d6ab1dc0SHong Zhang     }
1895d6ab1dc0SHong Zhang 
1896d6ab1dc0SHong Zhang     /* add received col data into lnk */
1897d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1898d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1899d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1900d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1901d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1902d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1903d6ab1dc0SHong Zhang       }
1904d6ab1dc0SHong Zhang     }
1905d6ab1dc0SHong Zhang     nnz = lnk[0];
1906d6ab1dc0SHong Zhang 
1907d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1908d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1909d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1910d6ab1dc0SHong Zhang       nspacedouble++;
1911d6ab1dc0SHong Zhang     }
1912d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1913d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1914d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
19152205254eSKarl Rupp 
1916d6ab1dc0SHong Zhang     current_space->array           += nnz;
1917d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1918d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
19192205254eSKarl Rupp 
1920d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1921d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1922d6ab1dc0SHong Zhang   }
19230d3441aeSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); //mpiexec -n 2 ./ex94 -f0 $D/tiny -f1 $D/tiny
1924d6ab1dc0SHong Zhang 
1925854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1926d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1927118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1928d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1929d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1930d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1931d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1932d6ab1dc0SHong Zhang 
1933d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1934d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19351795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1936d6ab1dc0SHong Zhang 
1937d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1938d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
193933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1940d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1941d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1942d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1943d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1944d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1945d6ab1dc0SHong Zhang     row  = i + rstart;
1946d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1947d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1948d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1949d6ab1dc0SHong Zhang   }
1950d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1951d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1952d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1953d6ab1dc0SHong Zhang 
1954d6ab1dc0SHong Zhang   merge->bi        = bi;
1955d6ab1dc0SHong Zhang   merge->bj        = bj;
1956d6ab1dc0SHong Zhang   merge->coi       = coi;
1957d6ab1dc0SHong Zhang   merge->coj       = coj;
1958d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1959d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1960d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1961d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1962d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1963d6ab1dc0SHong Zhang 
19646da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1965d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1966c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1967d6ab1dc0SHong Zhang 
1968d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1969d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19702205254eSKarl Rupp 
1971d6ab1dc0SHong Zhang   c->ptap     = ptap;
19720298fd71SBarry Smith   ptap->api   = NULL;
19730298fd71SBarry Smith   ptap->apj   = NULL;
1974d6ab1dc0SHong Zhang   ptap->merge = merge;
1975d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19760298fd71SBarry Smith   ptap->apa   = NULL;
1977f96d37fcSHong Zhang 
1978f96d37fcSHong Zhang   *C = Cmpi;
1979f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1980f96d37fcSHong Zhang   if (bi[pn] != 0) {
198157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
198257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1983f96d37fcSHong Zhang   } else {
1984f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1985f96d37fcSHong Zhang   }
1986f96d37fcSHong Zhang #endif
1987f96d37fcSHong Zhang   PetscFunctionReturn(0);
1988f96d37fcSHong Zhang }
1989