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