xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 3dfa136d8b789e6905fd401804d07a78b7147c92)
1 #define PETSCMAT_DLL
2 
3 /*
4   Defines projective product routines where A is a MPIAIJ matrix
5           C = P^T * A * P
6 */
7 
8 #include "src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
9 #include "src/mat/utils/freespace.h"
10 #include "src/mat/impls/aij/mpi/mpiaij.h"
11 #include "petscbt.h"
12 
13 EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
14 #undef __FUNCT__
15 #define __FUNCT__ "MatDestroy_MPIAIJ_MatPtAP"
16 PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_MatPtAP(Mat A)
17 {
18   PetscErrorCode       ierr;
19   Mat_Merge_SeqsToMPI  *merge;
20   PetscObjectContainer container;
21 
22   PetscFunctionBegin;
23   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
24   if (container) {
25     ierr = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
26     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
27     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
28     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
29     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
30     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
31     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
32     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
33     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
34     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
35     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
36 
37     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
38     ierr = PetscObjectCompose((PetscObject)A,"MatMergeSeqsToMPI",0);CHKERRQ(ierr);
39   }
40   ierr = merge->MatDestroy(A);CHKERRQ(ierr);
41   ierr = PetscFree(merge);CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
47 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) {
48   PetscErrorCode       ierr;
49   Mat_Merge_SeqsToMPI  *merge;
50   PetscObjectContainer container;
51 
52   PetscFunctionBegin;
53   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
54   if (container) {
55     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
56   } else {
57     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
58   }
59   ierr = (*merge->MatDuplicate)(A,op,M);CHKERRQ(ierr);
60   (*M)->ops->destroy   = merge->MatDestroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
61   (*M)->ops->duplicate = merge->MatDuplicate; /* =MatDuplicate_ MPIAIJ */
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
67 PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
68 {
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   if (!P->ops->ptapsymbolic_mpiaij) {
73     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name);
74   }
75   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
81 PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
82 {
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   if (!P->ops->ptapnumeric_mpiaij) {
87     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name);
88   }
89   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
95 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
96 {
97   PetscErrorCode       ierr;
98   Mat                  B_mpi;
99   Mat_MatMatMultMPI    *ap;
100   PetscObjectContainer container;
101   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
102   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
103   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
104   Mat_SeqAIJ           *p_loc,*p_oth;
105   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
106   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
107   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
108   PetscInt             am=A->rmap.n,pN=P->cmap.N,pn=P->cmap.n;
109   PetscBT              lnkbt;
110   MPI_Comm             comm=A->comm;
111   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
112   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
113   PetscInt             len,proc,*dnz,*onz,*owners;
114   PetscInt             nzi,*bi,*bj;
115   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
116   MPI_Request          *swaits,*rwaits;
117   MPI_Status           *sstatus,rstatus;
118   Mat_Merge_SeqsToMPI  *merge;
119   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon;
120   PetscMPIInt          j;
121 
122   PetscFunctionBegin;
123   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
124   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
125 
126   /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */
127   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
128   if (container) {
129     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
130     ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
131   }
132 
133   /* create the container 'Mat_MatMatMultMPI' and attach it to P */
134   ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr);
135   ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
136   ap->abnz_max = 0;
137 
138   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
139   ierr = PetscObjectContainerSetPointer(container,ap);CHKERRQ(ierr);
140   ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
141   ap->MatDestroy  = P->ops->destroy;
142   P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
143   ap->reuse       = MAT_INITIAL_MATRIX;
144   ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
145 
146   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
147   ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
148   /* get P_loc by taking all local rows of P */
149   ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr);
150 
151   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
152   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
153   pi_loc = p_loc->i; pj_loc = p_loc->j;
154   pi_oth = p_oth->i; pj_oth = p_oth->j;
155 
156   /* first, compute symbolic AP = A_loc*P */
157   /*--------------------------------------*/
158   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
159   ap->abi = api;
160   api[0] = 0;
161 
162   /* create and initialize a linked list */
163   nlnk = pN+1;
164   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
165 
166   /* Initial FreeSpace size is fill*nnz(A) */
167   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
168   current_space = free_space;
169 
170   for (i=0;i<am;i++) {
171     apnz = 0;
172     /* diagonal portion of A */
173     nzi = adi[i+1] - adi[i];
174     for (j=0; j<nzi; j++){
175       row = *adj++;
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 = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
180       apnz += nlnk;
181     }
182     /* off-diagonal portion of A */
183     nzi = aoi[i+1] - aoi[i];
184     for (j=0; j<nzi; j++){
185       row = *aoj++;
186       pnz = pi_oth[row+1] - pi_oth[row];
187       Jptr  = pj_oth + pi_oth[row];
188       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
189       apnz += nlnk;
190     }
191 
192     api[i+1] = api[i] + apnz;
193     if (ap->abnz_max < apnz) ap->abnz_max = apnz;
194 
195     /* if free space is not available, double the total space in the list */
196     if (current_space->local_remaining<apnz) {
197       ierr = PetscFreeSpaceGet(current_space->total_array_size,&current_space);CHKERRQ(ierr);
198     }
199 
200     /* Copy data into free space, then initialize lnk */
201     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
202     current_space->array           += apnz;
203     current_space->local_used      += apnz;
204     current_space->local_remaining -= apnz;
205   }
206   /* Allocate space for apj, initialize apj, and */
207   /* destroy list of free space and other temporary array(s) */
208   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr);
209   apj = ap->abj;
210   ierr = PetscFreeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr);
211 
212   /* determine symbolic Co=(p->B)^T*AP - send to others */
213   /*----------------------------------------------------*/
214   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
215 
216   /* then, compute symbolic Co = (p->B)^T*AP */
217   pon = (p->B)->cmap.n; /* total num of rows to be sent to other processors
218                          >= (num of nonzero rows of C_seq) - pn */
219   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
220   coi[0] = 0;
221 
222   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
223   nnz           = 3*pon*ap->abnz_max + 1;
224   ierr          = PetscFreeSpaceGet(nnz,&free_space);
225   current_space = free_space;
226 
227   for (i=0; i<pon; i++) {
228     nnz  = 0;
229     pnz = poti[i+1] - poti[i];
230     j     = pnz;
231     ptJ   = potj + poti[i+1];
232     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
233       j--; ptJ--;
234       row  = *ptJ; /* 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 = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
239       nnz += nlnk;
240     }
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(current_space->total_array_size,&current_space);CHKERRQ(ierr);
245     }
246 
247     /* Copy data into free space, and zero out denserows */
248     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
249     current_space->array           += nnz;
250     current_space->local_used      += nnz;
251     current_space->local_remaining -= nnz;
252     coi[i+1] = coi[i] + nnz;
253   }
254   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
255   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
256   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
257 
258   /* send j-array (coj) of Co to other processors */
259   /*----------------------------------------------*/
260   /* determine row ownership */
261   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
262   merge->rowmap.n = pn;
263   merge->rowmap.N = PETSC_DECIDE;
264   ierr = PetscMapInitialize(comm,&merge->rowmap);CHKERRQ(ierr);
265   owners = merge->rowmap.range;
266 
267   /* determine the number of messages to send, their lengths */
268   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
269   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
270   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
271   len_s = merge->len_s;
272   merge->nsend = 0;
273 
274   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
275   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
276 
277   proc = 0;
278   for (i=0; i<pon; i++){
279     while (prmap[i] >= owners[proc+1]) proc++;
280     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
281     len_s[proc] += coi[i+1] - coi[i];
282   }
283 
284   len   = 0;  /* max length of buf_si[] */
285   owners_co[0] = 0;
286   for (proc=0; proc<size; proc++){
287     owners_co[proc+1] = owners_co[proc] + len_si[proc];
288     if (len_si[proc]){
289       merge->nsend++;
290       len_si[proc] = 2*(len_si[proc] + 1);
291       len += len_si[proc];
292     }
293   }
294 
295   /* determine the number and length of messages to receive for coi and coj  */
296   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
297   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
298 
299   /* post the Irecv and Isend of coj */
300   ierr = PetscObjectGetNewTag((PetscObject)merge,&tag);CHKERRQ(ierr);
301   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
302 
303   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
304 
305   for (proc=0, k=0; proc<size; proc++){
306     if (!len_s[proc]) continue;
307     i = owners_co[proc];
308     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
309     k++;
310   }
311 
312   /* receives and sends of coj are complete */
313   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
314   i = merge->nrecv;
315   while (i--) {
316     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&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   /* send and recv coi */
322   /*-------------------*/
323   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
324 
325   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&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       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
343       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
344       nrows++;
345     }
346     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
347     k++;
348     buf_si += len_si[proc];
349   }
350   i = merge->nrecv;
351   while (i--) {
352     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
353   }
354   ierr = PetscFree(rwaits);CHKERRQ(ierr);
355   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
356 
357   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
358   for (i=0; i<merge->nrecv; i++){
359     ierr = PetscInfo3(A,"recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr);
360   }
361 
362   ierr = PetscFree(len_si);CHKERRQ(ierr);
363   ierr = PetscFree(len_ri);CHKERRQ(ierr);
364   ierr = PetscFree(swaits);CHKERRQ(ierr);
365   ierr = PetscFree(sstatus);CHKERRQ(ierr);
366   ierr = PetscFree(buf_s);CHKERRQ(ierr);
367 
368   /* compute the local portion of C (mpi mat) */
369   /*------------------------------------------*/
370   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
371 
372   /* allocate bi array and free space for accumulating nonzero column info */
373   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
374   bi[0] = 0;
375 
376   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
377   nnz           = 3*pn*ap->abnz_max + 1;
378   ierr          = PetscFreeSpaceGet(nnz,&free_space);
379   current_space = free_space;
380 
381   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
382   nextrow = buf_ri_k + merge->nrecv;
383   nextci  = nextrow + merge->nrecv;
384   for (k=0; k<merge->nrecv; k++){
385     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
386     nrows       = *buf_ri_k[k];
387     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
388     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
389   }
390   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
391   for (i=0; i<pn; i++) {
392     /* add pdt[i,:]*AP into lnk */
393     nnz = 0;
394     pnz  = pdti[i+1] - pdti[i];
395     j    = pnz;
396     ptJ  = pdtj + pdti[i+1];
397     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
398       j--; ptJ--;
399       row  = *ptJ; /* 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 = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
404       nnz += nlnk;
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 = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
412         nnz += nlnk;
413         nextrow[k]++; nextci[k]++;
414       }
415     }
416 
417     /* if free space is not available, make more free space */
418     if (current_space->local_remaining<nnz) {
419       ierr = PetscFreeSpaceGet(current_space->total_array_size,&current_space);CHKERRQ(ierr);
420     }
421     /* copy data into free space, then initialize lnk */
422     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
423     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
424     current_space->array           += nnz;
425     current_space->local_used      += nnz;
426     current_space->local_remaining -= nnz;
427     bi[i+1] = bi[i] + nnz;
428   }
429   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
430   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
431 
432   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
433   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
434   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
435 
436   /* create symbolic parallel matrix B_mpi */
437   /*---------------------------------------*/
438   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
439   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
440   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
441   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
442   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
443 
444   merge->bi            = bi;
445   merge->bj            = bj;
446   merge->coi           = coi;
447   merge->coj           = coj;
448   merge->buf_ri        = buf_ri;
449   merge->buf_rj        = buf_rj;
450   merge->owners_co     = owners_co;
451   merge->MatDestroy    = B_mpi->ops->destroy;
452   merge->MatDuplicate  = B_mpi->ops->duplicate;
453 
454   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
455   B_mpi->assembled     = PETSC_FALSE;
456   B_mpi->ops->destroy  = MatDestroy_MPIAIJ_MatPtAP;
457   B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
458 
459   /* attach the supporting struct to B_mpi for reuse */
460   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
461   ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr);
462   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
463   *C = B_mpi;
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
469 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
470 {
471   PetscErrorCode       ierr;
472   Mat_Merge_SeqsToMPI  *merge;
473   Mat_MatMatMultMPI    *ap;
474   PetscObjectContainer cont_merge,cont_ptap;
475   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
476   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
477   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
478   Mat_SeqAIJ           *p_loc,*p_oth;
479   PetscInt             *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0;
480   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
481   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
482   MatScalar            *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
483   PetscInt             am=A->rmap.n,cm=C->rmap.n,pon=(p->B)->cmap.n;
484   MPI_Comm             comm=C->comm;
485   PetscMPIInt          size,rank,taga,*len_s;
486   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
487   PetscInt             **buf_ri,**buf_rj;
488   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
489   MPI_Request          *s_waits,*r_waits;
490   MPI_Status           *status;
491   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
492   PetscInt             *api,*apj,*coi,*coj;
493   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap.rstart,pcend=P->cmap.rend;
494 
495   PetscFunctionBegin;
496   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
497   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
498 
499   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
500   if (cont_merge) {
501     ierr  = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
502   } else {
503     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
504   }
505 
506   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
507   if (cont_ptap) {
508     ierr  = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr);
509     if (ap->reuse == MAT_INITIAL_MATRIX){
510       ap->reuse = MAT_REUSE_MATRIX;
511     } else { /* update numerical values of P_oth and P_loc */
512       ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
513       ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr);
514     }
515   } else {
516     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
517   }
518 
519   /* get data from symbolic products */
520   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
521   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
522   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
523   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
524 
525   coi = merge->coi; coj = merge->coj;
526   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
527   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
528 
529   bi     = merge->bi; bj = merge->bj;
530   owners = merge->rowmap.range;
531   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
532   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
533 
534   /* get data from symbolic A*P */
535   ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
536   ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
537 
538   /* compute numeric C_seq=P_loc^T*A_loc*P */
539   api = ap->abi; apj = ap->abj;
540   for (i=0;i<am;i++) {
541     /* form i-th sparse row of A*P */
542     apnz = api[i+1] - api[i];
543     apJ  = apj + api[i];
544     /* diagonal portion of A */
545     anz  = adi[i+1] - adi[i];
546     for (j=0;j<anz;j++) {
547       row = *adj++;
548       pnz = pi_loc[row+1] - pi_loc[row];
549       pj  = pj_loc + pi_loc[row];
550       pa  = pa_loc + pi_loc[row];
551       nextp = 0;
552       for (k=0; nextp<pnz; k++) {
553         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
554           apa[k] += (*ada)*pa[nextp++];
555         }
556       }
557       flops += 2*pnz;
558       ada++;
559     }
560     /* off-diagonal portion of A */
561     anz  = aoi[i+1] - aoi[i];
562     for (j=0; j<anz; j++) {
563       row = *aoj++;
564       pnz = pi_oth[row+1] - pi_oth[row];
565       pj  = pj_oth + pi_oth[row];
566       pa  = pa_oth + pi_oth[row];
567       nextp = 0;
568       for (k=0; nextp<pnz; k++) {
569         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
570           apa[k] += (*aoa)*pa[nextp++];
571         }
572       }
573       flops += 2*pnz;
574       aoa++;
575     }
576 
577     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
578     pnz = pi_loc[i+1] - pi_loc[i];
579     for (j=0; j<pnz; j++) {
580       nextap = 0;
581       row    = *pJ++; /* global index */
582       if (row < pcstart || row >=pcend) { /* put the value into Co */
583         cj  = coj + coi[*poJ];
584         ca  = coa + coi[*poJ++];
585       } else {                            /* put the value into Cd */
586         cj   = bj + bi[*pdJ];
587         ca   = ba + bi[*pdJ++];
588       }
589       for (k=0; nextap<apnz; k++) {
590         if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
591       }
592       flops += 2*apnz;
593       pA++;
594     }
595 
596     /* zero the current row info for A*P */
597     ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
598   }
599   ierr = PetscFree(apa);CHKERRQ(ierr);
600 
601   /* send and recv matrix values */
602   /*-----------------------------*/
603   buf_ri = merge->buf_ri;
604   buf_rj = merge->buf_rj;
605   len_s  = merge->len_s;
606   ierr = PetscObjectGetNewTag((PetscObject)merge,&taga);CHKERRQ(ierr);
607   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
608 
609   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
610   for (proc=0,k=0; proc<size; proc++){
611     if (!len_s[proc]) continue;
612     i = merge->owners_co[proc];
613     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
614     k++;
615   }
616   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
617   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
618   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
619   ierr = PetscFree(status);CHKERRQ(ierr);
620 
621   ierr = PetscFree(s_waits);CHKERRQ(ierr);
622   ierr = PetscFree(r_waits);CHKERRQ(ierr);
623   ierr = PetscFree(coa);CHKERRQ(ierr);
624 
625   /* insert local and received values into C */
626   /*-----------------------------------------*/
627   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
628   nextrow   = buf_ri_k + merge->nrecv;
629   nextci = nextrow + merge->nrecv;
630 
631   for (k=0; k<merge->nrecv; k++){
632     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
633     nrows       = *(buf_ri_k[k]);
634     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
635     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
636   }
637 
638   for (i=0; i<cm; i++) {
639     row = owners[rank] + i; /* global row index of C_seq */
640     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
641     ba_i = ba + bi[i];
642     bnz  = bi[i+1] - bi[i];
643     /* add received vals into ba */
644     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
645       /* i-th row */
646       if (i == *nextrow[k]) {
647         cnz = *(nextci[k]+1) - *nextci[k];
648         cj  = buf_rj[k] + *(nextci[k]);
649         ca  = abuf_r[k] + *(nextci[k]);
650         nextcj = 0;
651         for (j=0; nextcj<cnz; j++){
652           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
653             ba_i[j] += ca[nextcj++];
654           }
655         }
656         nextrow[k]++; nextci[k]++;
657       }
658     }
659     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
660     flops += 2*cnz;
661   }
662   ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */
663   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
664   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
665 
666   ierr = PetscFree(ba);CHKERRQ(ierr);
667   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
668   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
669   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
670   PetscFunctionReturn(0);
671 }
672