xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 68ae421e09e7e0a9b8058c8fd6e7192caa558d65)
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 
12 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
13 #undef __FUNCT__
14 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
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     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
30     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
31     if (ptap->apa){ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
32     if (merge) {
33       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
34       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
35       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
36       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
37       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
38       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
39       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
40       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
41       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
42       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
43       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
44       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
45       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
46       ierr = merge->destroy(A);CHKERRQ(ierr);
47       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
48     }
49     ierr = PetscFree(ptap);CHKERRQ(ierr);
50   }
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
56 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
57 {
58   PetscErrorCode       ierr;
59   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
60   Mat_PtAPMPI          *ptap = a->ptap;
61   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
62 
63   PetscFunctionBegin;
64   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
65   (*M)->ops->destroy   = merge->destroy;
66   (*M)->ops->duplicate = merge->duplicate;
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
72 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
73 {
74   PetscErrorCode ierr;
75 
76   PetscFunctionBegin;
77   if (scall == MAT_INITIAL_MATRIX){
78     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
79     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
80     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
81   }
82   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
83   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
84   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
90 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
91 {
92   PetscErrorCode       ierr;
93   Mat                  Cmpi;
94   Mat_PtAPMPI          *ptap;
95   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
96   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
97   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
98   Mat_SeqAIJ           *p_loc,*p_oth;
99   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
100   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
101   PetscInt             *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
102   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
103   PetscBT              lnkbt;
104   MPI_Comm             comm=((PetscObject)A)->comm;
105   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
106   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
107   PetscInt             len,proc,*dnz,*onz,*owners;
108   PetscInt             nzi,*bi,*bj;
109   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
110   MPI_Request          *swaits,*rwaits;
111   MPI_Status           *sstatus,rstatus;
112   Mat_Merge_SeqsToMPI  *merge;
113   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
114   PetscReal            afill=1.0,afill_tmp;
115   PetscInt             rmax;
116 
117   PetscFunctionBegin;
118   /* check if matrix local sizes are compatible */
119   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
120     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);
121   }
122   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
123     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);
124   }
125 
126   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
127   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
128 
129   /* create struct Mat_PtAPMPI and attached it to C later */
130   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
131   ptap->reuse = MAT_INITIAL_MATRIX;
132 
133   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
134   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
135 
136   /* get P_loc by taking all local rows of P */
137   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
138 
139   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
140   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
141   pi_loc = p_loc->i; pj_loc = p_loc->j;
142   pi_oth = p_oth->i; pj_oth = p_oth->j;
143 
144   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
145   /*-------------------------------------------------------------------*/
146   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
147   api[0] = 0;
148 
149   /* create and initialize a linked list */
150   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
151 
152   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
153   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
154   current_space = free_space;
155 
156   for (i=0; i<am; i++) {
157     /* diagonal portion of A */
158     nzi = adi[i+1] - adi[i];
159     aj  = ad->j + adi[i];
160     for (j=0; j<nzi; j++){
161       row  = aj[j];
162       pnz  = pi_loc[row+1] - pi_loc[row];
163       Jptr = pj_loc + pi_loc[row];
164       /* add non-zero cols of P into the sorted linked list lnk */
165       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
166     }
167     /* off-diagonal portion of A */
168     nzi = aoi[i+1] - aoi[i];
169     aj  = ao->j + aoi[i];
170     for (j=0; j<nzi; j++){
171       row = aj[j];
172       pnz = pi_oth[row+1] - pi_oth[row];
173       Jptr = pj_oth + pi_oth[row];
174       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
175     }
176     apnz = lnk[0];
177     api[i+1] = api[i] + apnz;
178     if (ap_rmax < apnz) ap_rmax = apnz;
179 
180     /* if free space is not available, double the total space in the list */
181     if (current_space->local_remaining<apnz) {
182       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
183       nspacedouble++;
184     }
185 
186     /* Copy data into free space, then initialize lnk */
187     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
188     current_space->array           += apnz;
189     current_space->local_used      += apnz;
190     current_space->local_remaining -= apnz;
191   }
192 
193   /* Allocate space for apj, initialize apj, and */
194   /* destroy list of free space and other temporary array(s) */
195   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
196   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
197   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
198   if (afill_tmp > afill) afill = afill_tmp;
199 
200   /* determine symbolic Co=(p->B)^T*AP - send to others */
201   /*----------------------------------------------------*/
202   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
203 
204   /* then, compute symbolic Co = (p->B)^T*AP */
205   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
206                          >= (num of nonzero rows of C_seq) - pn */
207   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
208   coi[0] = 0;
209 
210   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
211   nnz           = fill*(poti[pon] + api[am]);
212   ierr          = PetscFreeSpaceGet(nnz,&free_space);
213   current_space = free_space;
214 
215   for (i=0; i<pon; i++) {
216     pnz = poti[i+1] - poti[i];
217     ptJ = potj + poti[i];
218     for (j=0; j<pnz; j++){
219       row  = ptJ[j]; /* row of AP == col of Pot */
220       apnz = api[row+1] - api[row];
221       Jptr = apj + api[row];
222       /* add non-zero cols of AP into the sorted linked list lnk */
223       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
224     }
225     nnz = lnk[0];
226 
227     /* If free space is not available, double the total space in the list */
228     if (current_space->local_remaining<nnz) {
229       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
230       nspacedouble++;
231     }
232 
233     /* Copy data into free space, and zero out denserows */
234     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
235     current_space->array           += nnz;
236     current_space->local_used      += nnz;
237     current_space->local_remaining -= nnz;
238     coi[i+1] = coi[i] + nnz;
239   }
240   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
241   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
242   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
243   if (afill_tmp > afill) afill = afill_tmp;
244   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
245 
246   /* send j-array (coj) of Co to other processors */
247   /*----------------------------------------------*/
248   /* determine row ownership */
249   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
250   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
251   merge->rowmap->n  = pn;
252   merge->rowmap->bs = 1;
253   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
254   owners = merge->rowmap->range;
255 
256   /* determine the number of messages to send, their lengths */
257   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
258   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
259   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
260   len_s = merge->len_s;
261   merge->nsend = 0;
262 
263   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
264   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
265 
266   proc = 0;
267   for (i=0; i<pon; i++){
268     while (prmap[i] >= owners[proc+1]) proc++;
269     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
270     len_s[proc] += coi[i+1] - coi[i];
271   }
272 
273   len   = 0;  /* max length of buf_si[] */
274   owners_co[0] = 0;
275   for (proc=0; proc<size; proc++){
276     owners_co[proc+1] = owners_co[proc] + len_si[proc];
277     if (len_si[proc]){
278       merge->nsend++;
279       len_si[proc] = 2*(len_si[proc] + 1);
280       len += len_si[proc];
281     }
282   }
283 
284   /* determine the number and length of messages to receive for coi and coj  */
285   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
286   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
287 
288   /* post the Irecv and Isend of coj */
289   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
290   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
291   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
292   for (proc=0, k=0; proc<size; proc++){
293     if (!len_s[proc]) continue;
294     i = owners_co[proc];
295     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
296     k++;
297   }
298 
299   /* receives and sends of coj are complete */
300   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
301   for (i=0; i<merge->nrecv; i++){
302     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
303   }
304   ierr = PetscFree(rwaits);CHKERRQ(ierr);
305   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
306 
307   /* send and recv coi */
308   /*-------------------*/
309   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
310   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
311   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
312   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
313   for (proc=0,k=0; proc<size; proc++){
314     if (!len_s[proc]) continue;
315     /* form outgoing message for i-structure:
316          buf_si[0]:                 nrows to be sent
317                [1:nrows]:           row index (global)
318                [nrows+1:2*nrows+1]: i-structure index
319     */
320     /*-------------------------------------------*/
321     nrows = len_si[proc]/2 - 1;
322     buf_si_i    = buf_si + nrows+1;
323     buf_si[0]   = nrows;
324     buf_si_i[0] = 0;
325     nrows = 0;
326     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
327       nzi = coi[i+1] - coi[i];
328       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
329       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
330       nrows++;
331     }
332     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
333     k++;
334     buf_si += len_si[proc];
335   }
336   i = merge->nrecv;
337   while (i--) {
338     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
339   }
340   ierr = PetscFree(rwaits);CHKERRQ(ierr);
341   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
342 
343   ierr = PetscFree(len_si);CHKERRQ(ierr);
344   ierr = PetscFree(len_ri);CHKERRQ(ierr);
345   ierr = PetscFree(swaits);CHKERRQ(ierr);
346   ierr = PetscFree(sstatus);CHKERRQ(ierr);
347   ierr = PetscFree(buf_s);CHKERRQ(ierr);
348 
349   /* compute the local portion of C (mpi mat) */
350   /*------------------------------------------*/
351   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
352 
353   /* allocate bi array and free space for accumulating nonzero column info */
354   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
355   bi[0] = 0;
356 
357   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
358   nnz           = fill*(pi_loc[pm] + api[am]);
359   ierr          = PetscFreeSpaceGet(nnz,&free_space);
360   current_space = free_space;
361 
362   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
363   for (k=0; k<merge->nrecv; k++){
364     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
365     nrows       = *buf_ri_k[k];
366     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
367     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
368   }
369   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
370   rmax = 0;
371   for (i=0; i<pn; i++) {
372     /* add pdt[i,:]*AP into lnk */
373     pnz = pdti[i+1] - pdti[i];
374     ptJ = pdtj + pdti[i];
375     for (j=0; j<pnz; j++){
376       row  = ptJ[j];  /* row of AP == col of Pt */
377       apnz = api[row+1] - api[row];
378       Jptr = apj + api[row];
379       /* add non-zero cols of AP into the sorted linked list lnk */
380       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
381     }
382 
383     /* add received col data into lnk */
384     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
385       if (i == *nextrow[k]) { /* i-th row */
386         nzi = *(nextci[k]+1) - *nextci[k];
387         Jptr  = buf_rj[k] + *nextci[k];
388         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
389         nextrow[k]++; nextci[k]++;
390       }
391     }
392     nnz = lnk[0];
393 
394     /* if free space is not available, make more free space */
395     if (current_space->local_remaining<nnz) {
396       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
397       nspacedouble++;
398     }
399     /* copy data into free space, then initialize lnk */
400     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
401     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
402     current_space->array           += nnz;
403     current_space->local_used      += nnz;
404     current_space->local_remaining -= nnz;
405     bi[i+1] = bi[i] + nnz;
406     if (nnz > rmax) rmax = nnz;
407   }
408   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
409   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
410 
411   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
412   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
413   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]+1);
414   if (afill_tmp > afill) afill = afill_tmp;
415   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
416 
417   /* create symbolic parallel matrix Cmpi */
418   /*--------------------------------------*/
419   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
420   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
421   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
422   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
423   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
424   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
425 
426   merge->bi            = bi;
427   merge->bj            = bj;
428   merge->coi           = coi;
429   merge->coj           = coj;
430   merge->buf_ri        = buf_ri;
431   merge->buf_rj        = buf_rj;
432   merge->owners_co     = owners_co;
433   merge->destroy       = Cmpi->ops->destroy;
434   merge->duplicate     = Cmpi->ops->duplicate;
435 
436   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
437   Cmpi->assembled      = PETSC_FALSE;
438   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
439   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
440 
441   /* attach the supporting struct to Cmpi for reuse */
442   c = (Mat_MPIAIJ*)Cmpi->data;
443   c->ptap        = ptap;
444   ptap->api      = api;
445   ptap->apj      = apj;
446   ptap->merge    = merge;
447   ptap->rmax     = ap_rmax;
448   *C             = Cmpi;
449 
450   /* flag 'scalable' determines which implementations to be used:
451        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
452        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
453   /* set default scalable */
454   ptap->scalable = PETSC_TRUE;
455   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,PETSC_NULL);CHKERRQ(ierr);
456   if (!ptap->scalable){  /* Do dense axpy */
457     ierr = PetscMalloc(pN*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
458   } else {
459     ierr = PetscMalloc((ap_rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
460   }
461 
462 #if defined(PETSC_USE_INFO)
463   if (bi[pn] != 0) {
464     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
465     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
466   } else {
467     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
468   }
469 #endif
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
475 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
476 {
477   PetscErrorCode       ierr;
478   Mat_Merge_SeqsToMPI  *merge;
479   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
480   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
481   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
482   Mat_SeqAIJ           *p_loc,*p_oth;
483   Mat_PtAPMPI          *ptap;
484   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
485   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
486   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
487   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
488   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
489   MPI_Comm             comm=((PetscObject)C)->comm;
490   PetscMPIInt          size,rank,taga,*len_s;
491   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
492   PetscInt             **buf_ri,**buf_rj;
493   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
494   MPI_Request          *s_waits,*r_waits;
495   MPI_Status           *status;
496   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
497   PetscInt             *api,*apj,*coi,*coj;
498   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
499   PetscBool            scalable;
500 
501   PetscFunctionBegin;
502   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
503   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
504 
505   ptap = c->ptap;
506   if (!ptap) SETERRQ(((PetscObject)C)->comm,PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
507   merge    = ptap->merge;
508   apa      = ptap->apa;
509   scalable = ptap->scalable;
510 
511   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
512   /*--------------------------------------------------*/
513   if (ptap->reuse == MAT_INITIAL_MATRIX){
514     ptap->reuse = MAT_REUSE_MATRIX;
515   } else { /* update numerical values of P_oth and P_loc */
516     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
517     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
518   }
519 
520   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
521   /*--------------------------------------------------------------*/
522   /* get data from symbolic products */
523   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
524   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
525   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
526   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
527 
528   coi = merge->coi; coj = merge->coj;
529   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
530   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
531 
532   bi     = merge->bi; bj = merge->bj;
533   owners = merge->rowmap->range;
534   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
535   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
536 
537   api = ptap->api; apj = ptap->apj;
538 
539   if (!scalable){  /* Do dense axpy */
540     /*--------------------------------------------------*/
541     /* apa (length of pN) stores dense row A[i,:]*P - nonscalable! */
542     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
543 
544     for (i=0; i<am; i++) {
545       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
546       /*------------------------------------------------------------*/
547       apJ = apj + api[i];
548 
549       /* diagonal portion of A */
550       anz = adi[i+1] - adi[i];
551       adj = ad->j + adi[i];
552       ada = ad->a + adi[i];
553       for (j=0; j<anz; j++) {
554         row = adj[j];
555         pnz = pi_loc[row+1] - pi_loc[row];
556         pj  = pj_loc + pi_loc[row];
557         pa  = pa_loc + pi_loc[row];
558 
559         /* perform dense axpy */
560         valtmp = ada[j];
561         for (k=0; k<pnz; k++){
562           apa[pj[k]] += valtmp*pa[k];
563         }
564         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
565       }
566 
567       /* off-diagonal portion of A */
568       anz = aoi[i+1] - aoi[i];
569       aoj = ao->j + aoi[i];
570       aoa = ao->a + aoi[i];
571       for (j=0; j<anz; j++) {
572         row = aoj[j];
573         pnz = pi_oth[row+1] - pi_oth[row];
574         pj  = pj_oth + pi_oth[row];
575         pa  = pa_oth + pi_oth[row];
576 
577         /* perform dense axpy */
578         valtmp = aoa[j];
579         for (k=0; k<pnz; k++){
580           apa[pj[k]] += valtmp*pa[k];
581         }
582         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
583       }
584 
585       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
586       /*--------------------------------------------------------------*/
587       apnz = api[i+1] - api[i];
588       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
589       pnz = po->i[i+1] - po->i[i];
590       poJ = po->j + po->i[i];
591       pA  = po->a + po->i[i];
592       for (j=0; j<pnz; j++){
593         row = poJ[j];
594         cnz = coi[row+1] - coi[row];
595         cj  = coj + coi[row];
596         ca  = coa + coi[row];
597         /* perform dense axpy */
598         valtmp = pA[j];
599         for (k=0; k<cnz; k++) {
600           ca[k] += valtmp*apa[cj[k]];
601         }
602         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
603       }
604 
605       /* put the value into Cd (diagonal part) */
606       pnz = pd->i[i+1] - pd->i[i];
607       pdJ = pd->j + pd->i[i];
608       pA  = pd->a + pd->i[i];
609       for (j=0; j<pnz; j++){
610         row = pdJ[j];
611         cnz = bi[row+1] - bi[row];
612         cj  = bj + bi[row];
613         ca  = ba + bi[row];
614         /* perform dense axpy */
615         valtmp = pA[j];
616         for (k=0; k<cnz; k++) {
617           ca[k] += valtmp*apa[cj[k]];
618         }
619         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
620       }
621 
622       /* zero the current row of A*P */
623       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
624     }
625   } else {/* Perform sparse axpy */
626     /*----------------------------------------------------*/
627     /* apa (length ap_rmax) stores sparse row A[i,:]*P */
628     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
629 
630     pA=pa_loc;
631     for (i=0; i<am; i++) {
632       /* form i-th sparse row of A*P */
633       apnz = api[i+1] - api[i];
634       apJ  = apj + api[i];
635       /* diagonal portion of A */
636       anz = adi[i+1] - adi[i];
637       adj = ad->j + adi[i];
638       ada = ad->a + adi[i];
639       for (j=0; j<anz; j++) {
640         row = adj[j];
641         pnz = pi_loc[row+1] - pi_loc[row];
642         pj  = pj_loc + pi_loc[row];
643         pa  = pa_loc + pi_loc[row];
644         valtmp = ada[j];
645         nextp  = 0;
646         for (k=0; nextp<pnz; k++) {
647           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
648             apa[k] += valtmp*pa[nextp++];
649           }
650         }
651         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
652       }
653       /* off-diagonal portion of A */
654       anz = aoi[i+1] - aoi[i];
655       aoj = ao->j + aoi[i];
656       aoa = ao->a + aoi[i];
657       for (j=0; j<anz; j++) {
658         row = aoj[j];
659         pnz = pi_oth[row+1] - pi_oth[row];
660         pj  = pj_oth + pi_oth[row];
661         pa  = pa_oth + pi_oth[row];
662         valtmp = aoa[j];
663         nextp  = 0;
664         for (k=0; nextp<pnz; k++) {
665           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
666             apa[k] += valtmp*pa[nextp++];
667           }
668         }
669         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
670       }
671 
672       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
673       /*--------------------------------------------------------------*/
674       pnz = pi_loc[i+1] - pi_loc[i];
675       pJ  = pj_loc + pi_loc[i];
676       for (j=0; j<pnz; j++) {
677         nextap = 0;
678         row    = pJ[j]; /* global index */
679         if (row < pcstart || row >=pcend) { /* put the value into Co */
680           row = *poJ;
681           cj  = coj + coi[row];
682           ca  = coa + coi[row]; poJ++;
683         } else {                            /* put the value into Cd */
684           row  = *pdJ;
685           cj   = bj + bi[row];
686           ca   = ba + bi[row]; pdJ++;
687         }
688         valtmp = pA[j];
689         for (k=0; nextap<apnz; k++) {
690           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
691         }
692         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
693       }
694       pA += pnz;
695       /* zero the current row info for A*P */
696       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
697     }
698   }
699 
700   /* 3) send and recv matrix values coa */
701   /*------------------------------------*/
702   buf_ri = merge->buf_ri;
703   buf_rj = merge->buf_rj;
704   len_s  = merge->len_s;
705   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
706   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
707 
708   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
709   for (proc=0,k=0; proc<size; proc++){
710     if (!len_s[proc]) continue;
711     i = merge->owners_co[proc];
712     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
713     k++;
714   }
715   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
716   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
717 
718   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
719   ierr = PetscFree(r_waits);CHKERRQ(ierr);
720   ierr = PetscFree(coa);CHKERRQ(ierr);
721 
722   /* 4) insert local Cseq and received values into Cmpi */
723   /*------------------------------------------------------*/
724   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
725   for (k=0; k<merge->nrecv; k++){
726     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
727     nrows       = *(buf_ri_k[k]);
728     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
729     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
730   }
731 
732   for (i=0; i<cm; i++) {
733     row = owners[rank] + i; /* global row index of C_seq */
734     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
735     ba_i = ba + bi[i];
736     bnz  = bi[i+1] - bi[i];
737     /* add received vals into ba */
738     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
739       /* i-th row */
740       if (i == *nextrow[k]) {
741         cnz = *(nextci[k]+1) - *nextci[k];
742         cj  = buf_rj[k] + *(nextci[k]);
743         ca  = abuf_r[k] + *(nextci[k]);
744         nextcj = 0;
745         for (j=0; nextcj<cnz; j++){
746           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
747             ba_i[j] += ca[nextcj++];
748           }
749         }
750         nextrow[k]++; nextci[k]++;
751         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
752       }
753     }
754     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
755   }
756   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
757   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758 
759   ierr = PetscFree(ba);CHKERRQ(ierr);
760   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
761   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
762   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
763 #if defined(PETSC_USE_INFO)
764   if (scalable){
765     ierr = PetscInfo(C,"Use scalable sparse axpy\n");CHKERRQ(ierr);
766   } else {
767     ierr = PetscInfo(C,"Use non-scalable dense axpy\n");CHKERRQ(ierr);
768   }
769 #endif
770   PetscFunctionReturn(0);
771 }
772