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