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