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