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