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