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