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