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