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