xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
1 
2 /*
3   Defines projective product routines where A is a MPIAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <../src/mat/impls/aij/mpi/mpiaij.h>
10 #include <petscbt.h>
11 #include <petsctime.h>
12 
13 /* #define PTAP_PROFILE */
14 
15 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
16 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
17 {
18   PetscErrorCode ierr;
19   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
20   Mat_PtAPMPI    *ptap=a->ptap;
21 
22   PetscFunctionBegin;
23   if (ptap) {
24     Mat_Merge_SeqsToMPI *merge=ptap->merge;
25     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
26     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
27     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
28     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
29     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
30     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
31     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
32     if (ptap->AP_loc) { /* used by alg_rap */
33       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
34       ierr = PetscFree(ap->i);CHKERRQ(ierr);
35       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
36       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
37     } else { /* used by alg_ptap */
38       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
39       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
40     }
41     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
42     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
43     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
44 
45     if (merge) { /* used by alg_ptap */
46       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
47       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
48       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
49       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
50       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
51       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
52       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
53       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
54       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
55       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
56       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
57       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
58       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
59       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
60     }
61 
62     ierr = ptap->destroy(A);CHKERRQ(ierr);
63     ierr = PetscFree(ptap);CHKERRQ(ierr);
64   }
65   PetscFunctionReturn(0);
66 }
67 
68 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
69 {
70   PetscErrorCode ierr;
71   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
72   Mat_PtAPMPI    *ptap  = a->ptap;
73 
74   PetscFunctionBegin;
75   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
76   (*M)->ops->destroy   = ptap->destroy;
77   (*M)->ops->duplicate = ptap->duplicate;
78   PetscFunctionReturn(0);
79 }
80 
81 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
82 {
83   PetscErrorCode ierr;
84   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
85   MPI_Comm       comm;
86 
87   PetscFunctionBegin;
88   /* check if matrix local sizes are compatible */
89   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
90   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
91   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
92 
93   ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
94   if (scall == MAT_INITIAL_MATRIX) {
95     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
96     if (rap) { /* do R=P^T locally, then C=R*A*P */
97       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
98     } else {       /* do P^T*A*P */
99       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
100     }
101     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
102   }
103   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
104   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
105   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 #if defined(PETSC_HAVE_HYPRE)
110 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
111 #endif
112 
113 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
114 {
115   PetscErrorCode    ierr;
116   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
117   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
118   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
119   Mat_PtAPMPI       *ptap = c->ptap;
120   Mat               AP_loc,C_loc,C_oth;
121   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
122   PetscScalar       *apa;
123   const PetscInt    *cols;
124   const PetscScalar *vals;
125 
126   PetscFunctionBegin;
127   ierr = MatZeroEntries(C);CHKERRQ(ierr);
128 
129   /* 1) get R = Pd^T,Ro = Po^T */
130   if (ptap->reuse == MAT_REUSE_MATRIX) {
131     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
132     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
133   }
134 
135   /* 2) get AP_loc */
136   AP_loc = ptap->AP_loc;
137   ap = (Mat_SeqAIJ*)AP_loc->data;
138 
139   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
140   /*-----------------------------------------------------*/
141   if (ptap->reuse == MAT_REUSE_MATRIX) {
142     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
143     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
144     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
145   }
146 
147   /* 2-2) compute numeric A_loc*P - dominating part */
148   /* ---------------------------------------------- */
149   /* get data from symbolic products */
150   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
151   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
152   api   = ap->i;
153   apj   = ap->j;
154   for (i=0; i<am; i++) {
155     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
156     apnz = api[i+1] - api[i];
157     apa = ap->a + api[i];
158     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
159     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
160     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
161   }
162 
163   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
164   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
165   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
166 
167   C_loc = ptap->C_loc;
168   C_oth = ptap->C_oth;
169 
170   /* add C_loc and Co to to C */
171   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
172 
173   /* C_loc -> C */
174   cm    = C_loc->rmap->N;
175   c_seq = (Mat_SeqAIJ*)C_loc->data;
176   cols = c_seq->j;
177   vals = c_seq->a;
178   for (i=0; i<cm; i++) {
179     ncols = c_seq->i[i+1] - c_seq->i[i];
180     row = rstart + i;
181     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
182     cols += ncols; vals += ncols;
183   }
184 
185   /* Co -> C, off-processor part */
186   cm = C_oth->rmap->N;
187   c_seq = (Mat_SeqAIJ*)C_oth->data;
188   cols = c_seq->j;
189   vals = c_seq->a;
190   for (i=0; i<cm; i++) {
191     ncols = c_seq->i[i+1] - c_seq->i[i];
192     row = p->garray[i];
193     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
194     cols += ncols; vals += ncols;
195   }
196   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198 
199   ptap->reuse = MAT_REUSE_MATRIX;
200   PetscFunctionReturn(0);
201 }
202 
203 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
204 {
205   PetscErrorCode      ierr;
206   Mat_PtAPMPI         *ptap;
207   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
208   MPI_Comm            comm;
209   PetscMPIInt         size,rank;
210   Mat                 Cmpi,P_loc,P_oth;
211   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
212   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
213   PetscInt            *lnk,i,k,pnz,row,nsend;
214   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
215   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
216   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
217   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
218   MPI_Request         *swaits,*rwaits;
219   MPI_Status          *sstatus,rstatus;
220   PetscLayout         rowmap;
221   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
222   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
223   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
224   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
225   PetscScalar         *apv;
226   PetscTable          ta;
227 #if defined(PETSC_USE_INFO)
228   PetscReal           apfill;
229 #endif
230 
231   PetscFunctionBegin;
232   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
233   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
234   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
235 
236   /* create symbolic parallel matrix Cmpi */
237   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
238   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
239 
240   /* create struct Mat_PtAPMPI and attached it to C later */
241   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
242   ptap->reuse = MAT_INITIAL_MATRIX;
243 
244   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
245   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
246   /* get P_loc by taking all local rows of P */
247   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
248 
249   ptap->P_loc = P_loc;
250   ptap->P_oth = P_oth;
251 
252   /* (0) compute Rd = Pd^T, Ro = Po^T  */
253   /* --------------------------------- */
254   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
255   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
256 
257   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
258   /* ----------------------------------------------------------------- */
259   p_loc  = (Mat_SeqAIJ*)P_loc->data;
260   p_oth  = (Mat_SeqAIJ*)P_oth->data;
261 
262   /* create and initialize a linked list */
263   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
264   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
265   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
266   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
267 
268   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
269 
270   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
271   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
272   current_space = free_space;
273   nspacedouble  = 0;
274 
275   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
276   api[0] = 0;
277   for (i=0; i<am; i++) {
278     /* diagonal portion: Ad[i,:]*P */
279     ai = ad->i; pi = p_loc->i;
280     nzi = ai[i+1] - ai[i];
281     aj  = ad->j + ai[i];
282     for (j=0; j<nzi; j++) {
283       row  = aj[j];
284       pnz  = pi[row+1] - pi[row];
285       Jptr = p_loc->j + pi[row];
286       /* add non-zero cols of P into the sorted linked list lnk */
287       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
288     }
289     /* off-diagonal portion: Ao[i,:]*P */
290     ai = ao->i; pi = p_oth->i;
291     nzi = ai[i+1] - ai[i];
292     aj  = ao->j + ai[i];
293     for (j=0; j<nzi; j++) {
294       row  = aj[j];
295       pnz  = pi[row+1] - pi[row];
296       Jptr = p_oth->j + pi[row];
297       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
298     }
299     apnz     = lnk[0];
300     api[i+1] = api[i] + apnz;
301 
302     /* if free space is not available, double the total space in the list */
303     if (current_space->local_remaining<apnz) {
304       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
305       nspacedouble++;
306     }
307 
308     /* Copy data into free space, then initialize lnk */
309     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
310 
311     current_space->array           += apnz;
312     current_space->local_used      += apnz;
313     current_space->local_remaining -= apnz;
314   }
315   /* Allocate space for apj and apv, initialize apj, and */
316   /* destroy list of free space and other temporary array(s) */
317   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
318   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
319   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
320 
321   /* Create AP_loc for reuse */
322   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
323 
324 #if defined(PETSC_USE_INFO)
325   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
326   ptap->AP_loc->info.mallocs           = nspacedouble;
327   ptap->AP_loc->info.fill_ratio_given  = fill;
328   ptap->AP_loc->info.fill_ratio_needed = apfill;
329 
330   if (api[am]) {
331     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
332     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
333   } else {
334     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
335   }
336 #endif
337 
338   /* (2-1) compute symbolic Co = Ro*AP_loc  */
339   /* ------------------------------------ */
340   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
341 
342   /* (3) send coj of C_oth to other processors  */
343   /* ------------------------------------------ */
344   /* determine row ownership */
345   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
346   rowmap->n  = pn;
347   rowmap->bs = 1;
348   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
349   owners = rowmap->range;
350 
351   /* determine the number of messages to send, their lengths */
352   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
353   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
354   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
355 
356   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
357   coi   = c_oth->i; coj = c_oth->j;
358   con   = ptap->C_oth->rmap->n;
359   proc  = 0;
360   for (i=0; i<con; i++) {
361     while (prmap[i] >= owners[proc+1]) proc++;
362     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
363     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
364   }
365 
366   len          = 0; /* max length of buf_si[], see (4) */
367   owners_co[0] = 0;
368   nsend        = 0;
369   for (proc=0; proc<size; proc++) {
370     owners_co[proc+1] = owners_co[proc] + len_si[proc];
371     if (len_s[proc]) {
372       nsend++;
373       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
374       len         += len_si[proc];
375     }
376   }
377 
378   /* determine the number and length of messages to receive for coi and coj  */
379   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
380   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
381 
382   /* post the Irecv and Isend of coj */
383   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
384   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
385   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
386   for (proc=0, k=0; proc<size; proc++) {
387     if (!len_s[proc]) continue;
388     i    = owners_co[proc];
389     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
390     k++;
391   }
392 
393   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
394   /* ---------------------------------------- */
395   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
396   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
397 
398   /* receives coj are complete */
399   for (i=0; i<nrecv; i++) {
400     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
401   }
402   ierr = PetscFree(rwaits);CHKERRQ(ierr);
403   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
404 
405   /* add received column indices into ta to update Crmax */
406   for (k=0; k<nrecv; k++) {/* k-th received message */
407     Jptr = buf_rj[k];
408     for (j=0; j<len_r[k]; j++) {
409       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
410     }
411   }
412   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
413   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
414 
415   /* (4) send and recv coi */
416   /*-----------------------*/
417   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
418   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
419   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
420   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
421   for (proc=0,k=0; proc<size; proc++) {
422     if (!len_s[proc]) continue;
423     /* form outgoing message for i-structure:
424          buf_si[0]:                 nrows to be sent
425                [1:nrows]:           row index (global)
426                [nrows+1:2*nrows+1]: i-structure index
427     */
428     /*-------------------------------------------*/
429     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
430     buf_si_i    = buf_si + nrows+1;
431     buf_si[0]   = nrows;
432     buf_si_i[0] = 0;
433     nrows       = 0;
434     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
435       nzi = coi[i+1] - coi[i];
436       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
437       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
438       nrows++;
439     }
440     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
441     k++;
442     buf_si += len_si[proc];
443   }
444   for (i=0; i<nrecv; i++) {
445     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
446   }
447   ierr = PetscFree(rwaits);CHKERRQ(ierr);
448   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
449 
450   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
451   ierr = PetscFree(len_ri);CHKERRQ(ierr);
452   ierr = PetscFree(swaits);CHKERRQ(ierr);
453   ierr = PetscFree(buf_s);CHKERRQ(ierr);
454 
455   /* (5) compute the local portion of Cmpi      */
456   /* ------------------------------------------ */
457   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
458   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
459   current_space = free_space;
460 
461   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
462   for (k=0; k<nrecv; k++) {
463     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
464     nrows       = *buf_ri_k[k];
465     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
466     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
467   }
468 
469   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
470   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
471   for (i=0; i<pn; i++) {
472     /* add C_loc into Cmpi */
473     nzi  = c_loc->i[i+1] - c_loc->i[i];
474     Jptr = c_loc->j + c_loc->i[i];
475     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
476 
477     /* add received col data into lnk */
478     for (k=0; k<nrecv; k++) { /* k-th received message */
479       if (i == *nextrow[k]) { /* i-th row */
480         nzi  = *(nextci[k]+1) - *nextci[k];
481         Jptr = buf_rj[k] + *nextci[k];
482         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
483         nextrow[k]++; nextci[k]++;
484       }
485     }
486     nzi = lnk[0];
487 
488     /* copy data into free space, then initialize lnk */
489     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
490     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
491   }
492   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
493   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
494   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
495 
496   /* local sizes and preallocation */
497   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
498   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
499   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
500   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
501 
502   /* members in merge */
503   ierr = PetscFree(id_r);CHKERRQ(ierr);
504   ierr = PetscFree(len_r);CHKERRQ(ierr);
505   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
506   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
507   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
508   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
509   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
510 
511   /* attach the supporting struct to Cmpi for reuse */
512   c = (Mat_MPIAIJ*)Cmpi->data;
513   c->ptap         = ptap;
514   ptap->duplicate = Cmpi->ops->duplicate;
515   ptap->destroy   = Cmpi->ops->destroy;
516 
517   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
518   Cmpi->assembled        = PETSC_FALSE;
519   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
520   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
521   *C                     = Cmpi;
522   PetscFunctionReturn(0);
523 }
524 
525 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
526 {
527   PetscErrorCode      ierr;
528   Mat_PtAPMPI         *ptap;
529   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
530   MPI_Comm            comm;
531   PetscMPIInt         size,rank;
532   Mat                 Cmpi;
533   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
534   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
535   PetscInt            *lnk,i,k,pnz,row,nsend;
536   PetscBT             lnkbt;
537   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
538   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
539   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
540   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
541   MPI_Request         *swaits,*rwaits;
542   MPI_Status          *sstatus,rstatus;
543   PetscLayout         rowmap;
544   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
545   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
546   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
547   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
548   PetscScalar         *apv;
549   PetscTable          ta;
550 #if defined(PETSC_HAVE_HYPRE)
551   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
552   PetscInt            nalg = 3;
553 #else
554   const char          *algTypes[2] = {"scalable","nonscalable"};
555   PetscInt            nalg = 2;
556 #endif
557   PetscInt            alg = 1; /* set default algorithm */
558 #if defined(PETSC_USE_INFO)
559   PetscReal           apfill;
560 #endif
561 #if defined(PTAP_PROFILE)
562   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
563 #endif
564   PetscBool           flg;
565 
566   PetscFunctionBegin;
567   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
568 
569   /* pick an algorithm */
570   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
571   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr);
572   ierr = PetscOptionsEnd();CHKERRQ(ierr);
573 
574   if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
575     MatInfo     Ainfo,Pinfo;
576     PetscInt    nz_local;
577     PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
578 
579     ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
580     ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
581     nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
582 
583     if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
584     ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
585 
586     if (alg_scalable) {
587       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
588       ierr = PetscInfo2(P,"Use scalable algorithm, pN %D, fill*nz_allocated %g\n",pN,fill*nz_local);CHKERRQ(ierr);
589     }
590   }
591   /* ierr = PetscInfo2(P,"Use algorithm %D, pN %D\n",alg,pN);CHKERRQ(ierr); */
592 
593   if (alg == 0) {
594     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
595     (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
596     PetscFunctionReturn(0);
597 
598 #if defined(PETSC_HAVE_HYPRE)
599   } else if (alg == 2) {
600     /* Use boomerAMGBuildCoarseOperator */
601     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
602     PetscFunctionReturn(0);
603 #endif
604   }
605 
606 #if defined(PTAP_PROFILE)
607   ierr = PetscTime(&t0);CHKERRQ(ierr);
608 #endif
609 
610   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
611   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
612 
613   /* create symbolic parallel matrix Cmpi */
614   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
615   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
616 
617   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
618   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
619 
620   /* create struct Mat_PtAPMPI and attached it to C later */
621   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
622   ptap->reuse = MAT_INITIAL_MATRIX;
623 
624   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
625   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
626   /* get P_loc by taking all local rows of P */
627   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
628 
629   /* (0) compute Rd = Pd^T, Ro = Po^T  */
630   /* --------------------------------- */
631   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
632   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
633 #if defined(PTAP_PROFILE)
634   ierr = PetscTime(&t11);CHKERRQ(ierr);
635 #endif
636 
637   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
638   /* ----------------------------------------------------------------- */
639   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
640   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
641 
642   /* create and initialize a linked list */
643   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
644   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
645   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
646   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
647   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
648 
649   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
650 
651   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
652   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
653   current_space = free_space;
654   nspacedouble  = 0;
655 
656   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
657   api[0] = 0;
658   for (i=0; i<am; i++) {
659     /* diagonal portion: Ad[i,:]*P */
660     ai = ad->i; pi = p_loc->i;
661     nzi = ai[i+1] - ai[i];
662     aj  = ad->j + ai[i];
663     for (j=0; j<nzi; j++) {
664       row  = aj[j];
665       pnz  = pi[row+1] - pi[row];
666       Jptr = p_loc->j + pi[row];
667       /* add non-zero cols of P into the sorted linked list lnk */
668       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
669     }
670     /* off-diagonal portion: Ao[i,:]*P */
671     ai = ao->i; pi = p_oth->i;
672     nzi = ai[i+1] - ai[i];
673     aj  = ao->j + ai[i];
674     for (j=0; j<nzi; j++) {
675       row  = aj[j];
676       pnz  = pi[row+1] - pi[row];
677       Jptr = p_oth->j + pi[row];
678       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
679     }
680     apnz     = lnk[0];
681     api[i+1] = api[i] + apnz;
682     if (ap_rmax < apnz) ap_rmax = apnz;
683 
684     /* if free space is not available, double the total space in the list */
685     if (current_space->local_remaining<apnz) {
686       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
687       nspacedouble++;
688     }
689 
690     /* Copy data into free space, then initialize lnk */
691     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
692 
693     current_space->array           += apnz;
694     current_space->local_used      += apnz;
695     current_space->local_remaining -= apnz;
696   }
697   /* Allocate space for apj and apv, initialize apj, and */
698   /* destroy list of free space and other temporary array(s) */
699   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
700   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
701   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
702 
703   /* Create AP_loc for reuse */
704   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
705 
706 #if defined(PETSC_USE_INFO)
707   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
708   ptap->AP_loc->info.mallocs           = nspacedouble;
709   ptap->AP_loc->info.fill_ratio_given  = fill;
710   ptap->AP_loc->info.fill_ratio_needed = apfill;
711 
712   if (api[am]) {
713     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
714     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
715   } else {
716     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
717   }
718 #endif
719 
720 #if defined(PTAP_PROFILE)
721   ierr = PetscTime(&t12);CHKERRQ(ierr);
722 #endif
723 
724   /* (2-1) compute symbolic Co = Ro*AP_loc  */
725   /* ------------------------------------ */
726   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
727 #if defined(PTAP_PROFILE)
728   ierr = PetscTime(&t1);CHKERRQ(ierr);
729 #endif
730 
731   /* (3) send coj of C_oth to other processors  */
732   /* ------------------------------------------ */
733   /* determine row ownership */
734   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
735   rowmap->n  = pn;
736   rowmap->bs = 1;
737   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
738   owners = rowmap->range;
739 
740   /* determine the number of messages to send, their lengths */
741   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
742   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
743   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
744 
745   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
746   coi   = c_oth->i; coj = c_oth->j;
747   con   = ptap->C_oth->rmap->n;
748   proc  = 0;
749   for (i=0; i<con; i++) {
750     while (prmap[i] >= owners[proc+1]) proc++;
751     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
752     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
753   }
754 
755   len          = 0; /* max length of buf_si[], see (4) */
756   owners_co[0] = 0;
757   nsend        = 0;
758   for (proc=0; proc<size; proc++) {
759     owners_co[proc+1] = owners_co[proc] + len_si[proc];
760     if (len_s[proc]) {
761       nsend++;
762       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
763       len         += len_si[proc];
764     }
765   }
766 
767   /* determine the number and length of messages to receive for coi and coj  */
768   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
769   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
770 
771   /* post the Irecv and Isend of coj */
772   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
773   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
774   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
775   for (proc=0, k=0; proc<size; proc++) {
776     if (!len_s[proc]) continue;
777     i    = owners_co[proc];
778     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
779     k++;
780   }
781 
782   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
783   /* ---------------------------------------- */
784   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
785   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
786 
787   /* receives coj are complete */
788   for (i=0; i<nrecv; i++) {
789     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
790   }
791   ierr = PetscFree(rwaits);CHKERRQ(ierr);
792   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
793 
794   /* add received column indices into ta to update Crmax */
795   for (k=0; k<nrecv; k++) {/* k-th received message */
796     Jptr = buf_rj[k];
797     for (j=0; j<len_r[k]; j++) {
798       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
799     }
800   }
801   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
802   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
803 
804   /* (4) send and recv coi */
805   /*-----------------------*/
806   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
807   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
808   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
809   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
810   for (proc=0,k=0; proc<size; proc++) {
811     if (!len_s[proc]) continue;
812     /* form outgoing message for i-structure:
813          buf_si[0]:                 nrows to be sent
814                [1:nrows]:           row index (global)
815                [nrows+1:2*nrows+1]: i-structure index
816     */
817     /*-------------------------------------------*/
818     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
819     buf_si_i    = buf_si + nrows+1;
820     buf_si[0]   = nrows;
821     buf_si_i[0] = 0;
822     nrows       = 0;
823     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
824       nzi = coi[i+1] - coi[i];
825       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
826       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
827       nrows++;
828     }
829     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
830     k++;
831     buf_si += len_si[proc];
832   }
833   for (i=0; i<nrecv; i++) {
834     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
835   }
836   ierr = PetscFree(rwaits);CHKERRQ(ierr);
837   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
838 
839   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
840   ierr = PetscFree(len_ri);CHKERRQ(ierr);
841   ierr = PetscFree(swaits);CHKERRQ(ierr);
842   ierr = PetscFree(buf_s);CHKERRQ(ierr);
843 #if defined(PTAP_PROFILE)
844   ierr = PetscTime(&t2);CHKERRQ(ierr);
845 #endif
846   /* (5) compute the local portion of Cmpi      */
847   /* ------------------------------------------ */
848   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
849   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
850   current_space = free_space;
851 
852   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
853   for (k=0; k<nrecv; k++) {
854     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
855     nrows       = *buf_ri_k[k];
856     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
857     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
858   }
859 
860   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
861   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
862   for (i=0; i<pn; i++) {
863     /* add C_loc into Cmpi */
864     nzi  = c_loc->i[i+1] - c_loc->i[i];
865     Jptr = c_loc->j + c_loc->i[i];
866     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
867 
868     /* add received col data into lnk */
869     for (k=0; k<nrecv; k++) { /* k-th received message */
870       if (i == *nextrow[k]) { /* i-th row */
871         nzi  = *(nextci[k]+1) - *nextci[k];
872         Jptr = buf_rj[k] + *nextci[k];
873         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
874         nextrow[k]++; nextci[k]++;
875       }
876     }
877     nzi = lnk[0];
878 
879     /* copy data into free space, then initialize lnk */
880     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
881     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
882   }
883   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
884   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
885   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
886 #if defined(PTAP_PROFILE)
887   ierr = PetscTime(&t3);CHKERRQ(ierr);
888 #endif
889 
890   /* local sizes and preallocation */
891   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
892   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
893   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
894   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
895 
896   /* members in merge */
897   ierr = PetscFree(id_r);CHKERRQ(ierr);
898   ierr = PetscFree(len_r);CHKERRQ(ierr);
899   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
900   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
901   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
902   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
903   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
904 
905   /* attach the supporting struct to Cmpi for reuse */
906   c = (Mat_MPIAIJ*)Cmpi->data;
907   c->ptap         = ptap;
908   ptap->duplicate = Cmpi->ops->duplicate;
909   ptap->destroy   = Cmpi->ops->destroy;
910 
911   if (alg == 1) {
912     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
913   }
914 
915   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
916   Cmpi->assembled        = PETSC_FALSE;
917   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
918   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
919   *C                     = Cmpi;
920 
921 #if defined(PTAP_PROFILE)
922   ierr = PetscTime(&t4);CHKERRQ(ierr);
923   if (rank == 1) {
924     printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0);
925   }
926 #endif
927   PetscFunctionReturn(0);
928 }
929 
930 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
931 {
932   PetscErrorCode    ierr;
933   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
934   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
935   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
936   Mat_PtAPMPI       *ptap = c->ptap;
937   Mat               AP_loc,C_loc,C_oth;
938   PetscInt          i,rstart,rend,cm,ncols,row;
939   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
940   PetscScalar       *apa;
941   const PetscInt    *cols;
942   const PetscScalar *vals;
943 #if defined(PTAP_PROFILE)
944   PetscMPIInt       rank;
945   MPI_Comm          comm;
946   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
947 #endif
948 
949   PetscFunctionBegin;
950   ierr = MatZeroEntries(C);CHKERRQ(ierr);
951 
952   /* 1) get R = Pd^T,Ro = Po^T */
953 #if defined(PTAP_PROFILE)
954   ierr = PetscTime(&t0);CHKERRQ(ierr);
955 #endif
956   if (ptap->reuse == MAT_REUSE_MATRIX) {
957     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
958     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
959   }
960 #if defined(PTAP_PROFILE)
961   ierr = PetscTime(&t1);CHKERRQ(ierr);
962   eR = t1 - t0;
963 #endif
964 
965   /* 2) get AP_loc */
966   AP_loc = ptap->AP_loc;
967   ap = (Mat_SeqAIJ*)AP_loc->data;
968 
969   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
970   /*-----------------------------------------------------*/
971   if (ptap->reuse == MAT_REUSE_MATRIX) {
972     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
973     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
974     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
975   }
976 
977   /* 2-2) compute numeric A_loc*P - dominating part */
978   /* ---------------------------------------------- */
979   /* get data from symbolic products */
980   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
981   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
982   apa   = ptap->apa;
983   api   = ap->i;
984   apj   = ap->j;
985   for (i=0; i<am; i++) {
986     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
987     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
988     apnz = api[i+1] - api[i];
989     for (j=0; j<apnz; j++) {
990       col = apj[j+api[i]];
991       ap->a[j+ap->i[i]] = apa[col];
992       apa[col] = 0.0;
993     }
994     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
995   }
996 #if defined(PTAP_PROFILE)
997   ierr = PetscTime(&t2);CHKERRQ(ierr);
998   eAP = t2 - t1;
999 #endif
1000 
1001   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1002   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1003   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
1004   C_loc = ptap->C_loc;
1005   C_oth = ptap->C_oth;
1006 #if defined(PTAP_PROFILE)
1007   ierr = PetscTime(&t3);CHKERRQ(ierr);
1008   eCseq = t3 - t2;
1009 #endif
1010 
1011   /* add C_loc and Co to to C */
1012   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1013 
1014   /* C_loc -> C */
1015   cm    = C_loc->rmap->N;
1016   c_seq = (Mat_SeqAIJ*)C_loc->data;
1017   cols = c_seq->j;
1018   vals = c_seq->a;
1019   for (i=0; i<cm; i++) {
1020     ncols = c_seq->i[i+1] - c_seq->i[i];
1021     row = rstart + i;
1022     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1023     cols += ncols; vals += ncols;
1024   }
1025 
1026   /* Co -> C, off-processor part */
1027   cm = C_oth->rmap->N;
1028   c_seq = (Mat_SeqAIJ*)C_oth->data;
1029   cols = c_seq->j;
1030   vals = c_seq->a;
1031   for (i=0; i<cm; i++) {
1032     ncols = c_seq->i[i+1] - c_seq->i[i];
1033     row = p->garray[i];
1034     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1035     cols += ncols; vals += ncols;
1036   }
1037   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1038   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1039 #if defined(PTAP_PROFILE)
1040   ierr = PetscTime(&t4);CHKERRQ(ierr);
1041   eCmpi = t4 - t3;
1042 
1043   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1044   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1045   if (rank==1) {
1046     ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr);
1047   }
1048 #endif
1049   ptap->reuse = MAT_REUSE_MATRIX;
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
1054 {
1055   PetscErrorCode      ierr;
1056   Mat                 Cmpi;
1057   Mat_PtAPMPI         *ptap;
1058   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1059   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1060   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1061   Mat_SeqAIJ          *p_loc,*p_oth;
1062   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
1063   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
1064   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1065   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
1066   PetscBT             lnkbt;
1067   MPI_Comm            comm;
1068   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
1069   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1070   PetscInt            len,proc,*dnz,*onz,*owners;
1071   PetscInt            nzi,*pti,*ptj;
1072   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1073   MPI_Request         *swaits,*rwaits;
1074   MPI_Status          *sstatus,rstatus;
1075   Mat_Merge_SeqsToMPI *merge;
1076   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1077   PetscReal           afill=1.0,afill_tmp;
1078 
1079   PetscFunctionBegin;
1080   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1081   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1082   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1083 
1084   /* create struct Mat_PtAPMPI and attached it to C later */
1085   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1086   ierr        = PetscNew(&merge);CHKERRQ(ierr);
1087   ptap->merge = merge;
1088   ptap->reuse = MAT_INITIAL_MATRIX;
1089 
1090   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1091   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1092 
1093   /* get P_loc by taking all local rows of P */
1094   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1095 
1096   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1097   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1098   pi_loc = p_loc->i; pj_loc = p_loc->j;
1099   pi_oth = p_oth->i; pj_oth = p_oth->j;
1100 
1101   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
1102   /*--------------------------------------------------------------------------*/
1103   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
1104   api[0] = 0;
1105 
1106   /* create and initialize a linked list */
1107   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1108 
1109   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
1110   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
1111   current_space = free_space;
1112 
1113   for (i=0; i<am; i++) {
1114     /* diagonal portion of A */
1115     nzi = adi[i+1] - adi[i];
1116     aj  = ad->j + adi[i];
1117     for (j=0; j<nzi; j++) {
1118       row  = aj[j];
1119       pnz  = pi_loc[row+1] - pi_loc[row];
1120       Jptr = pj_loc + pi_loc[row];
1121       /* add non-zero cols of P into the sorted linked list lnk */
1122       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1123     }
1124     /* off-diagonal portion of A */
1125     nzi = aoi[i+1] - aoi[i];
1126     aj  = ao->j + aoi[i];
1127     for (j=0; j<nzi; j++) {
1128       row  = aj[j];
1129       pnz  = pi_oth[row+1] - pi_oth[row];
1130       Jptr = pj_oth + pi_oth[row];
1131       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1132     }
1133     apnz     = lnk[0];
1134     api[i+1] = api[i] + apnz;
1135     if (ap_rmax < apnz) ap_rmax = apnz;
1136 
1137     /* if free space is not available, double the total space in the list */
1138     if (current_space->local_remaining<apnz) {
1139       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1140       nspacedouble++;
1141     }
1142 
1143     /* Copy data into free space, then initialize lnk */
1144     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1145 
1146     current_space->array           += apnz;
1147     current_space->local_used      += apnz;
1148     current_space->local_remaining -= apnz;
1149   }
1150 
1151   /* Allocate space for apj, initialize apj, and */
1152   /* destroy list of free space and other temporary array(s) */
1153   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
1154   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1155   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
1156   if (afill_tmp > afill) afill = afill_tmp;
1157 
1158   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
1159   /*-----------------------------------------------------------------*/
1160   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1161 
1162   /* then, compute symbolic Co = (p->B)^T*AP */
1163   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1164                          >= (num of nonzero rows of C_seq) - pn */
1165   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1166   coi[0] = 0;
1167 
1168   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
1169   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
1170   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1171   current_space = free_space;
1172 
1173   for (i=0; i<pon; i++) {
1174     pnz = poti[i+1] - poti[i];
1175     ptJ = potj + poti[i];
1176     for (j=0; j<pnz; j++) {
1177       row  = ptJ[j]; /* row of AP == col of Pot */
1178       apnz = api[row+1] - api[row];
1179       Jptr = apj + api[row];
1180       /* add non-zero cols of AP into the sorted linked list lnk */
1181       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1182     }
1183     nnz = lnk[0];
1184 
1185     /* If free space is not available, double the total space in the list */
1186     if (current_space->local_remaining<nnz) {
1187       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1188       nspacedouble++;
1189     }
1190 
1191     /* Copy data into free space, and zero out denserows */
1192     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1193 
1194     current_space->array           += nnz;
1195     current_space->local_used      += nnz;
1196     current_space->local_remaining -= nnz;
1197 
1198     coi[i+1] = coi[i] + nnz;
1199   }
1200 
1201   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
1202   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1203   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
1204   if (afill_tmp > afill) afill = afill_tmp;
1205   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1206 
1207   /* (3) send j-array (coj) of Co to other processors */
1208   /*--------------------------------------------------*/
1209   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
1210   len_s        = merge->len_s;
1211   merge->nsend = 0;
1212 
1213 
1214   /* determine row ownership */
1215   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1216   merge->rowmap->n  = pn;
1217   merge->rowmap->bs = 1;
1218 
1219   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1220   owners = merge->rowmap->range;
1221 
1222   /* determine the number of messages to send, their lengths */
1223   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
1224   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1225   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1226 
1227   proc = 0;
1228   for (i=0; i<pon; i++) {
1229     while (prmap[i] >= owners[proc+1]) proc++;
1230     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1231     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1232   }
1233 
1234   len          = 0; /* max length of buf_si[], see (4) */
1235   owners_co[0] = 0;
1236   for (proc=0; proc<size; proc++) {
1237     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1238     if (len_s[proc]) {
1239       merge->nsend++;
1240       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1241       len         += len_si[proc];
1242     }
1243   }
1244 
1245   /* determine the number and length of messages to receive for coi and coj  */
1246   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1247   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1248 
1249   /* post the Irecv and Isend of coj */
1250   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1251   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1252   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1253   for (proc=0, k=0; proc<size; proc++) {
1254     if (!len_s[proc]) continue;
1255     i    = owners_co[proc];
1256     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1257     k++;
1258   }
1259 
1260   /* receives and sends of coj are complete */
1261   for (i=0; i<merge->nrecv; i++) {
1262     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1263   }
1264   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1265   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1266 
1267   /* (4) send and recv coi */
1268   /*-----------------------*/
1269   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1270   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1271   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1272   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1273   for (proc=0,k=0; proc<size; proc++) {
1274     if (!len_s[proc]) continue;
1275     /* form outgoing message for i-structure:
1276          buf_si[0]:                 nrows to be sent
1277                [1:nrows]:           row index (global)
1278                [nrows+1:2*nrows+1]: i-structure index
1279     */
1280     /*-------------------------------------------*/
1281     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1282     buf_si_i    = buf_si + nrows+1;
1283     buf_si[0]   = nrows;
1284     buf_si_i[0] = 0;
1285     nrows       = 0;
1286     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1287       nzi = coi[i+1] - coi[i];
1288       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1289       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1290       nrows++;
1291     }
1292     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1293     k++;
1294     buf_si += len_si[proc];
1295   }
1296   i = merge->nrecv;
1297   while (i--) {
1298     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1299   }
1300   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1301   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1302 
1303   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
1304   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1305   ierr = PetscFree(swaits);CHKERRQ(ierr);
1306   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1307 
1308   /* (5) compute the local portion of C (mpi mat) */
1309   /*----------------------------------------------*/
1310   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1311 
1312   /* allocate pti array and free space for accumulating nonzero column info */
1313   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
1314   pti[0] = 0;
1315 
1316   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1317   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
1318   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1319   current_space = free_space;
1320 
1321   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1322   for (k=0; k<merge->nrecv; k++) {
1323     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1324     nrows       = *buf_ri_k[k];
1325     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1326     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1327   }
1328   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1329 
1330   for (i=0; i<pn; i++) {
1331     /* add pdt[i,:]*AP into lnk */
1332     pnz = pdti[i+1] - pdti[i];
1333     ptJ = pdtj + pdti[i];
1334     for (j=0; j<pnz; j++) {
1335       row  = ptJ[j];  /* row of AP == col of Pt */
1336       apnz = api[row+1] - api[row];
1337       Jptr = apj + api[row];
1338       /* add non-zero cols of AP into the sorted linked list lnk */
1339       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1340     }
1341 
1342     /* add received col data into lnk */
1343     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1344       if (i == *nextrow[k]) { /* i-th row */
1345         nzi  = *(nextci[k]+1) - *nextci[k];
1346         Jptr = buf_rj[k] + *nextci[k];
1347         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1348         nextrow[k]++; nextci[k]++;
1349       }
1350     }
1351     nnz = lnk[0];
1352 
1353     /* if free space is not available, make more free space */
1354     if (current_space->local_remaining<nnz) {
1355       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1356       nspacedouble++;
1357     }
1358     /* copy data into free space, then initialize lnk */
1359     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1360     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1361 
1362     current_space->array           += nnz;
1363     current_space->local_used      += nnz;
1364     current_space->local_remaining -= nnz;
1365 
1366     pti[i+1] = pti[i] + nnz;
1367   }
1368   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1369   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1370 
1371   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
1372   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
1373   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
1374   if (afill_tmp > afill) afill = afill_tmp;
1375   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1376 
1377   /* (6) create symbolic parallel matrix Cmpi */
1378   /*------------------------------------------*/
1379   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1380   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1381   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1382   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1383   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1384   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1385 
1386   merge->bi        = pti;      /* Cseq->i */
1387   merge->bj        = ptj;      /* Cseq->j */
1388   merge->coi       = coi;      /* Co->i   */
1389   merge->coj       = coj;      /* Co->j   */
1390   merge->buf_ri    = buf_ri;
1391   merge->buf_rj    = buf_rj;
1392   merge->owners_co = owners_co;
1393   merge->destroy   = Cmpi->ops->destroy;
1394 
1395   /* attach the supporting struct to Cmpi for reuse */
1396   c           = (Mat_MPIAIJ*)Cmpi->data;
1397   c->ptap     = ptap;
1398   ptap->api   = api;
1399   ptap->apj   = apj;
1400   ptap->duplicate = Cmpi->ops->duplicate;
1401   ptap->destroy   = Cmpi->ops->destroy;
1402 
1403   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1404   Cmpi->assembled        = PETSC_FALSE;
1405   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1406   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
1407   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
1408   *C                     = Cmpi;
1409 
1410   /* flag 'scalable' determines which implementations to be used:
1411        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
1412        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
1413   /* set default scalable */
1414   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
1415 
1416   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1417   if (!ptap->scalable) {  /* Do dense axpy */
1418     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1419   } else {
1420     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1421   }
1422 
1423 #if defined(PETSC_USE_INFO)
1424   if (pti[pn] != 0) {
1425     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1426     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1427   } else {
1428     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1429   }
1430 #endif
1431   PetscFunctionReturn(0);
1432 }
1433 
1434 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
1435 {
1436   PetscErrorCode      ierr;
1437   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1438   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1439   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1440   Mat_SeqAIJ          *p_loc,*p_oth;
1441   Mat_PtAPMPI         *ptap;
1442   Mat_Merge_SeqsToMPI *merge;
1443   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
1444   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
1445   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1446   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1447   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1448   MPI_Comm            comm;
1449   PetscMPIInt         size,rank,taga,*len_s;
1450   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1451   PetscInt            **buf_ri,**buf_rj;
1452   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1453   MPI_Request         *s_waits,*r_waits;
1454   MPI_Status          *status;
1455   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1456   PetscInt            *api,*apj,*coi,*coj;
1457   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
1458   PetscBool           scalable;
1459 #if defined(PTAP_PROFILE)
1460   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
1461 #endif
1462 
1463   PetscFunctionBegin;
1464   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1465 #if defined(PTAP_PROFILE)
1466   ierr = PetscTime(&t0);CHKERRQ(ierr);
1467 #endif
1468   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1469   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1470 
1471   ptap = c->ptap;
1472   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
1473   merge    = ptap->merge;
1474   apa      = ptap->apa;
1475   scalable = ptap->scalable;
1476 
1477   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
1478   /*-----------------------------------------------------*/
1479   if (ptap->reuse == MAT_INITIAL_MATRIX) {
1480     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1481     ptap->reuse = MAT_REUSE_MATRIX;
1482   } else { /* update numerical values of P_oth and P_loc */
1483     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1484     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1485   }
1486 #if defined(PTAP_PROFILE)
1487   ierr = PetscTime(&t1);CHKERRQ(ierr);
1488   eP = t1-t0;
1489 #endif
1490   /*
1491   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
1492          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
1493          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
1494          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
1495    */
1496 
1497   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1498   /*--------------------------------------------------------------*/
1499   /* get data from symbolic products */
1500   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1501   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1502   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
1503   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
1504 
1505   coi  = merge->coi; coj = merge->coj;
1506   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1507 
1508   bi     = merge->bi; bj = merge->bj;
1509   owners = merge->rowmap->range;
1510   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
1511 
1512   api = ptap->api; apj = ptap->apj;
1513 
1514   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1515     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
1516 #if defined(PTAP_PROFILE)
1517     ierr = PetscTime(&t1);CHKERRQ(ierr);
1518 #endif
1519     for (i=0; i<am; i++) {
1520 #if defined(PTAP_PROFILE)
1521       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1522 #endif
1523       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1524       /*------------------------------------------------------------*/
1525       apJ = apj + api[i];
1526 
1527       /* diagonal portion of A */
1528       anz = adi[i+1] - adi[i];
1529       adj = ad->j + adi[i];
1530       ada = ad->a + adi[i];
1531       for (j=0; j<anz; j++) {
1532         row = adj[j];
1533         pnz = pi_loc[row+1] - pi_loc[row];
1534         pj  = pj_loc + pi_loc[row];
1535         pa  = pa_loc + pi_loc[row];
1536 
1537         /* perform dense axpy */
1538         valtmp = ada[j];
1539         for (k=0; k<pnz; k++) {
1540           apa[pj[k]] += valtmp*pa[k];
1541         }
1542         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1543       }
1544 
1545       /* off-diagonal portion of A */
1546       anz = aoi[i+1] - aoi[i];
1547       aoj = ao->j + aoi[i];
1548       aoa = ao->a + aoi[i];
1549       for (j=0; j<anz; j++) {
1550         row = aoj[j];
1551         pnz = pi_oth[row+1] - pi_oth[row];
1552         pj  = pj_oth + pi_oth[row];
1553         pa  = pa_oth + pi_oth[row];
1554 
1555         /* perform dense axpy */
1556         valtmp = aoa[j];
1557         for (k=0; k<pnz; k++) {
1558           apa[pj[k]] += valtmp*pa[k];
1559         }
1560         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1561       }
1562 #if defined(PTAP_PROFILE)
1563       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1564       et2_AP += t2_1 - t2_0;
1565 #endif
1566 
1567       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1568       /*--------------------------------------------------------------*/
1569       apnz = api[i+1] - api[i];
1570       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1571       pnz = po->i[i+1] - po->i[i];
1572       poJ = po->j + po->i[i];
1573       pA  = po->a + po->i[i];
1574       for (j=0; j<pnz; j++) {
1575         row = poJ[j];
1576         cnz = coi[row+1] - coi[row];
1577         cj  = coj + coi[row];
1578         ca  = coa + coi[row];
1579         /* perform dense axpy */
1580         valtmp = pA[j];
1581         for (k=0; k<cnz; k++) {
1582           ca[k] += valtmp*apa[cj[k]];
1583         }
1584         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1585       }
1586       /* put the value into Cd (diagonal part) */
1587       pnz = pd->i[i+1] - pd->i[i];
1588       pdJ = pd->j + pd->i[i];
1589       pA  = pd->a + pd->i[i];
1590       for (j=0; j<pnz; j++) {
1591         row = pdJ[j];
1592         cnz = bi[row+1] - bi[row];
1593         cj  = bj + bi[row];
1594         ca  = ba + bi[row];
1595         /* perform dense axpy */
1596         valtmp = pA[j];
1597         for (k=0; k<cnz; k++) {
1598           ca[k] += valtmp*apa[cj[k]];
1599         }
1600         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1601       }
1602 
1603       /* zero the current row of A*P */
1604       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1605 #if defined(PTAP_PROFILE)
1606       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
1607       ePtAP += t2_2 - t2_1;
1608 #endif
1609     }
1610   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1611     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
1612     /*-----------------------------------------------------------------------------------------*/
1613     pA=pa_loc;
1614     for (i=0; i<am; i++) {
1615 #if defined(PTAP_PROFILE)
1616       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1617 #endif
1618       /* form i-th sparse row of A*P */
1619       apnz = api[i+1] - api[i];
1620       apJ  = apj + api[i];
1621       /* diagonal portion of A */
1622       anz = adi[i+1] - adi[i];
1623       adj = ad->j + adi[i];
1624       ada = ad->a + adi[i];
1625       for (j=0; j<anz; j++) {
1626         row    = adj[j];
1627         pnz    = pi_loc[row+1] - pi_loc[row];
1628         pj     = pj_loc + pi_loc[row];
1629         pa     = pa_loc + pi_loc[row];
1630         valtmp = ada[j];
1631         nextp  = 0;
1632         for (k=0; nextp<pnz; k++) {
1633           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1634             apa[k] += valtmp*pa[nextp++];
1635           }
1636         }
1637         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1638       }
1639       /* off-diagonal portion of A */
1640       anz = aoi[i+1] - aoi[i];
1641       aoj = ao->j + aoi[i];
1642       aoa = ao->a + aoi[i];
1643       for (j=0; j<anz; j++) {
1644         row    = aoj[j];
1645         pnz    = pi_oth[row+1] - pi_oth[row];
1646         pj     = pj_oth + pi_oth[row];
1647         pa     = pa_oth + pi_oth[row];
1648         valtmp = aoa[j];
1649         nextp  = 0;
1650         for (k=0; nextp<pnz; k++) {
1651           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1652             apa[k] += valtmp*pa[nextp++];
1653           }
1654         }
1655         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1656       }
1657 #if defined(PTAP_PROFILE)
1658       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1659       et2_AP += t2_1 - t2_0;
1660 #endif
1661 
1662       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1663       /*--------------------------------------------------------------*/
1664       pnz = pi_loc[i+1] - pi_loc[i];
1665       pJ  = pj_loc + pi_loc[i];
1666       for (j=0; j<pnz; j++) {
1667         nextap = 0;
1668         row    = pJ[j]; /* global index */
1669         if (row < pcstart || row >=pcend) { /* put the value into Co */
1670           row = *poJ;
1671           cj  = coj + coi[row];
1672           ca  = coa + coi[row]; poJ++;
1673         } else {                            /* put the value into Cd */
1674           row = *pdJ;
1675           cj  = bj + bi[row];
1676           ca  = ba + bi[row]; pdJ++;
1677         }
1678         valtmp = pA[j];
1679         for (k=0; nextap<apnz; k++) {
1680           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
1681         }
1682         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1683       }
1684       pA += pnz;
1685       /* zero the current row info for A*P */
1686       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1687 #if defined(PTAP_PROFILE)
1688       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
1689       ePtAP += t2_2 - t2_1;
1690 #endif
1691     }
1692   }
1693 #if defined(PTAP_PROFILE)
1694   ierr = PetscTime(&t2);CHKERRQ(ierr);
1695 #endif
1696 
1697   /* 3) send and recv matrix values coa */
1698   /*------------------------------------*/
1699   buf_ri = merge->buf_ri;
1700   buf_rj = merge->buf_rj;
1701   len_s  = merge->len_s;
1702   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1703   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1704 
1705   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1706   for (proc=0,k=0; proc<size; proc++) {
1707     if (!len_s[proc]) continue;
1708     i    = merge->owners_co[proc];
1709     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1710     k++;
1711   }
1712   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1713   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1714 
1715   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1716   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1717   ierr = PetscFree(coa);CHKERRQ(ierr);
1718 #if defined(PTAP_PROFILE)
1719   ierr = PetscTime(&t3);CHKERRQ(ierr);
1720 #endif
1721 
1722   /* 4) insert local Cseq and received values into Cmpi */
1723   /*------------------------------------------------------*/
1724   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1725   for (k=0; k<merge->nrecv; k++) {
1726     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1727     nrows       = *(buf_ri_k[k]);
1728     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1729     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1730   }
1731 
1732   for (i=0; i<cm; i++) {
1733     row  = owners[rank] + i; /* global row index of C_seq */
1734     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1735     ba_i = ba + bi[i];
1736     bnz  = bi[i+1] - bi[i];
1737     /* add received vals into ba */
1738     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1739       /* i-th row */
1740       if (i == *nextrow[k]) {
1741         cnz    = *(nextci[k]+1) - *nextci[k];
1742         cj     = buf_rj[k] + *(nextci[k]);
1743         ca     = abuf_r[k] + *(nextci[k]);
1744         nextcj = 0;
1745         for (j=0; nextcj<cnz; j++) {
1746           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1747             ba_i[j] += ca[nextcj++];
1748           }
1749         }
1750         nextrow[k]++; nextci[k]++;
1751         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1752       }
1753     }
1754     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1755   }
1756   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1757   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1758 
1759   ierr = PetscFree(ba);CHKERRQ(ierr);
1760   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1761   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1762   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1763 #if defined(PTAP_PROFILE)
1764   ierr = PetscTime(&t4);CHKERRQ(ierr);
1765   if (rank==1) {
1766     ierr = PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
1767   }
1768 #endif
1769   PetscFunctionReturn(0);
1770 }
1771