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