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