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