xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 9ce11a7cd8203bfcb72e5062b2420c4de75b86d7)
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 #undef __FUNCT__
17 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
18 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
19 {
20   PetscErrorCode ierr;
21   Mat_MPIAIJ     *a   =(Mat_MPIAIJ*)A->data;
22   Mat_PtAPMPI    *ptap=a->ptap;
23 
24   PetscFunctionBegin;
25   if (ptap) {
26     Mat_Merge_SeqsToMPI *merge=ptap->merge;
27     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
28     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
29     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
30     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
31     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
32 
33     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
34     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
35     ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
36 
37     if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
38     if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
39     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
40     if (merge) {
41       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
42       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
43       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
44       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
45       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
46       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
47       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
48       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
49       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
50       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
51       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
52       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
53       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
54       ierr = merge->destroy(A);CHKERRQ(ierr);
55       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
56     }
57     ierr = PetscFree(ptap);CHKERRQ(ierr);
58   }
59   PetscFunctionReturn(0);
60 }
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
64 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
65 {
66   PetscErrorCode      ierr;
67   Mat_MPIAIJ          *a     = (Mat_MPIAIJ*)A->data;
68   Mat_PtAPMPI         *ptap  = a->ptap;
69   Mat_Merge_SeqsToMPI *merge = ptap->merge;
70 
71   PetscFunctionBegin;
72   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
73 
74   (*M)->ops->destroy   = merge->destroy;
75   (*M)->ops->duplicate = merge->duplicate;
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
81 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
82 {
83   PetscErrorCode ierr;
84   PetscBool      newalg=PETSC_FALSE;
85 
86   PetscFunctionBegin;
87   ierr = PetscOptionsGetBool(NULL,"-matptap_new",&newalg,NULL);CHKERRQ(ierr);
88   if (scall == MAT_INITIAL_MATRIX) {
89     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
90     if (newalg) {
91       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(A,P,fill,C);CHKERRQ(ierr);
92     } else {
93       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
94     }
95     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
96   }
97   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
98   if (newalg) {
99     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_new(A,P,*C);CHKERRQ(ierr);
100   } else {
101     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
102   }
103   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_new"
109 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(Mat A,Mat P,PetscReal fill,Mat *C)
110 {
111   PetscErrorCode ierr;
112   Mat_PtAPMPI    *ptap;
113   Mat_MPIAIJ     *c,*p=(Mat_MPIAIJ*)P->data;
114   Mat            AP;
115 
116   PetscFunctionBegin;
117   ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); //rm later !!!
118   c    = (Mat_MPIAIJ*)(*C)->data;
119   ptap = c->ptap;
120   ptap->reuse = MAT_INITIAL_MATRIX;
121   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
122   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
123 
124   ierr = MatMatMult(A,P,MAT_INITIAL_MATRIX,2.0,&AP);CHKERRQ(ierr);
125   ierr = MatMPIAIJGetLocalMat(AP,MAT_INITIAL_MATRIX,&ptap->AP_loc);CHKERRQ(ierr);
126   ierr = MatDestroy(&AP);CHKERRQ(ierr);
127 
128   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
129   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
130 
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_new"
136 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C)
137 {
138   PetscErrorCode    ierr;
139   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
140   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
141   Mat_PtAPMPI       *ptap = c->ptap;
142   Mat               AP_loc,C_loc,C_oth;
143   PetscInt          i,rstart,rend,cm,ncols,row;
144   PetscMPIInt       rank;
145   MPI_Comm          comm;
146   const PetscInt    *cols;
147   const PetscScalar *vals;
148   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
149 
150   PetscFunctionBegin;
151   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
152   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
153 
154   ierr = MatZeroEntries(C);CHKERRQ(ierr);
155 
156   /* 1) get R = Pd^T,Ro = Po^T */
157   ierr = PetscTime(&t0);CHKERRQ(ierr);
158   ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
159   ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
160   ierr = PetscTime(&t1);CHKERRQ(ierr);
161   eR = t1 - t0;
162 
163   /* 2) get AP_loc */
164   AP_loc = ptap->AP_loc;
165   Mat_SeqAIJ    *ap=(Mat_SeqAIJ*)AP_loc->data;
166 
167   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
168   /*-----------------------------------------------------*/
169   if (ptap->reuse == MAT_INITIAL_MATRIX) {
170     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
171     ptap->reuse = MAT_REUSE_MATRIX;
172   } else { /* update numerical values of P_oth and P_loc */
173     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
174     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
175   }
176 
177   /* 2-2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
178   /*--------------------------------------------------------------*/
179   /* get data from symbolic products */
180   Mat_SeqAIJ  *p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
181   Mat_SeqAIJ  *p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
182   PetscInt    *api,*apj,am = A->rmap->n,j,col,apnz;
183   PetscScalar *apa = ptap->apa;
184 
185   api = ptap->api; apj = ptap->apj;
186   for (i=0; i<am; i++) {
187     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
188     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
189     apnz = api[i+1] - api[i];
190     for (j=0; j<apnz; j++) {
191       col = apj[j+api[i]];
192       ap->a[j+ap->i[i]] = apa[col];
193       apa[col] = 0.0;
194     }
195   }
196 
197   ierr = PetscTime(&t2);CHKERRQ(ierr);
198   eAP = t2 - t1;
199 
200   /* 3) C_loc = R*AP_loc, Co = Ro*AP_loc */
201   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
202   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
203   C_loc = ptap->C_loc;
204   C_oth = ptap->C_oth;
205   //printf("[%d] Co %d, %d\n", rank,Co->rmap->N,Co->cmap->N);
206   ierr = PetscTime(&t3);CHKERRQ(ierr);
207   eCseq = t3 - t2;
208 
209   /* add C_loc and Co to to C */
210   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
211 
212   /* C_loc -> C */
213   cm = C_loc->rmap->N;
214   Mat_SeqAIJ *c_seq;
215   c_seq = (Mat_SeqAIJ*)C_loc->data;
216   for (i=0; i<cm; i++) {
217     ncols = c_seq->i[i+1] - c_seq->i[i];
218     row = rstart + i;
219     cols = c_seq->j + c_seq->i[i];
220     vals = c_seq->a + c_seq->i[i];
221     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
222   }
223 
224   /* Co -> C, off-processor part */
225   //printf("[%d] p->B %d, %d\n",rank,p->B->rmap->N,p->B->cmap->N);
226   cm = C_oth->rmap->N;
227   c_seq = (Mat_SeqAIJ*)C_oth->data;
228   for (i=0; i<cm; i++) {
229     ncols = c_seq->i[i+1] - c_seq->i[i];
230     row = p->garray[i];
231     cols = c_seq->j + c_seq->i[i];
232     vals = c_seq->a + c_seq->i[i];
233     //printf("[%d] row[%d] = %d\n",rank,i,row);
234     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
235   }
236   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238   ierr = PetscTime(&t4);CHKERRQ(ierr);
239   eCmpi = t4 - t3;
240 
241   ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
242   ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
243 
244   if (rank==1) {
245     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);
246   }
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
252 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
253 {
254   PetscErrorCode      ierr;
255   Mat                 Cmpi;
256   Mat_PtAPMPI         *ptap;
257   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
258   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
259   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
260   Mat_SeqAIJ          *p_loc,*p_oth;
261   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
262   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
263   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
264   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
265   PetscBT             lnkbt;
266   MPI_Comm            comm;
267   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
268   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
269   PetscInt            len,proc,*dnz,*onz,*owners;
270   PetscInt            nzi,*pti,*ptj;
271   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
272   MPI_Request         *swaits,*rwaits;
273   MPI_Status          *sstatus,rstatus;
274   Mat_Merge_SeqsToMPI *merge;
275   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
276   PetscReal           afill=1.0,afill_tmp;
277   PetscInt            rmax;
278 
279   PetscFunctionBegin;
280   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
281 
282   /* check if matrix local sizes are compatible */
283   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
284     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);
285   }
286   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
287     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);
288   }
289 
290   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
291   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
292 
293   /* create struct Mat_PtAPMPI and attached it to C later */
294   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
295   ierr        = PetscNew(&merge);CHKERRQ(ierr);
296   ptap->merge = merge;
297   ptap->reuse = MAT_INITIAL_MATRIX;
298 
299   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
300   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
301 
302   /* get P_loc by taking all local rows of P */
303   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
304 
305   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
306   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
307   pi_loc = p_loc->i; pj_loc = p_loc->j;
308   pi_oth = p_oth->i; pj_oth = p_oth->j;
309 
310   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
311   /*--------------------------------------------------------------------------*/
312   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
313   api[0] = 0;
314 
315   /* create and initialize a linked list */
316   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
317 
318   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
319   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
320 
321   current_space = free_space;
322 
323   for (i=0; i<am; i++) {
324     /* diagonal portion of A */
325     nzi = adi[i+1] - adi[i];
326     aj  = ad->j + adi[i];
327     for (j=0; j<nzi; j++) {
328       row  = aj[j];
329       pnz  = pi_loc[row+1] - pi_loc[row];
330       Jptr = pj_loc + pi_loc[row];
331       /* add non-zero cols of P into the sorted linked list lnk */
332       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
333     }
334     /* off-diagonal portion of A */
335     nzi = aoi[i+1] - aoi[i];
336     aj  = ao->j + aoi[i];
337     for (j=0; j<nzi; j++) {
338       row  = aj[j];
339       pnz  = pi_oth[row+1] - pi_oth[row];
340       Jptr = pj_oth + pi_oth[row];
341       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
342     }
343     apnz     = lnk[0];
344     api[i+1] = api[i] + apnz;
345     if (ap_rmax < apnz) ap_rmax = apnz;
346 
347     /* if free space is not available, double the total space in the list */
348     if (current_space->local_remaining<apnz) {
349       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
350       nspacedouble++;
351     }
352 
353     /* Copy data into free space, then initialize lnk */
354     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
355 
356     current_space->array           += apnz;
357     current_space->local_used      += apnz;
358     current_space->local_remaining -= apnz;
359   }
360 
361   /* Allocate space for apj, initialize apj, and */
362   /* destroy list of free space and other temporary array(s) */
363   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
364   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
365   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
366   if (afill_tmp > afill) afill = afill_tmp;
367 
368   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
369   /*-----------------------------------------------------------------*/
370   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
371 
372   /* then, compute symbolic Co = (p->B)^T*AP */
373   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
374                          >= (num of nonzero rows of C_seq) - pn */
375   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
376   coi[0] = 0;
377 
378   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
379   nnz           = fill*(poti[pon] + api[am]);
380   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
381   current_space = free_space;
382 
383   for (i=0; i<pon; i++) {
384     pnz = poti[i+1] - poti[i];
385     ptJ = potj + poti[i];
386     for (j=0; j<pnz; j++) {
387       row  = ptJ[j]; /* row of AP == col of Pot */
388       apnz = api[row+1] - api[row];
389       Jptr = apj + api[row];
390       /* add non-zero cols of AP into the sorted linked list lnk */
391       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
392     }
393     nnz = lnk[0];
394 
395     /* If free space is not available, double the total space in the list */
396     if (current_space->local_remaining<nnz) {
397       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
398       nspacedouble++;
399     }
400 
401     /* Copy data into free space, and zero out denserows */
402     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
403 
404     current_space->array           += nnz;
405     current_space->local_used      += nnz;
406     current_space->local_remaining -= nnz;
407 
408     coi[i+1] = coi[i] + nnz;
409   }
410 
411   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
412   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
413   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
414   if (afill_tmp > afill) afill = afill_tmp;
415   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
416 
417   /* (3) send j-array (coj) of Co to other processors */
418   /*--------------------------------------------------*/
419   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
420   len_s        = merge->len_s;
421   merge->nsend = 0;
422 
423 
424   /* determine row ownership */
425   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
426   merge->rowmap->n  = pn;
427   merge->rowmap->bs = 1;
428 
429   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
430   owners = merge->rowmap->range;
431 
432   /* determine the number of messages to send, their lengths */
433   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
434   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
435   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
436 
437   proc = 0;
438   for (i=0; i<pon; i++) {
439     while (prmap[i] >= owners[proc+1]) proc++;
440     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
441     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
442   }
443 
444   len          = 0; /* max length of buf_si[], see (4) */
445   owners_co[0] = 0;
446   for (proc=0; proc<size; proc++) {
447     owners_co[proc+1] = owners_co[proc] + len_si[proc];
448     if (len_s[proc]) {
449       merge->nsend++;
450       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
451       len         += len_si[proc];
452     }
453   }
454 
455   /* determine the number and length of messages to receive for coi and coj  */
456   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
457   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
458 
459   /* post the Irecv and Isend of coj */
460   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
461   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
462   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
463   for (proc=0, k=0; proc<size; proc++) {
464     if (!len_s[proc]) continue;
465     i    = owners_co[proc];
466     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
467     k++;
468   }
469 
470   /* receives and sends of coj are complete */
471   for (i=0; i<merge->nrecv; i++) {
472     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
473   }
474   ierr = PetscFree(rwaits);CHKERRQ(ierr);
475   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
476 
477   /* (4) send and recv coi */
478   /*-----------------------*/
479   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
480   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
481   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
482   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
483   for (proc=0,k=0; proc<size; proc++) {
484     if (!len_s[proc]) continue;
485     /* form outgoing message for i-structure:
486          buf_si[0]:                 nrows to be sent
487                [1:nrows]:           row index (global)
488                [nrows+1:2*nrows+1]: i-structure index
489     */
490     /*-------------------------------------------*/
491     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
492     buf_si_i    = buf_si + nrows+1;
493     buf_si[0]   = nrows;
494     buf_si_i[0] = 0;
495     nrows       = 0;
496     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
497       nzi = coi[i+1] - coi[i];
498       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
499       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
500       nrows++;
501     }
502     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
503     k++;
504     buf_si += len_si[proc];
505   }
506   i = merge->nrecv;
507   while (i--) {
508     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
509   }
510   ierr = PetscFree(rwaits);CHKERRQ(ierr);
511   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
512 
513   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
514   ierr = PetscFree(len_ri);CHKERRQ(ierr);
515   ierr = PetscFree(swaits);CHKERRQ(ierr);
516   ierr = PetscFree(buf_s);CHKERRQ(ierr);
517 
518   /* (5) compute the local portion of C (mpi mat) */
519   /*----------------------------------------------*/
520   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
521 
522   /* allocate pti array and free space for accumulating nonzero column info */
523   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
524   pti[0] = 0;
525 
526   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
527   nnz           = fill*(pi_loc[pm] + api[am]);
528   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
529   current_space = free_space;
530 
531   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
532   for (k=0; k<merge->nrecv; k++) {
533     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
534     nrows       = *buf_ri_k[k];
535     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
536     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
537   }
538   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
539   rmax = 0;
540   for (i=0; i<pn; i++) {
541     /* add pdt[i,:]*AP into lnk */
542     pnz = pdti[i+1] - pdti[i];
543     ptJ = pdtj + pdti[i];
544     for (j=0; j<pnz; j++) {
545       row  = ptJ[j];  /* row of AP == col of Pt */
546       apnz = api[row+1] - api[row];
547       Jptr = apj + api[row];
548       /* add non-zero cols of AP into the sorted linked list lnk */
549       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
550     }
551 
552     /* add received col data into lnk */
553     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
554       if (i == *nextrow[k]) { /* i-th row */
555         nzi  = *(nextci[k]+1) - *nextci[k];
556         Jptr = buf_rj[k] + *nextci[k];
557         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
558         nextrow[k]++; nextci[k]++;
559       }
560     }
561     nnz = lnk[0];
562 
563     /* if free space is not available, make more free space */
564     if (current_space->local_remaining<nnz) {
565       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
566       nspacedouble++;
567     }
568     /* copy data into free space, then initialize lnk */
569     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
570     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
571 
572     current_space->array           += nnz;
573     current_space->local_used      += nnz;
574     current_space->local_remaining -= nnz;
575 
576     pti[i+1] = pti[i] + nnz;
577     if (nnz > rmax) rmax = nnz;
578   }
579   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
580   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
581 
582   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
583   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
584   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
585   if (afill_tmp > afill) afill = afill_tmp;
586   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
587 
588   /* (6) create symbolic parallel matrix Cmpi */
589   /*------------------------------------------*/
590   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
591   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
592   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
593   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
594   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
595   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
596 
597   merge->bi        = pti;      /* Cseq->i */
598   merge->bj        = ptj;      /* Cseq->j */
599   merge->coi       = coi;      /* Co->i   */
600   merge->coj       = coj;      /* Co->j   */
601   merge->buf_ri    = buf_ri;
602   merge->buf_rj    = buf_rj;
603   merge->owners_co = owners_co;
604   merge->destroy   = Cmpi->ops->destroy;
605   merge->duplicate = Cmpi->ops->duplicate;
606 
607   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
608   Cmpi->assembled      = PETSC_FALSE;
609   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
610   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
611 
612   /* attach the supporting struct to Cmpi for reuse */
613   c           = (Mat_MPIAIJ*)Cmpi->data;
614   c->ptap     = ptap;
615   ptap->api   = api;
616   ptap->apj   = apj;
617   ptap->rmax  = ap_rmax;
618   *C          = Cmpi;
619 
620   /* flag 'scalable' determines which implementations to be used:
621        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
622        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
623   /* set default scalable */
624   ptap->scalable = PETSC_FALSE; //PETSC_TRUE;
625 
626   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
627   if (!ptap->scalable) {  /* Do dense axpy */
628     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
629   } else {
630     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
631   }
632 
633 #if defined(PETSC_USE_INFO)
634   if (pti[pn] != 0) {
635     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
636     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
637   } else {
638     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
639   }
640 #endif
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
646 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
647 {
648   PetscErrorCode      ierr;
649   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
650   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
651   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
652   Mat_SeqAIJ          *p_loc,*p_oth;
653   Mat_PtAPMPI         *ptap;
654   Mat_Merge_SeqsToMPI *merge;
655   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
656   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
657   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
658   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
659   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
660   MPI_Comm            comm;
661   PetscMPIInt         size,rank,taga,*len_s;
662   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
663   PetscInt            **buf_ri,**buf_rj;
664   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
665   MPI_Request         *s_waits,*r_waits;
666   MPI_Status          *status;
667   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
668   PetscInt            *api,*apj,*coi,*coj;
669   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
670   PetscBool           scalable;
671 #if defined(PTAP_PROFILE)
672   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
673 #endif
674 
675   PetscFunctionBegin;
676   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
677 #if defined(PTAP_PROFILE)
678   ierr = PetscTime(&t0);CHKERRQ(ierr);
679 #endif
680   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
681   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
682 
683   ptap = c->ptap;
684   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");
685   merge    = ptap->merge;
686   apa      = ptap->apa;
687   scalable = ptap->scalable;
688 
689   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
690   /*-----------------------------------------------------*/
691   if (ptap->reuse == MAT_INITIAL_MATRIX) {
692     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
693     ptap->reuse = MAT_REUSE_MATRIX;
694   } else { /* update numerical values of P_oth and P_loc */
695     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
696     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
697   }
698 #if defined(PTAP_PROFILE)
699   ierr = PetscTime(&t1);CHKERRQ(ierr);
700   eP = t1-t0;
701 #endif
702   /*
703   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
704          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
705          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
706          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
707    */
708 
709   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
710   /*--------------------------------------------------------------*/
711   /* get data from symbolic products */
712   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
713   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
714   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
715   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
716 
717   coi  = merge->coi; coj = merge->coj;
718   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
719 
720   bi     = merge->bi; bj = merge->bj;
721   owners = merge->rowmap->range;
722   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
723 
724   api = ptap->api; apj = ptap->apj;
725 
726   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
727     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
728 #if 0
729     /* ------ 10x slower -------------- */
730     /*==================================*/
731     Mat         R = ptap->R;
732     Mat_SeqAIJ  *r = (Mat_SeqAIJ*)R->data;
733     PetscInt    *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N;
734     PetscScalar *ra=r->a,tmp,cdense[pN];
735 
736     ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr);
737     for (i=0; i<cm; i++) { /* each row of C or R */
738       rnz = ri[i+1] - ri[i];
739 
740       for (j=0; j<rnz; j++) { /* each nz of R */
741         arow = rj[ri[i] + j];
742 
743         /* diagonal portion of A */
744         anz  = ad->i[arow+1] - ad->i[arow];
745         for (k=0; k<anz; k++) { /* each nz of Ad */
746           tmp  = ra[ri[i] + j]*ad->a[ad->i[arow] + k];
747           prow = ad->j[ad->i[arow] + k];
748           pnz  = pi_loc[prow+1] - pi_loc[prow];
749 
750           for (l=0; l<pnz; l++) { /* each nz of P_loc */
751             pcol = pj_loc[pi_loc[prow] + l];
752             cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l];
753           }
754         }
755 
756         /* off-diagonal portion of A */
757         anz  = ao->i[arow+1] - ao->i[arow];
758         for (k=0; k<anz; k++) { /* each nz of Ao */
759           tmp  = ra[ri[i] + j]*ao->a[ao->i[arow] + k];
760           prow = ao->j[ao->i[arow] + k];
761           pnz  = pi_oth[prow+1] - pi_oth[prow];
762 
763           for (l=0; l<pnz; l++) { /* each nz of P_oth */
764             pcol = pj_oth[pi_oth[prow] + l];
765             cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l];
766           }
767         }
768 
769       } //for (j=0; j<rnz; j++)
770 
771       /* copy cdense[] into ca; zero cdense[] */
772       cnz = bi[i+1] - bi[i];
773       cj  = bj + bi[i];
774       ca  = ba + bi[i];
775       for (j=0; j<cnz; j++) {
776         ca[j] += cdense[cj[j]];
777         cdense[cj[j]] = 0.0;
778       }
779 #if 0
780       if (rank == 0) {
781         printf("[%d] row %d: ",rank,i);
782         for (j=0; j<pN; j++) printf(" %g,",cdense[j]);
783         printf("\n");
784       }
785       for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[]
786 #endif
787     } //for (i=0; i<cm; i++) {
788 #endif
789 
790     //==========================================
791 
792     ierr = PetscTime(&t1);CHKERRQ(ierr);
793     for (i=0; i<am; i++) {
794 #if defined(PTAP_PROFILE)
795       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
796 #endif
797       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
798       /*------------------------------------------------------------*/
799       apJ = apj + api[i];
800 
801       /* diagonal portion of A */
802       anz = adi[i+1] - adi[i];
803       adj = ad->j + adi[i];
804       ada = ad->a + adi[i];
805       for (j=0; j<anz; j++) {
806         row = adj[j];
807         pnz = pi_loc[row+1] - pi_loc[row];
808         pj  = pj_loc + pi_loc[row];
809         pa  = pa_loc + pi_loc[row];
810 
811         /* perform dense axpy */
812         valtmp = ada[j];
813         for (k=0; k<pnz; k++) {
814           apa[pj[k]] += valtmp*pa[k];
815         }
816         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
817       }
818 
819       /* off-diagonal portion of A */
820       anz = aoi[i+1] - aoi[i];
821       aoj = ao->j + aoi[i];
822       aoa = ao->a + aoi[i];
823       for (j=0; j<anz; j++) {
824         row = aoj[j];
825         pnz = pi_oth[row+1] - pi_oth[row];
826         pj  = pj_oth + pi_oth[row];
827         pa  = pa_oth + pi_oth[row];
828 
829         /* perform dense axpy */
830         valtmp = aoa[j];
831         for (k=0; k<pnz; k++) {
832           apa[pj[k]] += valtmp*pa[k];
833         }
834         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
835       }
836 #if defined(PTAP_PROFILE)
837       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
838       et2_AP += t2_1 - t2_0;
839 #endif
840 
841       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
842       /*--------------------------------------------------------------*/
843       apnz = api[i+1] - api[i];
844       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
845       pnz = po->i[i+1] - po->i[i];
846       poJ = po->j + po->i[i];
847       pA  = po->a + po->i[i];
848       for (j=0; j<pnz; j++) {
849         row = poJ[j];
850         cnz = coi[row+1] - coi[row];
851         cj  = coj + coi[row];
852         ca  = coa + coi[row];
853         /* perform dense axpy */
854         valtmp = pA[j];
855         for (k=0; k<cnz; k++) {
856           ca[k] += valtmp*apa[cj[k]];
857         }
858         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
859       }
860 #if 1
861       /* put the value into Cd (diagonal part) */
862       pnz = pd->i[i+1] - pd->i[i];
863       pdJ = pd->j + pd->i[i];
864       pA  = pd->a + pd->i[i];
865       for (j=0; j<pnz; j++) {
866         row = pdJ[j];
867         cnz = bi[row+1] - bi[row];
868         cj  = bj + bi[row];
869         ca  = ba + bi[row];
870         /* perform dense axpy */
871         valtmp = pA[j];
872         for (k=0; k<cnz; k++) {
873           ca[k] += valtmp*apa[cj[k]];
874         }
875         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
876       }
877 #endif
878       /* zero the current row of A*P */
879       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
880 #if defined(PTAP_PROFILE)
881       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
882       ePtAP += t2_2 - t2_1;
883 #endif
884     }
885 
886     if (rank == 100) {
887     for (row=0; row<cm; row++) {
888       printf("[%d] row %d: ",rank,row);
889       cnz = bi[row+1] - bi[row];
890       for (j=0; j<cnz; j++) printf(" %g,",ba[bi[row]+j]);
891       printf("\n");
892     }
893     }
894 
895   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
896     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
897     /*-----------------------------------------------------------------------------------------*/
898     pA=pa_loc;
899     for (i=0; i<am; i++) {
900 #if defined(PTAP_PROFILE)
901       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
902 #endif
903       /* form i-th sparse row of A*P */
904       apnz = api[i+1] - api[i];
905       apJ  = apj + api[i];
906       /* diagonal portion of A */
907       anz = adi[i+1] - adi[i];
908       adj = ad->j + adi[i];
909       ada = ad->a + adi[i];
910       for (j=0; j<anz; j++) {
911         row    = adj[j];
912         pnz    = pi_loc[row+1] - pi_loc[row];
913         pj     = pj_loc + pi_loc[row];
914         pa     = pa_loc + pi_loc[row];
915         valtmp = ada[j];
916         nextp  = 0;
917         for (k=0; nextp<pnz; k++) {
918           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
919             apa[k] += valtmp*pa[nextp++];
920           }
921         }
922         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
923       }
924       /* off-diagonal portion of A */
925       anz = aoi[i+1] - aoi[i];
926       aoj = ao->j + aoi[i];
927       aoa = ao->a + aoi[i];
928       for (j=0; j<anz; j++) {
929         row    = aoj[j];
930         pnz    = pi_oth[row+1] - pi_oth[row];
931         pj     = pj_oth + pi_oth[row];
932         pa     = pa_oth + pi_oth[row];
933         valtmp = aoa[j];
934         nextp  = 0;
935         for (k=0; nextp<pnz; k++) {
936           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
937             apa[k] += valtmp*pa[nextp++];
938           }
939         }
940         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
941       }
942 #if defined(PTAP_PROFILE)
943       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
944       et2_AP += t2_1 - t2_0;
945 #endif
946 
947       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
948       /*--------------------------------------------------------------*/
949       pnz = pi_loc[i+1] - pi_loc[i];
950       pJ  = pj_loc + pi_loc[i];
951       for (j=0; j<pnz; j++) {
952         nextap = 0;
953         row    = pJ[j]; /* global index */
954         if (row < pcstart || row >=pcend) { /* put the value into Co */
955           row = *poJ;
956           cj  = coj + coi[row];
957           ca  = coa + coi[row]; poJ++;
958         } else {                            /* put the value into Cd */
959           row = *pdJ;
960           cj  = bj + bi[row];
961           ca  = ba + bi[row]; pdJ++;
962         }
963         valtmp = pA[j];
964         for (k=0; nextap<apnz; k++) {
965           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
966         }
967         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
968       }
969       pA += pnz;
970       /* zero the current row info for A*P */
971       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
972 #if defined(PTAP_PROFILE)
973       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
974       ePtAP += t2_2 - t2_1;
975 #endif
976     }
977   }
978 #if defined(PTAP_PROFILE)
979   ierr = PetscTime(&t2);CHKERRQ(ierr);
980 #endif
981 
982   /* 3) send and recv matrix values coa */
983   /*------------------------------------*/
984   buf_ri = merge->buf_ri;
985   buf_rj = merge->buf_rj;
986   len_s  = merge->len_s;
987   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
988   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
989 
990   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
991   for (proc=0,k=0; proc<size; proc++) {
992     if (!len_s[proc]) continue;
993     i    = merge->owners_co[proc];
994     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
995     k++;
996   }
997   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
998   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
999 
1000   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1001   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1002   ierr = PetscFree(coa);CHKERRQ(ierr);
1003 #if defined(PTAP_PROFILE)
1004   ierr = PetscTime(&t3);CHKERRQ(ierr);
1005 #endif
1006 
1007   /* 4) insert local Cseq and received values into Cmpi */
1008   /*------------------------------------------------------*/
1009   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1010   for (k=0; k<merge->nrecv; k++) {
1011     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1012     nrows       = *(buf_ri_k[k]);
1013     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1014     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1015   }
1016 
1017   for (i=0; i<cm; i++) {
1018     row  = owners[rank] + i; /* global row index of C_seq */
1019     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1020     ba_i = ba + bi[i];
1021     bnz  = bi[i+1] - bi[i];
1022     /* add received vals into ba */
1023     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1024       /* i-th row */
1025       if (i == *nextrow[k]) {
1026         cnz    = *(nextci[k]+1) - *nextci[k];
1027         cj     = buf_rj[k] + *(nextci[k]);
1028         ca     = abuf_r[k] + *(nextci[k]);
1029         nextcj = 0;
1030         for (j=0; nextcj<cnz; j++) {
1031           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1032             ba_i[j] += ca[nextcj++];
1033           }
1034         }
1035         nextrow[k]++; nextci[k]++;
1036         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1037       }
1038     }
1039     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1040   }
1041   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1042   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1043 
1044   ierr = PetscFree(ba);CHKERRQ(ierr);
1045   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1046   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1047   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1048 #if defined(PTAP_PROFILE)
1049   ierr = PetscTime(&t4);CHKERRQ(ierr);
1050   if (rank==1) {
1051     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);
1052   }
1053 #endif
1054   PetscFunctionReturn(0);
1055 }
1056