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